BlackoilWellModelNldd.hpp
Go to the documentation of this file.
1/*
2 Copyright 2016 SINTEF ICT, Applied Mathematics.
3 Copyright 2016 - 2017 Statoil ASA.
4 Copyright 2017 Dr. Blatt - HPC-Simulation-Software & Services
5 Copyright 2016 - 2018 IRIS AS
6
7 This file is part of the Open Porous Media project (OPM).
8
9 OPM is free software: you can redistribute it and/or modify
10 it under the terms of the GNU General Public License as published by
11 the Free Software Foundation, either version 3 of the License, or
12 (at your option) any later version.
13
14 OPM is distributed in the hope that it will be useful,
15 but WITHOUT ANY WARRANTY; without even the implied warranty of
16 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
17 GNU General Public License for more details.
18
19 You should have received a copy of the GNU General Public License
20 along with OPM. If not, see <http://www.gnu.org/licenses/>.
21*/
22#ifndef OPM_BLACKOILWELLMODEL_NLDD_HEADER_INCLUDED
23#define OPM_BLACKOILWELLMODEL_NLDD_HEADER_INCLUDED
24
25#include <opm/grid/utility/SparseTable.hpp>
26
28
30
32
33#include <map>
34#include <string>
35#include <vector>
36
37namespace Opm {
38
39template<typename Scalar, typename IndexTraits>
41{
42public:
43 std::vector<Scalar> getPrimaryVarsDomain(const int domainIdx) const;
44 void setPrimaryVarsDomain(const int domainIdx, const std::vector<Scalar>& vars);
45
47 { return well_local_cells_; }
48
49 const std::map<std::string, int>& well_domain() const
50 { return well_domain_; }
51
52protected:
54 : genWellModel_(model)
55 {}
56
57 void calcDomains(const std::vector<const SubDomainIndices*>& domains);
58
59private:
60 void logDomains() const;
61
62 void findWellDomains(const std::vector<const SubDomainIndices*>& domains);
63
64 void calcLocalIndices(const std::vector<const SubDomainIndices*>& domains);
65
67
68 // Keep track of the domain of each well
69 std::map<std::string, int> well_domain_{};
70
71 // Store the local index of the wells perforated cells in the domain
72 SparseTable<int> well_local_cells_{};
73};
74
76template<typename TypeTag>
78 public BlackoilWellModelNlddGeneric<GetPropType<TypeTag, Properties::Scalar>,
79 typename GetPropType<TypeTag, Properties::FluidSystem>::IndexTraitsType>
80{
81public:
82 // --------- Types ---------
88 using IndexTraits = typename FluidSystem::IndexTraitsType;
90
92
95 , wellModel_(model)
96 {}
97
99 const BVector& weights,
100 const bool use_well_weights,
101 const int domainIndex) const;
102
103 // prototype for assemble function for ASPIN solveLocal()
104 // will try to merge back to assemble() when done prototyping
105 void assemble(const int iterationIdx,
106 const double dt,
107 const Domain& domain);
108
109 void updateWellControls(DeferredLogger& deferred_logger,
110 const Domain& domain);
111
112 void setupDomains(const std::vector<Domain>& domains);
113
114 // Check if well equations are converged locally.
116 const std::vector<Scalar>& B_avg,
117 DeferredLogger& local_deferredLogger) const;
118
119 // using the solution x to recover the solution xw for wells and applying
120 // xw to update Well State
122 const int domainIdx);
123
124 // Get number of wells on this rank
125 int numLocalWells() const
126 {
127 return wellModel_.numLocalWells();
128 }
129
130 // Get number of wells on this rank
131 int numLocalWellsEnd() const
132 {
133 return wellModel_.numLocalWellsEnd();
134 }
135
136private:
137 BlackoilWellModel<TypeTag>& wellModel_;
138
139 void assembleWellEq(const double dt,
140 const Domain& domain,
141 DeferredLogger& deferred_logger);
142
143 // These members are used to avoid reallocation in specific functions
144 // instead of using local variables.
145 // Their state is not relevant between function calls, so they can
146 // (and must) be mutable, as the functions using them are const.
147 mutable BVector x_local_;
148
149};
150
151} // namespace Opm
152
154
155#endif
Class for handling the blackoil well model.
Definition: BlackoilWellModelGeneric.hpp:95
Definition: BlackoilWellModelNldd.hpp:41
const SparseTable< int > & well_local_cells() const
Definition: BlackoilWellModelNldd.hpp:46
BlackoilWellModelNlddGeneric(BlackoilWellModelGeneric< Scalar, IndexTraits > &model)
Definition: BlackoilWellModelNldd.hpp:53
std::vector< Scalar > getPrimaryVarsDomain(const int domainIdx) const
void calcDomains(const std::vector< const SubDomainIndices * > &domains)
const std::map< std::string, int > & well_domain() const
Definition: BlackoilWellModelNldd.hpp:49
void setPrimaryVarsDomain(const int domainIdx, const std::vector< Scalar > &vars)
Class for handling the blackoil well model in a NLDD solver.
Definition: BlackoilWellModelNldd.hpp:80
ConvergenceReport getWellConvergence(const Domain &domain, const std::vector< Scalar > &B_avg, DeferredLogger &local_deferredLogger) const
Definition: BlackoilWellModelNldd_impl.hpp:132
GetPropType< TypeTag, Properties::Grid > Grid
Definition: BlackoilWellModelNldd.hpp:83
void assemble(const int iterationIdx, const double dt, const Domain &domain)
Definition: BlackoilWellModelNldd_impl.hpp:39
int numLocalWellsEnd() const
Definition: BlackoilWellModelNldd.hpp:131
GetPropType< TypeTag, Properties::Indices > Indices
Definition: BlackoilWellModelNldd.hpp:89
SubDomain< Grid > Domain
Definition: BlackoilWellModelNldd.hpp:91
typename FluidSystem::IndexTraitsType IndexTraits
Definition: BlackoilWellModelNldd.hpp:88
void recoverWellSolutionAndUpdateWellState(const BVector &x, const int domainIdx)
Definition: BlackoilWellModelNldd_impl.hpp:101
GetPropType< TypeTag, Properties::Scalar > Scalar
Definition: BlackoilWellModelNldd.hpp:84
BlackoilWellModelNldd(BlackoilWellModel< TypeTag > &model)
Definition: BlackoilWellModelNldd.hpp:93
typename BlackoilWellModel< TypeTag >::BVector BVector
Definition: BlackoilWellModelNldd.hpp:86
void setupDomains(const std::vector< Domain > &domains)
Definition: BlackoilWellModelNldd_impl.hpp:205
typename BlackoilWellModel< TypeTag >::PressureMatrix PressureMatrix
Definition: BlackoilWellModelNldd.hpp:85
void addWellPressureEquations(PressureMatrix &jacobian, const BVector &weights, const bool use_well_weights, const int domainIndex) const
Definition: BlackoilWellModelNldd_impl.hpp:76
GetPropType< TypeTag, Properties::FluidSystem > FluidSystem
Definition: BlackoilWellModelNldd.hpp:87
int numLocalWells() const
Definition: BlackoilWellModelNldd.hpp:125
void updateWellControls(DeferredLogger &deferred_logger, const Domain &domain)
Definition: BlackoilWellModelNldd_impl.hpp:178
Class for handling the blackoil well model.
Definition: BlackoilWellModel.hpp:102
Dune::BCRSMatrix< Opm::MatrixBlock< Scalar, 1, 1 > > PressureMatrix
Definition: BlackoilWellModel.hpp:299
Dune::BlockVector< VectorBlockType > BVector
Definition: BlackoilWellModel.hpp:133
Definition: ConvergenceReport.hpp:38
Definition: DeferredLogger.hpp:57
Definition: blackoilboundaryratevector.hh:39
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property (equivalent to old macro GET_PROP_TYPE(....
Definition: propertysystem.hh:233
Definition: SubDomain.hpp:85