2012-10-04 08:15:32 -05:00
|
|
|
/*
|
|
|
|
Copyright 2012 SINTEF ICT, Applied Mathematics.
|
|
|
|
|
|
|
|
This file is part of the Open Porous Media project (OPM).
|
|
|
|
|
|
|
|
OPM is free software: you can redistribute it and/or modify
|
|
|
|
it under the terms of the GNU General Public License as published by
|
|
|
|
the Free Software Foundation, either version 3 of the License, or
|
|
|
|
(at your option) any later version.
|
|
|
|
|
|
|
|
OPM is distributed in the hope that it will be useful,
|
|
|
|
but WITHOUT ANY WARRANTY; without even the implied warranty of
|
|
|
|
MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
|
|
|
|
GNU General Public License for more details.
|
|
|
|
|
|
|
|
You should have received a copy of the GNU General Public License
|
|
|
|
along with OPM. If not, see <http://www.gnu.org/licenses/>.
|
|
|
|
*/
|
|
|
|
|
2014-11-18 18:44:54 -06:00
|
|
|
#include <config.h>
|
|
|
|
|
2012-10-04 08:15:32 -05:00
|
|
|
#include <opm/polymer/PolymerInflow.hpp>
|
2013-03-18 07:10:32 -05:00
|
|
|
#include <opm/core/wells.h>
|
2014-04-26 04:09:43 -05:00
|
|
|
#include <opm/parser/eclipse/Deck/Deck.hpp>
|
2014-11-20 02:32:02 -06:00
|
|
|
#include <opm/parser/eclipse/EclipseState/Schedule/Well.hpp>
|
|
|
|
#include <opm/parser/eclipse/EclipseState/EclipseState.hpp>
|
|
|
|
#include <opm/parser/eclipse/EclipseState/Schedule/Schedule.hpp>
|
|
|
|
#include <opm/parser/eclipse/EclipseState/Schedule/ScheduleEnums.hpp>
|
|
|
|
#include <opm/parser/eclipse/EclipseState/Schedule/WellPolymerProperties.hpp>
|
|
|
|
#include <opm/parser/eclipse/EclipseState/Schedule/WellInjectionProperties.hpp>
|
2012-10-04 08:15:32 -05:00
|
|
|
#include <map>
|
2014-11-20 02:32:02 -06:00
|
|
|
#include <unordered_map>
|
|
|
|
#include <memory>
|
2012-10-04 08:15:32 -05:00
|
|
|
namespace Opm
|
|
|
|
{
|
|
|
|
|
|
|
|
// ---------- Methods of PolymerInflowBasic ----------
|
|
|
|
|
|
|
|
/// Constructor.
|
|
|
|
/// @param[in] starttime Start time of injection in seconds.
|
|
|
|
/// @param[in] endtime End time of injection in seconds.
|
|
|
|
/// @param[in] amount Amount to be injected per second.
|
|
|
|
PolymerInflowBasic::PolymerInflowBasic(const double starttime,
|
|
|
|
const double endtime,
|
|
|
|
const double amount)
|
|
|
|
: stime_(starttime), etime_(endtime), amount_(amount)
|
|
|
|
{
|
|
|
|
}
|
|
|
|
|
|
|
|
void PolymerInflowBasic::getInflowValues(const double step_start,
|
|
|
|
const double step_end,
|
2012-10-04 08:59:15 -05:00
|
|
|
std::vector<double>& poly_inflow_c) const
|
2012-10-04 08:15:32 -05:00
|
|
|
{
|
|
|
|
const double eps = 1e-5*(step_end - step_start);
|
|
|
|
if (step_start + eps >= stime_ && step_end - eps <= etime_) {
|
|
|
|
std::fill(poly_inflow_c.begin(), poly_inflow_c.end(), amount_);
|
|
|
|
} else if (step_start + eps <= etime_ && step_end - eps >= stime_) {
|
2013-09-03 08:00:30 -05:00
|
|
|
OPM_MESSAGE("Warning: polymer injection set to change inside timestep. Using value at start of step.");
|
2012-10-04 08:15:32 -05:00
|
|
|
std::fill(poly_inflow_c.begin(), poly_inflow_c.end(), amount_);
|
|
|
|
} else {
|
|
|
|
std::fill(poly_inflow_c.begin(), poly_inflow_c.end(), 0.0);
|
|
|
|
}
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
|
|
// ---------- Methods of PolymerInflowFromDeck ----------
|
|
|
|
|
|
|
|
|
2014-11-20 02:32:02 -06:00
|
|
|
void
|
2016-01-06 07:52:33 -06:00
|
|
|
PolymerInflowFromDeck::setInflowValues(Opm::EclipseStateConstPtr eclipseState,
|
2014-11-20 02:32:02 -06:00
|
|
|
size_t currentStep)
|
|
|
|
{
|
2015-01-26 22:24:26 -06:00
|
|
|
ScheduleConstPtr schedule = eclipseState->getSchedule();
|
2016-01-06 07:52:33 -06:00
|
|
|
for (const auto& well : schedule->getWells(currentStep)) {
|
2016-08-02 22:10:37 -05:00
|
|
|
if (well->isInjector(currentStep)) {
|
|
|
|
WellInjectionProperties injection = well->getInjectionProperties(currentStep);
|
|
|
|
if (injection.injectorType == WellInjector::WATER) {
|
|
|
|
WellPolymerProperties polymer = well->getPolymerProperties(currentStep);
|
|
|
|
wellPolymerRate_.insert(std::make_pair(well->name(), polymer.m_polymerConcentration));
|
|
|
|
} else {
|
|
|
|
OPM_THROW(std::logic_error, "For polymer injector you must have a water injector");
|
|
|
|
}
|
2014-11-20 02:32:02 -06:00
|
|
|
}
|
|
|
|
}
|
|
|
|
}
|
|
|
|
|
|
|
|
/// Constructor.
|
|
|
|
/// @param[in] deck Input deck expected to contain WPOLYMER.
|
2016-01-06 07:52:33 -06:00
|
|
|
PolymerInflowFromDeck::PolymerInflowFromDeck(Opm::EclipseStateConstPtr eclipseState,
|
2014-11-20 02:32:02 -06:00
|
|
|
const Wells& wells,
|
|
|
|
const int num_cells,
|
|
|
|
size_t currentStep)
|
|
|
|
: sparse_inflow_(num_cells)
|
|
|
|
{
|
2016-01-06 07:52:33 -06:00
|
|
|
setInflowValues(eclipseState, currentStep);
|
2014-11-20 02:32:02 -06:00
|
|
|
|
|
|
|
std::unordered_map<std::string, double>::const_iterator map_it;
|
|
|
|
// Extract concentrations and put into cell->concentration map.
|
|
|
|
std::map<int, double> perfcell_conc;
|
2015-05-20 01:57:26 -05:00
|
|
|
for (map_it = wellPolymerRate_.begin(); map_it != wellPolymerRate_.end(); ++map_it) {
|
2014-11-20 02:32:02 -06:00
|
|
|
// Only use well name and polymer concentration.
|
|
|
|
// That is, we ignore salt concentration and group
|
|
|
|
// names.
|
|
|
|
int wix = 0;
|
|
|
|
for (; wix < wells.number_of_wells; ++wix) {
|
2015-05-20 01:57:26 -05:00
|
|
|
if (wellPolymerRate_.count(wells.name[wix]) > 0) {
|
2014-11-20 02:32:02 -06:00
|
|
|
break;
|
|
|
|
}
|
|
|
|
}
|
2015-05-19 23:35:04 -05:00
|
|
|
if (wix == wells.number_of_wells) {
|
2015-05-20 01:57:26 -05:00
|
|
|
OPM_THROW(std::runtime_error, "Could not find a match for well "
|
|
|
|
<< map_it->first
|
|
|
|
<< " from WPOLYMER.");
|
|
|
|
|
2015-05-19 23:35:04 -05:00
|
|
|
}
|
2014-11-20 02:32:02 -06:00
|
|
|
for (int j = wells.well_connpos[wix]; j < wells.well_connpos[wix+1]; ++j) {
|
|
|
|
const int perf_cell = wells.well_cells[j];
|
|
|
|
perfcell_conc[perf_cell] = map_it->second;
|
|
|
|
}
|
|
|
|
}
|
|
|
|
|
|
|
|
// Build sparse vector from map.
|
|
|
|
std::map<int, double>::const_iterator it = perfcell_conc.begin();
|
|
|
|
for (; it != perfcell_conc.end(); ++it) {
|
|
|
|
sparse_inflow_.addElement(it->second, it->first);
|
|
|
|
}
|
|
|
|
}
|
2012-10-04 08:15:32 -05:00
|
|
|
|
|
|
|
void PolymerInflowFromDeck::getInflowValues(const double /*step_start*/,
|
|
|
|
const double /*step_end*/,
|
2012-10-04 08:59:15 -05:00
|
|
|
std::vector<double>& poly_inflow_c) const
|
2012-10-04 08:15:32 -05:00
|
|
|
{
|
|
|
|
// This method does not depend on the given time,
|
|
|
|
// instead one would have a new epoch (and create a new
|
|
|
|
// instance) for each change in WPOLYMER.
|
|
|
|
std::fill(poly_inflow_c.begin(), poly_inflow_c.end(), 0.0);
|
|
|
|
const int nnz = sparse_inflow_.nonzeroSize();
|
|
|
|
for (int i = 0; i < nnz; ++i) {
|
|
|
|
poly_inflow_c[sparse_inflow_.nonzeroIndex(i)] = sparse_inflow_.nonzeroElement(i) ;
|
|
|
|
}
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
|
|
} // namespace Opm
|