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 double dt,
40 const Domain& domain)
41{
42 OPM_TIMEBLOCK(assemble);
43 // NLDD domain well assembly for local solves.
44 //
45 // Call chain: BlackoilWellModelNldd::assemble()
46 // -> BlackoilWellModelNldd::assembleWellEq(dt, domain)
47 // -> WellInterface::assembleWellEq()
48 // -> WellInterface::prepareWellBeforeAssembling()
49 // -> WellInterface::assembleWellEqWithoutIteration()
50 //
51 // Key assumptions and design decisions:
52 //
53 // 1. Global well initialization (BlackoilWellModel::prepareTimeStep, including the
54 // optional initial well solve) has already been called via
55 // BlackoilWellModel::assemble() during the initial global linearization.
56 // This is ensured by NewtonIterationContext::needsTimestepInit().
57 //
58 // 2. Inner well iterations (BlackoilWellModel::iterateWellEquations) are skipped
59 // during local solves via NewtonIterationContext::shouldRunInnerWellIterations()
60 // (inLocalSolve=true).
61 // The NLDD local Newton loop already iterates wells to convergence.
62 //
63 // 3. No MPI collective operations (do_mpi_gather=false).
64 //
65 // 4. Only individual well constraints are checked (see
66 // BlackoilWellModelNldd::updateWellControls).
67 // Group controls, network balancing, gas lift, and NUPCOL/VREP updates
68 // are NOT performed during domain solves. This means:
69 // - Wells can switch between individual controls (BHP, THP, rate limits)
70 // to respect physical safety constraints.
71 // - Wells cannot be switched to/from GRUP control. Individual constraint
72 // checks (WellInterface::checkIndividualConstraints) never produce GRUP
73 // as a result; only WellInterface::checkGroupConstraints can do that,
74 // and it is not called here.
75 // - Group production targets and VREP/REIN injection targets use stale
76 // group state from the last global assembly.
77 // - Well control switches during domain solves are not logged to
78 // well_control_log_ (gated by inLocalSolve), preventing local
79 // oscillations from exhausting the global oscillation budget.
80 //
81 // The final Newton step after NLDD domain solves re-synchronizes everything
82 // via a global BlackoilWellModel::assemble() which performs the full control
83 // update pipeline: BlackoilWellModelGeneric::updateAndCommunicateGroupData,
84 // BlackoilWellModelGeneric::updateGroupControls,
85 // group+individual well constraint checks, network balancing, and gas lift.
86 // Importantly, the global step also checks individual constraints (after
87 // group constraints), so if a domain solve switched a well to e.g. BHP
88 // because of a genuine BHP limit violation, the global step will reach
89 // the same conclusion -- the two levels are consistent for individual
90 // constraints.
91
92 // Use do_mpi_gather=false to avoid MPI collective operations in domain solves.
93 auto loggerGuard = wellModel_.groupStateHelper().pushLogger(/*do_mpi_gather=*/false);
94
95 this->updateWellControls(domain);
96 this->assembleWellEq(dt, domain);
97
98 // Update cellRates_ with current contributions from wells in this domain for reservoir linearization
99 wellModel_.updateCellRatesForDomain(domain.index, this->well_domain());
100}
101
102template<typename TypeTag>
103void
105assembleWellEq(const double dt,
106 const Domain& domain)
107{
108 OPM_TIMEBLOCK(assembleWellEq);
109 for (const auto& well : wellModel_.localNonshutWells()) {
110 if (this->well_domain().at(well->name()) == domain.index) {
111 well->assembleWellEq(wellModel_.simulator(),
112 dt,
113 wellModel_.groupStateHelper(),
114 wellModel_.wellState());
115 }
116 }
117}
118
119template<typename TypeTag>
120void
123 const BVector& /*weights*/,
124 const bool /*use_well_weights*/,
125 const int /*domainIndex*/) const
126{
127 throw std::logic_error("CPRW is not yet implemented for NLDD subdomains");
128 // To fix this function, rdofs should be the size of the domain, and the nw should be the number of wells in the domain
129 // int nw = this->numLocalWellsEnd(); // should number of wells in the domain
130 // int rdofs = local_num_cells_; // should be the size of the domain
131 // for ( int i = 0; i < nw; i++ ) {
132 // int wdof = rdofs + i;
133 // jacobian[wdof][wdof] = 1.0;// better scaling ?
134 // }
135
136 // for ( const auto& well : well_container_ ) {
137 // if (well_domain_.at(well->name()) == domainIndex) {
138 // weights should be the size of the domain
139 // well->addWellPressureEquations(jacobian, weights, pressureVarIndex, use_well_weights, this->wellState());
140 // }
141 // }
142}
143
144template<typename TypeTag>
145void
148 const int domainIdx)
149{
150 // Note: no point in trying to do a parallel gathering
151 // try/catch here, as this function is not called in
152 // parallel but for each individual domain of each rank.
153 // Use do_mpi_gather=false to avoid MPI collective operations.
154 auto loggerGuard = wellModel_.groupStateHelper().pushLogger(/*do_mpi_gather=*/false);
155 for (const auto& well : wellModel_.localNonshutWells()) {
156 if (this->well_domain().at(well->name()) == domainIdx) {
157 const auto& cells = well->cells();
158 x_local_.resize(cells.size());
159
160 for (size_t i = 0; i < cells.size(); ++i) {
161 x_local_[i] = x[cells[i]];
162 }
163 well->recoverWellSolutionAndUpdateWellState(wellModel_.simulator(),
164 x_local_,
165 wellModel_.groupStateHelper(),
166 wellModel_.wellState());
167 }
168 }
169}
170
171template<typename TypeTag>
174getWellConvergence(const Domain& domain,
175 const std::vector<Scalar>& B_avg,
176 DeferredLogger& local_deferredLogger) const
177{
178 const auto& iterCtx = wellModel_.simulator().problem().iterationContext();
179 const bool relax_tolerance = iterCtx.shouldRelax(wellModel_.numStrictIterations() + 1);
180
181 ConvergenceReport report;
182 {
183 // Use do_mpi_gather=false to avoid MPI collective operations in domain solves.
184 auto loggerGuard = wellModel_.groupStateHelper().pushLogger(/*do_mpi_gather=*/false);
185
186 for (const auto& well : wellModel_.localNonshutWells()) {
187 if ((this->well_domain().at(well->name()) == domain.index)) {
188 if (well->isOperableAndSolvable() || well->wellIsStopped()) {
189 report += well->getWellConvergence(wellModel_.groupStateHelper(),
190 B_avg,
191 relax_tolerance);
192 } else {
193 ConvergenceReport xreport;
194 using CR = ConvergenceReport;
195 xreport.setWellFailed({CR::WellFailure::Type::Unsolvable,
196 CR::Severity::Normal, -1, well->name()});
197 report += xreport;
198 }
199 }
200 }
201 } // loggerGuard goes out of scope here, before the OpmLog::debug() calls below
202
203 // Log debug messages for NaN or too large residuals.
204 if (wellModel_.terminalOutput()) {
205 for (const auto& f : report.wellFailures()) {
206 if (f.severity() == ConvergenceReport::Severity::NotANumber) {
207 local_deferredLogger.debug("NaN residual found with phase " +
208 std::to_string(f.phase()) +
209 " for well " + f.wellName());
210 } else if (f.severity() == ConvergenceReport::Severity::TooLarge) {
211 local_deferredLogger.debug("Too large residual found with phase " +
212 std::to_string(f.phase()) +
213 " for well " + f.wellName());
214 }
215 }
216 }
217 return report;
218}
219
220template<typename TypeTag>
221void
223updateWellControls(const Domain& domain)
224{
225 OPM_TIMEBLOCK(updateWellControls);
226 if (!wellModel_.wellsActive()) {
227 return;
228 }
229
230 // Only individual well constraints (BHP, THP, rate limits) are checked.
231 // Group constraints, network balancing, and NUPCOL/VREP updates require
232 // MPI collectives and consistent global group state, which are not
233 // available during domain-local solves. The final global Newton step
234 // handles all group-level control decisions.
235 //
236 // Individual constraints are physical safety limits that are essential
237 // for well convergence (e.g., BHP minimum prevents infeasible operating
238 // points). Skipping them would cause domain convergence failures.
239 for (const auto& well : wellModel_.localNonshutWells()) {
240 if (this->well_domain().at(well->name()) == domain.index) {
242 well->updateWellControl(wellModel_.simulator(),
243 mode,
244 wellModel_.groupStateHelper(),
245 wellModel_.wellState());
246 }
247 }
248}
249
250template <typename TypeTag>
251void
253setupDomains(const std::vector<Domain>& domains)
254{
255 std::vector<const SubDomainIndices*> genDomains;
256 std::ranges::transform(domains, std::back_inserter(genDomains),
257 [](const auto& domain)
258 { return static_cast<const SubDomainIndices*>(&domain); });
259 this->calcDomains(genDomains);
260}
261
262} // namespace Opm
263
264#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:174
void recoverWellSolutionAndUpdateWellState(const BVector &x, const int domainIdx)
Definition: BlackoilWellModelNldd_impl.hpp:147
void assemble(const double dt, const Domain &domain)
Definition: BlackoilWellModelNldd_impl.hpp:39
typename BlackoilWellModel< TypeTag >::BVector BVector
Definition: BlackoilWellModelNldd.hpp:86
void setupDomains(const std::vector< Domain > &domains)
Definition: BlackoilWellModelNldd_impl.hpp:253
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:122
void updateWellControls(const Domain &domain)
Definition: BlackoilWellModelNldd_impl.hpp:223
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 GroupStateHelperType &groupStateHelper, WellStateType &well_state)
Definition: WellInterface_impl.hpp:190
Definition: blackoilbioeffectsmodules.hh:45
std::string to_string(const ConvergenceReport::ReservoirFailure::Type t)
Definition: SubDomain.hpp:61
int index
Definition: SubDomain.hpp:64
Definition: SubDomain.hpp:85