MultisegmentWellPrimaryVariables.hpp
Go to the documentation of this file.
1/*
2 Copyright 2017 SINTEF Digital, Mathematics and Cybernetics.
3 Copyright 2017 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
21
22#ifndef OPM_MULTISEGMENTWELL_PRIMARY_VARIABLES_HEADER_INCLUDED
23#define OPM_MULTISEGMENTWELL_PRIMARY_VARIABLES_HEADER_INCLUDED
24
25#include <opm/material/densead/Evaluation.hpp>
26
28#include <opm/input/eclipse/Schedule/SummaryState.hpp>
29#include <opm/input/eclipse/Units/Units.hpp>
30
31#include <array>
32#include <cstddef>
33#include <vector>
34
35namespace Opm
36{
37
38class DeferredLogger;
39template <typename Scalar, typename IndexTraits> class MultisegmentWellGeneric;
40template<class FluidSystem, class Indices> class WellInterfaceIndices;
41template<typename Scalar, typename IndexTraits> class WellState;
42
43template<class FluidSystem, class Indices>
45{
46public:
47 // TODO: for now, not considering the polymer, solvent and so on to simplify the development process.
48
49 // TODO: we need to have order for the primary variables and also the order for the well equations.
50 // sometimes, they are similar, while sometimes, they can have very different forms.
51
52 // Table showing the primary variable indices, depending on what phases are present:
53 //
54 // WOG OG WG WO W/O/G (single phase)
55 // WQTotal 0 0 0 0 0
56 // WFrac 1 -1000 -1000 1 -1000
57 // GFrac 2 1 1 -1000 -1000
58 // Spres 3 2 2 2 1
59
60 static constexpr bool has_wfrac_variable = Indices::waterEnabled && Indices::oilEnabled;
61 static constexpr bool has_gfrac_variable = Indices::gasEnabled && Indices::numPhases > 1;
62
63 static constexpr int WQTotal = 0;
64 static constexpr int WFrac = has_wfrac_variable ? 1 : -1000;
65 static constexpr int GFrac = has_gfrac_variable ? has_wfrac_variable + 1 : -1000;
66 static constexpr int SPres = has_wfrac_variable + has_gfrac_variable + 1;
67
68 // the number of well equations TODO: it should have a more general strategy for it
69 static constexpr int numWellEq = Indices::numPhases + 1;
70
71 using Scalar = typename FluidSystem::Scalar;
72 using IndexTraits = typename FluidSystem::IndexTraitsType;
73 using EvalWell = DenseAd::Evaluation<Scalar, /*size=*/Indices::numEq + numWellEq>;
74
77
79 : well_(well)
80 {}
81
83 void resize(const int numSegments);
84
86 void update(const WellState<Scalar, IndexTraits>& well_state,
87 const bool stop_or_zero_rate_target);
88
90 void updateNewton(const BVectorWell& dwells,
91 const Scalar relaxation_factor,
92 const Scalar DFLimit,
93 const bool stop_or_zero_rate_target,
94 const Scalar max_pressure_change);
95
98 const Scalar rho,
100 const SummaryState& summary_state,
101 DeferredLogger& deferred_logger) const;
102
106 const int compIdx) const;
107
111 const int compIdx) const;
112
115 const int seg_upwind,
116 const int comp_idx) const;
117
120
122 EvalWell getSegmentPressure(const int seg) const;
123
126 const int comp_idx) const;
127
129 EvalWell getQs(const int comp_idx) const;
130
133
135 const std::array<EvalWell, numWellEq>& eval(const int idx) const
136 { return evaluation_[idx]; }
137
139 const std::array<Scalar, numWellEq>& value(const int idx) const
140 { return value_[idx]; }
141
143 void setValue(const int idx, const std::array<Scalar, numWellEq>& val)
144 { value_[idx] = val; }
145
147 void outputLowLimitPressureSegments(DeferredLogger& deferred_logger) const;
148
149private:
151 void setEvaluationsFromValues();
152
154 void processFractions(const int seg);
155
157 EvalWell volumeFraction(const int seg,
158 const int compIdx) const;
159
162 std::vector<std::array<Scalar, numWellEq>> value_;
163
166 std::vector<std::array<EvalWell, numWellEq>> evaluation_;
167
169
170 static constexpr double bhp_lower_limit = 1. * unit::barsa - 1. * unit::Pascal;
171 static constexpr double seg_pres_lower_limit = 0.;
172};
173
174}
175
176#endif // OPM_MULTISEGMENTWELL_PRIMARY_VARIABLES_HEADER_INCLUDED
Definition: DeferredLogger.hpp:57
Dune::BlockVector< VectorBlockWellType > BVectorWell
Definition: MultisegmentWellEquations.hpp:60
Definition: MultisegmentWellGeneric.hpp:39
Definition: MultisegmentWellPrimaryVariables.hpp:45
static constexpr int WFrac
Definition: MultisegmentWellPrimaryVariables.hpp:64
void resize(const int numSegments)
Resize values and evaluations.
void copyToWellState(const MultisegmentWellGeneric< Scalar, IndexTraits > &mswell, const Scalar rho, WellState< Scalar, IndexTraits > &well_state, const SummaryState &summary_state, DeferredLogger &deferred_logger) const
Copy values to well state.
EvalWell getBhp() const
Get bottomhole pressure.
DenseAd::Evaluation< Scalar, Indices::numEq+numWellEq > EvalWell
Definition: MultisegmentWellPrimaryVariables.hpp:73
static constexpr int numWellEq
Definition: MultisegmentWellPrimaryVariables.hpp:69
static constexpr bool has_gfrac_variable
Definition: MultisegmentWellPrimaryVariables.hpp:61
void update(const WellState< Scalar, IndexTraits > &well_state, const bool stop_or_zero_rate_target)
Copy values from well state.
EvalWell getQs(const int comp_idx) const
Returns scaled rate for a component.
void outputLowLimitPressureSegments(DeferredLogger &deferred_logger) const
output the segments with pressure close to lower pressure limit for debugging purpose
EvalWell getSegmentRateUpwinding(const int seg, const int seg_upwind, const int comp_idx) const
Returns upwinding rate for a component in a segment.
typename FluidSystem::Scalar Scalar
Definition: MultisegmentWellPrimaryVariables.hpp:71
EvalWell getWQTotal() const
Get WQTotal.
EvalWell getSegmentPressure(const int seg) const
Get pressure for a segment.
const std::array< EvalWell, numWellEq > & eval(const int idx) const
Returns a const ref to an array of evaluations.
Definition: MultisegmentWellPrimaryVariables.hpp:135
typename FluidSystem::IndexTraitsType IndexTraits
Definition: MultisegmentWellPrimaryVariables.hpp:72
static constexpr bool has_wfrac_variable
Definition: MultisegmentWellPrimaryVariables.hpp:60
EvalWell volumeFractionScaled(const int seg, const int compIdx) const
Returns scaled volume fraction for a component in a segment.
MultisegmentWellPrimaryVariables(const WellInterfaceIndices< FluidSystem, Indices > &well)
Definition: MultisegmentWellPrimaryVariables.hpp:78
void setValue(const int idx, const std::array< Scalar, numWellEq > &val)
Set a value array. Note that this does not also set the corresponding evaluation.
Definition: MultisegmentWellPrimaryVariables.hpp:143
static constexpr int SPres
Definition: MultisegmentWellPrimaryVariables.hpp:66
static constexpr int GFrac
Definition: MultisegmentWellPrimaryVariables.hpp:65
typename Equations::BVectorWell BVectorWell
Definition: MultisegmentWellPrimaryVariables.hpp:76
void updateNewton(const BVectorWell &dwells, const Scalar relaxation_factor, const Scalar DFLimit, const bool stop_or_zero_rate_target, const Scalar max_pressure_change)
Update values from newton update vector.
static constexpr int WQTotal
Definition: MultisegmentWellPrimaryVariables.hpp:63
EvalWell getSegmentRate(const int seg, const int comp_idx) const
Get rate for a component in a segment.
const std::array< Scalar, numWellEq > & value(const int idx) const
Returns a value array.
Definition: MultisegmentWellPrimaryVariables.hpp:139
EvalWell surfaceVolumeFraction(const int seg, const int compIdx) const
Returns surface volume fraction for a component in a segment.
Definition: WellInterfaceIndices.hpp:34
Definition: WellState.hpp:66
Definition: blackoilboundaryratevector.hh:39