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