MultisegmentWellEval.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_EVAL_HEADER_INCLUDED
23#define OPM_MULTISEGMENTWELL_EVAL_HEADER_INCLUDED
24
29
30#include <opm/material/densead/Evaluation.hpp>
31
32#include <utility>
33#include <vector>
34
35namespace Opm
36{
37
38class ConvergenceReport;
39class Schedule;
40class SummaryState;
41
42template<class FluidSystem, class Indices> class WellInterfaceIndices;
43template<class Scalar> class WellState;
44
45template<typename FluidSystem, typename Indices>
46class MultisegmentWellEval : public MultisegmentWellGeneric<typename FluidSystem::Scalar>
47{
48protected:
49 using Scalar = typename FluidSystem::Scalar;
51 static constexpr int numWellEq = PrimaryVariables::numWellEq;
52 static constexpr int SPres = PrimaryVariables::SPres;
53 static constexpr int WQTotal = PrimaryVariables::WQTotal;
54
57
58 using BVector = typename Equations::BVector;
60
61 // TODO: for more efficient implementation, we should have EvalReservoir, EvalWell, and EvalRerservoirAndWell
62 // EvalR (Eval), EvalW, EvalRW
63 // TODO: for now, we only use one type to save some implementation efforts, while improve later.
65 using Eval = DenseAd::Evaluation<Scalar, /*size=*/Indices::numEq>;
66
67public:
69 const Equations& linSys() const
70 { return linSys_; }
71
72protected:
74
75 void initMatrixAndVectors(const int num_cells);
76
77 void assembleDefaultPressureEq(const int seg,
78 WellState<Scalar>& well_state,
79 const bool use_average_density);
80
81 // assemble pressure equation for ICD segments
82 void assembleICDPressureEq(const int seg,
83 const UnitSystem& unit_system,
84 WellState<Scalar>& well_state,
85 const SummaryState& summary_state,
86 const bool use_average_density,
87 DeferredLogger& deferred_logger);
88
90 WellState<Scalar>& well_state,
91 const bool use_average_density);
92
93
94 void assemblePressureEq(const int seg,
95 const UnitSystem& unit_system,
96 WellState<Scalar>& well_state,
97 const SummaryState& summary_state,
98 const bool use_average_density,
99 DeferredLogger& deferred_logger);
100
103 const std::vector<double>& B_avg,
104 DeferredLogger& deferred_logger,
105 const double max_residual_allowed,
106 const double tolerance_wells,
107 const double relaxed_inner_tolerance_flow_ms_well,
108 const double tolerance_pressure_ms_wells,
109 const double relaxed_inner_tolerance_pressure_ms_well,
110 const bool relax_tolerance,
111 const bool well_is_stopped) const;
112
113 std::pair<bool, std::vector<Scalar> >
114 getFiniteWellResiduals(const std::vector<Scalar>& B_avg,
115 DeferredLogger& deferred_logger) const;
116
117 double getControlTolerance(const WellState<Scalar>& well_state,
118 const double tolerance_wells,
119 const double tolerance_pressure_ms_wells,
120 DeferredLogger& deferred_logger) const;
121
123 const std::vector<double>& residuals,
124 const double tolerance_wells,
125 const double tolerance_pressure_ms_wells,
126 DeferredLogger& deferred_logger) const;
127
129 WellState<Scalar>& well_state);
130
132 const UnitSystem& unit_system) const;
133
134 // convert a Eval from reservoir to contain the derivative related to wells
135 EvalWell extendEval(const Eval& in) const;
136
138
142
143 // depth difference between perforations and the perforated grid cells
144 std::vector<double> cell_perforation_depth_diffs_;
145 // pressure correction due to the different depth of the perforation and
146 // center depth of the grid block
148};
149
150}
151
152#endif // OPM_MULTISEGMENTWELL_GENERIC_HEADER_INCLUDED
Definition: ConvergenceReport.hpp:38
Definition: DeferredLogger.hpp:57
Dune::BlockVector< VectorBlockType > BVector
Definition: MultisegmentWellEquations.hpp:60
Dune::BlockVector< VectorBlockWellType > BVectorWell
Definition: MultisegmentWellEquations.hpp:57
Definition: MultisegmentWellEval.hpp:47
ConvergenceReport getWellConvergence(const WellState< Scalar > &well_state, const std::vector< double > &B_avg, DeferredLogger &deferred_logger, const double max_residual_allowed, const double tolerance_wells, const double relaxed_inner_tolerance_flow_ms_well, const double tolerance_pressure_ms_wells, const double relaxed_inner_tolerance_pressure_ms_well, const bool relax_tolerance, const bool well_is_stopped) const
check whether the well equations get converged for this well
std::vector< double > cell_perforation_depth_diffs_
Definition: MultisegmentWellEval.hpp:144
double getResidualMeasureValue(const WellState< Scalar > &well_state, const std::vector< double > &residuals, const double tolerance_wells, const double tolerance_pressure_ms_wells, DeferredLogger &deferred_logger) const
static constexpr int numWellEq
Definition: MultisegmentWellEval.hpp:51
double getControlTolerance(const WellState< Scalar > &well_state, const double tolerance_wells, const double tolerance_pressure_ms_wells, DeferredLogger &deferred_logger) const
EvalWell extendEval(const Eval &in) const
static constexpr int SPres
Definition: MultisegmentWellEval.hpp:52
static constexpr int WQTotal
Definition: MultisegmentWellEval.hpp:53
const WellInterfaceIndices< FluidSystem, Indices > & baseif_
Definition: MultisegmentWellEval.hpp:137
void assembleICDPressureEq(const int seg, const UnitSystem &unit_system, WellState< Scalar > &well_state, const SummaryState &summary_state, const bool use_average_density, DeferredLogger &deferred_logger)
std::pair< bool, std::vector< Scalar > > getFiniteWellResiduals(const std::vector< Scalar > &B_avg, DeferredLogger &deferred_logger) const
Equations linSys_
The equation system.
Definition: MultisegmentWellEval.hpp:139
typename PrimaryVariables::EvalWell EvalWell
Definition: MultisegmentWellEval.hpp:64
void assembleAccelerationAndHydroPressureLosses(const int seg, WellState< Scalar > &well_state, const bool use_average_density)
typename Equations::BVectorWell BVectorWell
Definition: MultisegmentWellEval.hpp:59
const Equations & linSys() const
Returns a const reference to equation system.
Definition: MultisegmentWellEval.hpp:69
void assemblePressureEq(const int seg, const UnitSystem &unit_system, WellState< Scalar > &well_state, const SummaryState &summary_state, const bool use_average_density, DeferredLogger &deferred_logger)
PrimaryVariables primary_variables_
The primary variables.
Definition: MultisegmentWellEval.hpp:140
MultisegmentWellEval(WellInterfaceIndices< FluidSystem, Indices > &baseif)
void assembleDefaultPressureEq(const int seg, WellState< Scalar > &well_state, const bool use_average_density)
MSWSegments segments_
Segment properties.
Definition: MultisegmentWellEval.hpp:141
EvalWell pressureDropAutoICD(const int seg, const UnitSystem &unit_system) const
DenseAd::Evaluation< Scalar, Indices::numEq > Eval
Definition: MultisegmentWellEval.hpp:65
std::vector< double > cell_perforation_pressure_diffs_
Definition: MultisegmentWellEval.hpp:147
void initMatrixAndVectors(const int num_cells)
typename FluidSystem::Scalar Scalar
Definition: MultisegmentWellEval.hpp:49
void assembleAccelerationPressureLoss(const int seg, WellState< Scalar > &well_state)
typename Equations::BVector BVector
Definition: MultisegmentWellEval.hpp:58
Definition: MultisegmentWellGeneric.hpp:42
Definition: MultisegmentWellPrimaryVariables.hpp:45
DenseAd::Evaluation< double, Indices::numEq+numWellEq > EvalWell
Definition: MultisegmentWellPrimaryVariables.hpp:79
static constexpr int numWellEq
Definition: MultisegmentWellPrimaryVariables.hpp:76
static constexpr int SPres
Definition: MultisegmentWellPrimaryVariables.hpp:73
static constexpr int WQTotal
Definition: MultisegmentWellPrimaryVariables.hpp:70
Definition: MultisegmentWellSegments.hpp:45
Definition: WellInterfaceIndices.hpp:35
Definition: WellState.hpp:62
Definition: BlackoilPhases.hpp:27