opm-simulators
Loading...
Searching...
No Matches
MultisegmentWellEval.hpp
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
25#include <opm/simulators/wells/MultisegmentWellEquations.hpp>
26#include <opm/simulators/wells/MultisegmentWellGeneric.hpp>
27#include <opm/simulators/wells/MultisegmentWellPrimaryVariables.hpp>
28#include <opm/simulators/wells/MultisegmentWellSegments.hpp>
29#include <opm/simulators/wells/ParallelWellInfo.hpp>
30
31#include <opm/material/densead/Evaluation.hpp>
32
33#include <utility>
34#include <vector>
35
36namespace Opm {
37
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;
61 using BVectorWell = typename Equations::BVectorWell;
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.
66 using EvalWell = typename PrimaryVariables::EvalWell;
67 using Eval = DenseAd::Evaluation<Scalar, /*size=*/Indices::numEq>;
68
69public:
71 const Equations& linSys() const
72 { return linSys_; }
73 const ParallelWellInfo<Scalar>& pw_info_;
74
75protected:
77
78 void initMatrixAndVectors();
79
80 void assembleDefaultPressureEq(const int seg,
82 const bool use_average_density);
83
84 // assemble pressure equation for ICD segments
85 void assembleICDPressureEq(const int seg,
86 const UnitSystem& unit_system,
88 const SummaryState& summary_state,
89 const bool use_average_density,
90 DeferredLogger& deferred_logger);
91
92 void assembleAccelerationAndHydroPressureLosses(const int seg,
94 const bool use_average_density);
95
96
97 void assemblePressureEq(const int seg,
98 const UnitSystem& unit_system,
100 const SummaryState& summary_state,
101 const bool use_average_density,
102 DeferredLogger& deferred_logger);
103
106 const std::vector<Scalar>& B_avg,
107 DeferredLogger& deferred_logger,
108 const Scalar max_residual_allowed,
109 const Scalar tolerance_wells,
110 const Scalar relaxed_inner_tolerance_flow_ms_well,
111 const Scalar tolerance_pressure_ms_wells,
112 const Scalar relaxed_inner_tolerance_pressure_ms_well,
113 const bool relax_tolerance,
114 const bool well_is_stopped) const;
115
116 std::pair<bool, std::vector<Scalar> >
117 getFiniteWellResiduals(const std::vector<Scalar>& B_avg,
118 DeferredLogger& deferred_logger) const;
119
120 Scalar getControlTolerance(const WellState<Scalar, IndexTraits>& well_state,
121 const Scalar tolerance_wells,
122 const Scalar tolerance_pressure_ms_wells,
123 DeferredLogger& deferred_logger) const;
124
125 Scalar getResidualMeasureValue(const WellState<Scalar, IndexTraits>& well_state,
126 const std::vector<Scalar>& residuals,
127 const Scalar tolerance_wells,
128 const Scalar tolerance_pressure_ms_wells,
129 DeferredLogger& deferred_logger) const;
130
131 void assembleAccelerationPressureLoss(const int seg,
133
134 EvalWell pressureDropAutoICD(const int seg,
135 const UnitSystem& unit_system) const;
136
137 // convert a Eval from reservoir to contain the derivative related to wells
138 EvalWell extendEval(const Eval& in) const;
139
141
142 Equations linSys_;
143 PrimaryVariables primary_variables_;
144 MSWSegments segments_;
145
146 // depth difference between perforations and the perforated grid cells
147 std::vector<Scalar> cell_perforation_depth_diffs_;
148 // pressure correction due to the different depth of the perforation and
149 // center depth of the grid block
150 std::vector<Scalar> cell_perforation_pressure_diffs_;
151};
152
153}
154
155#endif // OPM_MULTISEGMENTWELL_GENERIC_HEADER_INCLUDED
Represents the convergence status of the whole simulator, to make it possible to query and store the ...
Definition ConvergenceReport.hpp:38
Definition DeferredLogger.hpp:57
Definition MultisegmentWellEquations.hpp:52
Definition MultisegmentWellEval.hpp:48
const Equations & linSys() const
Returns a const reference to equation system.
Definition MultisegmentWellEval.hpp:71
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
Definition MultisegmentWellEval.cpp:82
Definition MultisegmentWellPrimaryVariables.hpp:45
Definition MultisegmentWellSegments.hpp:45
Class encapsulating some information about parallel wells.
Definition ParallelWellInfo.hpp:198
Definition WellInterfaceIndices.hpp:34
The state of a set of wells, tailored for use by the fully implicit blackoil simulator.
Definition WellState.hpp:66
This file contains a set of helper functions used by VFPProd / VFPInj.
Definition blackoilbioeffectsmodules.hh:43