BlackoilWellModelNldd_impl.hpp
Go to the documentation of this file.
1/*
2 Copyright 2016 - 2019 SINTEF Digital, Mathematics & Cybernetics.
3 Copyright 2016 - 2018 Equinor ASA.
4 Copyright 2017 Dr. Blatt - HPC-Simulation-Software & Services
5 Copyright 2016 - 2018 Norce 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
23#ifndef OPM_BLACKOILWELLMODEL_NLDD_IMPL_HEADER_INCLUDED
24#define OPM_BLACKOILWELLMODEL_NLDD_IMPL_HEADER_INCLUDED
25
26// Improve IDE experience
27#ifndef OPM_BLACKOILWELLMODEL_NLDD_HEADER_INCLUDED
28#include <config.h>
30#endif
31
32#include <algorithm>
33
34namespace Opm {
35
36template<typename TypeTag>
37void
39assemble(const int /*iterationIdx*/,
40 const double dt,
41 const Domain& domain)
42{
43 OPM_TIMEBLOCK(assemble);
44 // We assume that calculateExplicitQuantities() and
45 // prepareTimeStep() have been called already for the entire
46 // well model, so we do not need to do it here (when
47 // iterationIdx is 0).
48
49 DeferredLogger local_deferredLogger;
50 this->updateWellControls(local_deferredLogger, domain);
51 this->assembleWellEq(dt, domain, local_deferredLogger);
52
53 // Update cellRates_ with current contributions from wells in this domain for reservoir linearization
54 wellModel_.updateCellRatesForDomain(domain.index, this->well_domain());
55}
56
57template<typename TypeTag>
58void
60assembleWellEq(const double dt,
61 const Domain& domain,
62 DeferredLogger& deferred_logger)
63{
64 OPM_TIMEBLOCK(assembleWellEq);
65 for (const auto& well : wellModel_.localNonshutWells()) {
66 if (this->well_domain().at(well->name()) == domain.index) {
67 well->assembleWellEq(wellModel_.simulator(),
68 dt,
69 wellModel_.wgHelper(),
70 wellModel_.wellState(),
71 deferred_logger);
72 }
73 }
74}
75
76template<typename TypeTag>
77void
80 const BVector& /*weights*/,
81 const bool /*use_well_weights*/,
82 const int /*domainIndex*/) const
83{
84 throw std::logic_error("CPRW is not yet implemented for NLDD subdomains");
85 // To fix this function, rdofs should be the size of the domain, and the nw should be the number of wells in the domain
86 // int nw = this->numLocalWellsEnd(); // should number of wells in the domain
87 // int rdofs = local_num_cells_; // should be the size of the domain
88 // for ( int i = 0; i < nw; i++ ) {
89 // int wdof = rdofs + i;
90 // jacobian[wdof][wdof] = 1.0;// better scaling ?
91 // }
92
93 // for ( const auto& well : well_container_ ) {
94 // if (well_domain_.at(well->name()) == domainIndex) {
95 // weights should be the size of the domain
96 // well->addWellPressureEquations(jacobian, weights, pressureVarIndex, use_well_weights, this->wellState());
97 // }
98 // }
99}
100
101template<typename TypeTag>
102void
105 const int domainIdx)
106{
107 // Note: no point in trying to do a parallel gathering
108 // try/catch here, as this function is not called in
109 // parallel but for each individual domain of each rank.
110 DeferredLogger local_deferredLogger;
111 for (const auto& well : wellModel_.localNonshutWells()) {
112 if (this->well_domain().at(well->name()) == domainIdx) {
113 const auto& cells = well->cells();
114 x_local_.resize(cells.size());
115
116 for (size_t i = 0; i < cells.size(); ++i) {
117 x_local_[i] = x[cells[i]];
118 }
119 well->recoverWellSolutionAndUpdateWellState(wellModel_.simulator(),
120 x_local_,
121 wellModel_.wellState(),
122 local_deferredLogger);
123 }
124 }
125 // TODO: avoid losing the logging information that could
126 // be stored in the local_deferredlogger in a parallel case.
127 if (wellModel_.terminalOutput()) {
128 local_deferredLogger.logMessages();
129 }
130}
131
132template<typename TypeTag>
135getWellConvergence(const Domain& domain,
136 const std::vector<Scalar>& B_avg,
137 DeferredLogger& local_deferredLogger) const
138{
139 const int iterationIdx = wellModel_.simulator().model().newtonMethod().numIterations();
140 const bool relax_tolerance = iterationIdx > wellModel_.numStrictIterations();
141
142 ConvergenceReport report;
143 for (const auto& well : wellModel_.localNonshutWells()) {
144 if ((this->well_domain().at(well->name()) == domain.index)) {
145 if (well->isOperableAndSolvable() || well->wellIsStopped()) {
146 report += well->getWellConvergence(wellModel_.simulator(),
147 wellModel_.wellState(),
148 B_avg,
149 local_deferredLogger,
150 relax_tolerance);
151 } else {
152 ConvergenceReport xreport;
153 using CR = ConvergenceReport;
154 xreport.setWellFailed({CR::WellFailure::Type::Unsolvable,
155 CR::Severity::Normal, -1, well->name()});
156 report += xreport;
157 }
158 }
159 }
160
161 // Log debug messages for NaN or too large residuals.
162 if (wellModel_.terminalOutput()) {
163 for (const auto& f : report.wellFailures()) {
164 if (f.severity() == ConvergenceReport::Severity::NotANumber) {
165 local_deferredLogger.debug("NaN residual found with phase " +
166 std::to_string(f.phase()) +
167 " for well " + f.wellName());
168 } else if (f.severity() == ConvergenceReport::Severity::TooLarge) {
169 local_deferredLogger.debug("Too large residual found with phase " +
170 std::to_string(f.phase()) +
171 " for well " + f.wellName());
172 }
173 }
174 }
175 return report;
176}
177
178template<typename TypeTag>
179void
181updateWellControls(DeferredLogger& deferred_logger,
182 const Domain& domain)
183{
184 OPM_TIMEBLOCK(updateWellControls);
185 if (!wellModel_.wellsActive()) {
186 return;
187 }
188
189 // TODO: decide on and implement an approach to handling of
190 // group controls, network and similar for domain solves.
191
192 // Check only individual well constraints and communicate.
193 for (const auto& well : wellModel_.localNonshutWells()) {
194 if (this->well_domain().at(well->name()) == domain.index) {
196 well->updateWellControl(wellModel_.simulator(),
197 mode,
198 wellModel_.wgHelper(),
199 wellModel_.wellState(),
200 deferred_logger);
201 }
202 }
203}
204
205template <typename TypeTag>
206void
208setupDomains(const std::vector<Domain>& domains)
209{
210 std::vector<const SubDomainIndices*> genDomains;
211 std::transform(domains.begin(), domains.end(),
212 std::back_inserter(genDomains),
213 [](const auto& domain)
214 { return static_cast<const SubDomainIndices*>(&domain); });
215 this->calcDomains(genDomains);
216}
217
218} // namespace Opm
219
220#endif
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:135
void assemble(const int iterationIdx, const double dt, const Domain &domain)
Definition: BlackoilWellModelNldd_impl.hpp:39
void recoverWellSolutionAndUpdateWellState(const BVector &x, const int domainIdx)
Definition: BlackoilWellModelNldd_impl.hpp:104
typename BlackoilWellModel< TypeTag >::BVector BVector
Definition: BlackoilWellModelNldd.hpp:86
void setupDomains(const std::vector< Domain > &domains)
Definition: BlackoilWellModelNldd_impl.hpp:208
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:79
void updateWellControls(DeferredLogger &deferred_logger, const Domain &domain)
Definition: BlackoilWellModelNldd_impl.hpp:181
Definition: ConvergenceReport.hpp:38
void setWellFailed(const WellFailure &wf)
Definition: ConvergenceReport.hpp:272
const std::vector< WellFailure > & wellFailures() const
Definition: ConvergenceReport.hpp:380
Definition: DeferredLogger.hpp:57
void debug(const std::string &tag, const std::string &message)
Definition: WellInterface.hpp:77
bool updateWellControl(const Simulator &simulator, const IndividualOrGroup iog, const WellGroupHelperType &wgHelper, WellStateType &well_state, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:189
Definition: blackoilbioeffectsmodules.hh:43
std::string to_string(const ConvergenceReport::ReservoirFailure::Type t)
int index
Definition: SubDomain.hpp:64
Definition: SubDomain.hpp:85