SimulatorFullyImplicitCompressiblePolymer_impl.hpp
Go to the documentation of this file.
1/*
2 Copyright 2014 SINTEF ICT, Applied Mathematics.
3 Copyright 2014 STATOIL ASA.
4
5 This file is part of the Open Porous Media project (OPM).
6
7 OPM is free software: you can redistribute it and/or modify
8 it under the terms of the GNU General Public License as published by
9 the Free Software Foundation, either version 3 of the License, or
10 (at your option) any later version.
11
12 OPM is distributed in the hope that it will be useful,
13 but WITHOUT ANY WARRANTY; without even the implied warranty of
14 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
15 GNU General Public License for more details.
16
17 You should have received a copy of the GNU General Public License
18 along with OPM. If not, see <http://www.gnu.org/licenses/>.
19*/
20#ifndef OPM_SIMULATORFULLYIMPLICITCOMPRESSIBLEPOLYMER_IMPL_HEADER_INCLUDED
21#define OPM_SIMULATORFULLYIMPLICITCOMPRESSIBLEPOLYMER_IMPL_HEADER_INCLUDED
22
23namespace Opm
24{
25
27template <class GridT>
29SimulatorFullyImplicitCompressiblePolymer(const parameter::ParameterGroup& param,
30 const GridT& grid,
31 DerivedGeology& geo,
32 BlackoilPropsAdInterface& props,
33 const PolymerPropsAd& polymer_props,
34 const RockCompressibility* rock_comp_props,
35 std::shared_ptr<EclipseState> eclipse_state,
36 BlackoilOutputWriter& output_writer,
37 Opm::DeckConstPtr& deck,
38 NewtonIterationBlackoilInterface& linsolver,
39 const double* gravity)
40: BaseType(param,
41 grid,
42 geo,
43 props,
44 rock_comp_props,
45 linsolver,
46 gravity,
47 /*disgas=*/false,
48 /*vapoil=*/false,
49 eclipse_state,
50 output_writer,
51 /*threshold_pressures_by_face=*/std::vector<double>())
52 , deck_(deck)
53 , polymer_props_(polymer_props)
54
55{
56}
57
58template <class GridT>
60createSolver(const Wells* wells)
61 -> std::unique_ptr<Solver>
62{
63 return std::unique_ptr<Solver>(new Solver(BaseType::grid_,
64 BaseType::props_,
65 BaseType::geo_,
66 BaseType::rock_comp_props_,
67 polymer_props_,
68 *wells,
69 BaseType::solver_));
70}
71
72template <class GridT>
74handleAdditionalWellInflow(SimulatorTimer& timer,
75 WellsManager& wells_manager,
76 typename BaseType::WellState& well_state,
77 const Wells* wells)
78{
79 // compute polymer inflow
80 std::unique_ptr<PolymerInflowInterface> polymer_inflow_ptr;
81 if (deck_->hasKeyword("WPOLYMER")) {
82 if (wells_manager.c_wells() == 0) {
83 OPM_THROW(std::runtime_error, "Cannot control polymer injection via WPOLYMER without wells.");
84 }
85 polymer_inflow_ptr.reset(new PolymerInflowFromDeck(deck_, BaseType::eclipse_state_, *wells, Opm::UgGridHelpers::numCells(BaseType::grid_), timer.currentStepNum()));
86 } else {
87 polymer_inflow_ptr.reset(new PolymerInflowBasic(0.0*Opm::unit::day,
88 1.0*Opm::unit::day,
89 0.0));
90 }
91 std::vector<double> polymer_inflow_c(Opm::UgGridHelpers::numCells(BaseType::grid_));
92 polymer_inflow_ptr->getInflowValues(timer.simulationTimeElapsed(),
93 timer.simulationTimeElapsed() + timer.currentStepLength(),
94 polymer_inflow_c);
95 well_state.polymerInflow() = polymer_inflow_c;
96}
97
98
99} // namespace Opm
100
101#endif // OPM_SIMULATORFULLYIMPLICITCOMPRESSIBLEPOLYMER_HEADER_INCLUDED
Basic polymer injection behaviour class. This class gives all injectors the same polymer concentratio...
Definition: PolymerInflow.hpp:60
Polymer injection behaviour class using deck WPOLYMER. This class reads the accumulated WPOLYMER line...
Definition: PolymerInflow.hpp:89
Definition: PolymerPropsAd.hpp:33
void handleAdditionalWellInflow(SimulatorTimer &timer, WellsManager &wells_manager, typename BaseType::WellState &well_state, const Wells *wells)
Definition: SimulatorFullyImplicitCompressiblePolymer_impl.hpp:74
SimulatorFullyImplicitCompressiblePolymer(const parameter::ParameterGroup &param, const GridT &grid, DerivedGeology &geo, BlackoilPropsAdInterface &props, const PolymerPropsAd &polymer_props, const RockCompressibility *rock_comp_props, std::shared_ptr< EclipseState > eclipse_state, BlackoilOutputWriter &output_writer, Opm::DeckConstPtr &deck, NewtonIterationBlackoilInterface &linsolver, const double *gravity)
Initialise from parameters and objects to observe.
Definition: SimulatorFullyImplicitCompressiblePolymer_impl.hpp:29
std::unique_ptr< Solver > createSolver(const Wells *wells)
Definition: SimulatorFullyImplicitCompressiblePolymer_impl.hpp:60
Definition: CompressibleTpfaPolymer.hpp:33