BlackoilModel_impl.hpp
Go to the documentation of this file.
1/*
2 Copyright 2013, 2015 SINTEF ICT, Applied Mathematics.
3 Copyright 2014, 2015 Dr. Blatt - HPC-Simulation-Software & Services
4 Copyright 2014, 2015 Statoil ASA.
5 Copyright 2015 NTNU
6 Copyright 2015, 2016, 2017 IRIS AS
7
8 This file is part of the Open Porous Media project (OPM).
9
10 OPM is free software: you can redistribute it and/or modify
11 it under the terms of the GNU General Public License as published by
12 the Free Software Foundation, either version 3 of the License, or
13 (at your option) any later version.
14
15 OPM is distributed in the hope that it will be useful,
16 but WITHOUT ANY WARRANTY; without even the implied warranty of
17 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
18 GNU General Public License for more details.
19
20 You should have received a copy of the GNU General Public License
21 along with OPM. If not, see <http://www.gnu.org/licenses/>.
22*/
23
24#ifndef OPM_BLACKOILMODEL_IMPL_HEADER_INCLUDED
25#define OPM_BLACKOILMODEL_IMPL_HEADER_INCLUDED
26
27#ifndef OPM_BLACKOILMODEL_HEADER_INCLUDED
28#include <config.h>
30#endif
31
32#include <dune/common/timer.hh>
33
34#include <opm/common/ErrorMacros.hpp>
35#include <opm/common/OpmLog/OpmLog.hpp>
36
38
39#include <algorithm>
40#include <cmath>
41#include <filesystem>
42#include <functional>
43#include <iomanip>
44#include <limits>
45#include <memory>
46#include <numeric>
47#include <sstream>
48#include <stdexcept>
49#include <string>
50#include <string_view>
51#include <tuple>
52#include <utility>
53#include <vector>
54
55#include <fmt/format.h>
56
57namespace {
58 template <typename TypeTag>
59 std::string_view
60 make_string(const typename Opm::BlackoilModel<TypeTag>::DebugFlags f)
61 {
63
64 switch (f) {
65 case F::STRICT: return "Strict";
66 case F::RELAXED: return "Relaxed";
67 case F::TUNINGDP: return "TuningDP";
68 }
69
70 return "< ??? >";
71 }
72} // Anonymous namespace
73
74namespace Opm {
75
76template <class TypeTag>
78BlackoilModel(Simulator& simulator,
79 const ModelParameters& param,
81 const bool terminal_output)
82 : simulator_(simulator)
83 , grid_(simulator_.vanguard().grid())
84 , param_( param )
85 , well_model_ (well_model)
86 , terminal_output_ (terminal_output)
87 , current_relaxation_(1.0)
88 , dx_old_(simulator_.model().numGridDof())
89 , conv_monitor_(param_.monitor_params_)
90{
91 // compute global sum of number of cells
93 convergence_reports_.reserve(300); // Often insufficient, but avoids frequent moves.
94 // TODO: remember to fix!
95 if (param_.nonlinear_solver_ == "nldd") {
96 if (terminal_output) {
97 OpmLog::info("Using Non-Linear Domain Decomposition solver (nldd).");
98 }
99 nlddSolver_ = std::make_unique<BlackoilModelNldd<TypeTag>>(*this);
100 } else if (param_.nonlinear_solver_ == "newton") {
101 if (terminal_output) {
102 OpmLog::info("Using Newton nonlinear solver.");
103 }
104 } else {
105 OPM_THROW(std::runtime_error, "Unknown nonlinear solver option: " +
107 }
108}
109
110template <class TypeTag>
114{
115 OPM_TIMEFUNCTION();
117 Dune::Timer perfTimer;
118 perfTimer.start();
119 // update the solution variables in the model
120 int lastStepFailed = timer.lastStepFailed();
121 if (grid_.comm().size() > 1 && grid_.comm().max(lastStepFailed) != grid_.comm().min(lastStepFailed)) {
122 OPM_THROW(std::runtime_error, "Misalignment of the parallel simulation run in prepareStep " +
123 "- the previous step succeeded on some ranks but failed on others.");
124 }
125 if (lastStepFailed) {
126 simulator_.model().updateFailed();
127 }
128 else {
129 simulator_.model().advanceTimeLevel();
130 }
131
132 // Set the timestep size and episode index for the model explicitly.
133 // The model needs to know the report step/episode index because of
134 // timing dependent data despite the fact that flow uses its own time stepper.
135 // (The length of the episode does not matter, though.)
136 simulator_.setTime(timer.simulationTimeElapsed());
137 simulator_.setTimeStepSize(timer.currentStepLength());
138
139 simulator_.problem().resetIterationForNewTimestep();
140
141 simulator_.problem().beginTimeStep();
142
143 unsigned numDof = simulator_.model().numGridDof();
144 wasSwitched_.resize(numDof);
145 std::fill(wasSwitched_.begin(), wasSwitched_.end(), false);
146
147 if (param_.update_equations_scaling_) {
148 OpmLog::error("Equation scaling not supported");
149 //updateEquationsScaling();
150 }
151
152 if (hasNlddSolver()) {
153 nlddSolver_->prepareStep();
154 }
155
156 report.pre_post_time += perfTimer.stop();
157
158 auto getIdx = [](unsigned phaseIdx) -> int
159 {
160 if (FluidSystem::phaseIsActive(phaseIdx)) {
161 const unsigned sIdx = FluidSystem::solventComponentIndex(phaseIdx);
162 return FluidSystem::canonicalToActiveCompIdx(sIdx);
163 }
164
165 return -1;
166 };
167 const auto& schedule = simulator_.vanguard().schedule();
168 auto& rst_conv = simulator_.problem().eclWriter().mutableOutputModule().getConv();
169 rst_conv.init(simulator_.vanguard().globalNumCells(),
170 schedule[timer.reportStepNum()].rst_config(),
171 {getIdx(FluidSystem::oilPhaseIdx),
172 getIdx(FluidSystem::gasPhaseIdx),
173 getIdx(FluidSystem::waterPhaseIdx),
174 contiPolymerEqIdx,
175 contiBrineEqIdx,
176 contiSolventEqIdx});
177
178 return report;
179}
180
181template <class TypeTag>
182void
185 const int minIter,
186 const int maxIter,
187 const SimulatorTimerInterface& timer)
188{
189 // ----------- Set up reports and timer -----------
190 failureReport_ = SimulatorReportSingle();
191 Dune::Timer perfTimer;
192
193 perfTimer.start();
194 report.total_linearizations = 1;
195
196 // ----------- Assemble -----------
197 try {
198 report += assembleReservoir(timer);
199 report.assemble_time += perfTimer.stop();
200 // Mark timestep initialized after assembleReservoir(), because the well model's
201 // assemble() also uses needsTimestepInit() to trigger prepareTimeStep().
202 simulator_.problem().markTimestepInitialized();
203 }
204 catch (...) {
205 report.assemble_time += perfTimer.stop();
206 failureReport_ += report;
207 throw; // continue throwing the stick
208 }
209
210 // ----------- Check if converged -----------
211 std::vector<Scalar> residual_norms;
212 perfTimer.reset();
213 perfTimer.start();
214 // the step is not considered converged until at least minIter iterations is done
215 {
216 auto convrep = getConvergence(timer, maxIter, residual_norms);
217 report.converged = convrep.converged() &&
218 simulator_.problem().iterationContext().iteration() >= minIter;
219 ConvergenceReport::Severity severity = convrep.severityOfWorstFailure();
220 convergence_reports_.back().report.push_back(std::move(convrep));
221
222 // Throw if any NaN or too large residual found.
224 failureReport_ += report;
225 OPM_THROW_PROBLEM(NumericalProblem, "NaN residual found!");
226 } else if (severity == ConvergenceReport::Severity::TooLarge) {
227 failureReport_ += report;
228 OPM_THROW_NOLOG(NumericalProblem, "Too large residual found!");
230 failureReport_ += report;
231 OPM_THROW_PROBLEM(ConvergenceMonitorFailure,
232 fmt::format(
233 "Total penalty count exceeded cut-off-limit of {}",
234 param_.monitor_params_.cutoff_
235 ));
236 }
237 }
238 report.update_time += perfTimer.stop();
239 residual_norms_history_.push_back(residual_norms);
240}
241
242template <class TypeTag>
243template <class NonlinearSolverType>
247 NonlinearSolverType& nonlinear_solver)
248{
249 // Model-level timestep initialization (once per timestep).
250 // markTimestepInitialized() is called later in initialLinearization(),
251 // after assembleReservoir() has triggered the well model's prepareTimeStep().
252 if (simulator_.problem().iterationContext().needsTimestepInit()) {
253 residual_norms_history_.clear();
254 conv_monitor_.reset();
255 current_relaxation_ = 1.0;
256 dx_old_ = 0.0;
257 convergence_reports_.push_back({timer.reportStepNum(), timer.currentStepNum(), {}});
258 convergence_reports_.back().report.reserve(11);
259 }
260
262 if (this->param_.nonlinear_solver_ != "nldd")
263 {
264 result = this->nonlinearIterationNewton(timer, nonlinear_solver);
265 }
266 else {
267 result = this->nlddSolver_->nonlinearIterationNldd(timer, nonlinear_solver);
268 }
269
270 auto& rst_conv = simulator_.problem().eclWriter().mutableOutputModule().getConv();
271 rst_conv.update(simulator_.model().linearizer().residual());
272
273 simulator_.problem().advanceIteration();
274 return result;
275}
276
277template <class TypeTag>
278template <class NonlinearSolverType>
282 NonlinearSolverType& nonlinear_solver)
283{
284 OPM_TIMEFUNCTION();
285
286 // ----------- Set up reports and timer -----------
288 Dune::Timer perfTimer;
289
290 this->initialLinearization(report,
291 this->param_.newton_min_iter_,
292 this->param_.newton_max_iter_,
293 timer);
294
295 // ----------- If not converged, solve linear system and do Newton update -----------
296 if (!report.converged) {
297 perfTimer.reset();
298 perfTimer.start();
299 report.total_newton_iterations = 1;
300
301 // Compute the nonlinear update.
302 unsigned nc = simulator_.model().numGridDof();
303 BVector x(nc);
304
305 // Solve the linear system.
306 linear_solve_setup_time_ = 0.0;
307 try {
308 // Apply the Schur complement of the well model to
309 // the reservoir linearized equations.
310 // Note that linearize may throw for MSwells.
311 wellModel().linearize(simulator().model().linearizer().jacobian(),
312 simulator().model().linearizer().residual());
313
314 // ---- Solve linear system ----
315 solveJacobianSystem(x);
316
317 report.linear_solve_setup_time += linear_solve_setup_time_;
318 report.linear_solve_time += perfTimer.stop();
319 report.total_linear_iterations += linearIterationsLastSolve();
320 }
321 catch (...) {
322 report.linear_solve_setup_time += linear_solve_setup_time_;
323 report.linear_solve_time += perfTimer.stop();
324 report.total_linear_iterations += linearIterationsLastSolve();
325
326 failureReport_ += report;
327 throw; // re-throw up
328 }
329
330 perfTimer.reset();
331 perfTimer.start();
332
333 // handling well state update before oscillation treatment is a decision based
334 // on observation to avoid some big performance degeneration under some circumstances.
335 // there is no theorectical explanation which way is better for sure.
336 wellModel().postSolve(x);
337
338 if (param_.use_update_stabilization_) {
339 // Stabilize the nonlinear update.
340 bool isOscillate = false;
341 bool isStagnate = false;
342 nonlinear_solver.detectOscillations(residual_norms_history_,
343 residual_norms_history_.size() - 1,
344 isOscillate,
345 isStagnate);
346 if (isOscillate) {
347 current_relaxation_ -= nonlinear_solver.relaxIncrement();
348 current_relaxation_ = std::max(current_relaxation_,
349 nonlinear_solver.relaxMax());
350 if (terminalOutputEnabled()) {
351 std::string msg = " Oscillating behavior detected: Relaxation set to "
352 + std::to_string(current_relaxation_);
353 OpmLog::info(msg);
354 }
355 }
356 nonlinear_solver.stabilizeNonlinearUpdate(x, dx_old_, current_relaxation_);
357 }
358
359 // ---- Newton update ----
360 // Apply the update, with considering model-dependent limitations and
361 // chopping of the update.
362 updateSolution(x);
363
364 report.update_time += perfTimer.stop();
365 }
366
367 return report;
368}
369
370template <class TypeTag>
374{
375 // -------- Mass balance equations --------
376 simulator_.problem().beginIteration();
377 simulator_.model().linearizer().linearizeDomain();
378 simulator_.problem().endIteration();
379 return wellModel().lastReport();
380}
381
382template <class TypeTag>
385relativeChange() const
386{
387 Scalar resultDelta = 0.0;
388 Scalar resultDenom = 0.0;
389
390 const auto& elemMapper = simulator_.model().elementMapper();
391 const auto& gridView = simulator_.gridView();
392 for (const auto& elem : elements(gridView, Dune::Partitions::interior)) {
393 unsigned globalElemIdx = elemMapper.index(elem);
394 const auto& priVarsNew = simulator_.model().solution(/*timeIdx=*/0)[globalElemIdx];
395
396 Scalar pressureNew;
397 pressureNew = priVarsNew[Indices::pressureSwitchIdx];
398
399 Scalar saturationsNew[FluidSystem::numPhases] = { 0.0 };
400 Scalar oilSaturationNew = 1.0;
401 if (FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx) &&
402 FluidSystem::numActivePhases() > 1 &&
403 priVarsNew.primaryVarsMeaningWater() == PrimaryVariables::WaterMeaning::Sw)
404 {
405 saturationsNew[FluidSystem::waterPhaseIdx] = priVarsNew[Indices::waterSwitchIdx];
406 oilSaturationNew -= saturationsNew[FluidSystem::waterPhaseIdx];
407 }
408
409 if (FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx) &&
410 FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx) &&
411 priVarsNew.primaryVarsMeaningGas() == PrimaryVariables::GasMeaning::Sg)
412 {
413 assert(Indices::compositionSwitchIdx >= 0);
414 saturationsNew[FluidSystem::gasPhaseIdx] = priVarsNew[Indices::compositionSwitchIdx];
415 oilSaturationNew -= saturationsNew[FluidSystem::gasPhaseIdx];
416 }
417
418 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx)) {
419 saturationsNew[FluidSystem::oilPhaseIdx] = oilSaturationNew;
420 }
421
422 const auto& priVarsOld = simulator_.model().solution(/*timeIdx=*/1)[globalElemIdx];
423
424 Scalar pressureOld;
425 pressureOld = priVarsOld[Indices::pressureSwitchIdx];
426
427 Scalar saturationsOld[FluidSystem::numPhases] = { 0.0 };
428 Scalar oilSaturationOld = 1.0;
429
430 // NB fix me! adding pressures changes to saturation changes does not make sense
431 Scalar tmp = pressureNew - pressureOld;
432 resultDelta += tmp*tmp;
433 resultDenom += pressureNew*pressureNew;
434
435 if (FluidSystem::numActivePhases() > 1) {
436 if (priVarsOld.primaryVarsMeaningWater() == PrimaryVariables::WaterMeaning::Sw) {
437 saturationsOld[FluidSystem::waterPhaseIdx] =
438 priVarsOld[Indices::waterSwitchIdx];
439 oilSaturationOld -= saturationsOld[FluidSystem::waterPhaseIdx];
440 }
441
442 if (priVarsOld.primaryVarsMeaningGas() == PrimaryVariables::GasMeaning::Sg)
443 {
444 assert(Indices::compositionSwitchIdx >= 0 );
445 saturationsOld[FluidSystem::gasPhaseIdx] =
446 priVarsOld[Indices::compositionSwitchIdx];
447 oilSaturationOld -= saturationsOld[FluidSystem::gasPhaseIdx];
448 }
449
450 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx)) {
451 saturationsOld[FluidSystem::oilPhaseIdx] = oilSaturationOld;
452 }
453 for (unsigned phaseIdx = 0; phaseIdx < FluidSystem::numPhases; ++ phaseIdx) {
454 Scalar tmpSat = saturationsNew[phaseIdx] - saturationsOld[phaseIdx];
455 resultDelta += tmpSat*tmpSat;
456 resultDenom += saturationsNew[phaseIdx]*saturationsNew[phaseIdx];
457 assert(std::isfinite(resultDelta));
458 assert(std::isfinite(resultDenom));
459 }
460 }
461 }
462
463 resultDelta = gridView.comm().sum(resultDelta);
464 resultDenom = gridView.comm().sum(resultDenom);
465
466 return resultDenom > 0.0 ? resultDelta / resultDenom : 0.0;
467}
468
469template <class TypeTag>
470void
473{
474 auto& jacobian = simulator_.model().linearizer().jacobian().istlMatrix();
475 auto& residual = simulator_.model().linearizer().residual();
476 auto& linSolver = simulator_.model().newtonMethod().linearSolver();
477
478 const int numSolvers = linSolver.numAvailableSolvers();
479 if (numSolvers > 1 && (linSolver.getSolveCount() % 100 == 0)) {
480 if (terminal_output_) {
481 OpmLog::debug("\nRunning speed test for comparing available linear solvers.");
482 }
483
484 Dune::Timer perfTimer;
485 std::vector<double> times(numSolvers);
486 std::vector<double> setupTimes(numSolvers);
487
488 x = 0.0;
489 std::vector<BVector> x_trial(numSolvers, x);
490 for (int solver = 0; solver < numSolvers; ++solver) {
491 BVector x0(x);
492 linSolver.setActiveSolver(solver);
493 perfTimer.start();
494 linSolver.prepare(jacobian, residual);
495 setupTimes[solver] = perfTimer.stop();
496 perfTimer.reset();
497 linSolver.setResidual(residual);
498 perfTimer.start();
499 linSolver.solve(x_trial[solver]);
500 times[solver] = perfTimer.stop();
501 perfTimer.reset();
502 if (terminal_output_) {
503 OpmLog::debug(fmt::format(fmt::runtime("Solver time {}: {}"), solver, times[solver]));
504 }
505 }
506
507 int fastest_solver = std::ranges::min_element(times) - times.begin();
508 // Use timing on rank 0 to determine fastest, must be consistent across ranks.
509 grid_.comm().broadcast(&fastest_solver, 1, 0);
510 linear_solve_setup_time_ = setupTimes[fastest_solver];
511 x = x_trial[fastest_solver];
512 linSolver.setActiveSolver(fastest_solver);
513 }
514 else {
515 // set initial guess
516 x = 0.0;
517
518 Dune::Timer perfTimer;
519 perfTimer.start();
520 linSolver.prepare(jacobian, residual);
521 linear_solve_setup_time_ = perfTimer.stop();
522 linSolver.setResidual(residual);
523 // actually, the error needs to be calculated after setResidual in order to
524 // account for parallelization properly. since the residual of ECFV
525 // discretizations does not need to be synchronized across processes to be
526 // consistent, this is not relevant for OPM-flow...
527 linSolver.solve(x);
528 }
529}
530
531template <class TypeTag>
532void
534updateSolution(const BVector& dx)
535{
536 OPM_TIMEBLOCK(updateSolution);
537 // Prepare to store solution update for convergence check
538 if (this->param_.tolerance_max_dp_ > 0.0 || this->param_.tolerance_max_ds_ > 0.0
539 || this->param_.tolerance_max_drs_ > 0.0 || this->param_.tolerance_max_drv_ > 0.0) {
540 prepareStoringSolutionUpdate();
541 }
542
543 auto& newtonMethod = simulator_.model().newtonMethod();
544 SolutionVector& solution = simulator_.model().solution(/*timeIdx=*/0);
545
546 newtonMethod.update_(/*nextSolution=*/solution,
547 /*curSolution=*/solution,
548 /*update=*/dx,
549 /*resid=*/dx); // the update routines of the black
550 // oil model do not care about the
551 // residual
552
553 // if the solution is updated, the intensive quantities need to be recalculated
554 {
555 OPM_TIMEBLOCK(invalidateAndUpdateIntensiveQuantities);
556 simulator_.model().invalidateAndUpdateIntensiveQuantities(/*timeIdx=*/0);
557 }
558
559 // Store solution update
560 if (this->param_.tolerance_max_dp_ > 0.0 || this->param_.tolerance_max_ds_ > 0.0
561 || this->param_.tolerance_max_drs_ > 0.0 || this->param_.tolerance_max_drv_ > 0.0) {
562 storeSolutionUpdate(dx);
563 }
564}
565
566template <class TypeTag>
567void
570{
571 // Init. solution update vector
572 unsigned nc = simulator_.model().numGridDof();
573 solUpd_.resize(nc);
574
575 const auto& elemMapper = simulator_.model().elementMapper();
576 const auto& gridView = simulator_.gridView();
577 for (const auto& elem : elements(gridView, Dune::Partitions::interior)) {
578 // Copy solution vector to transfer primary variables meaning
579 unsigned globalElemIdx = elemMapper.index(elem);
580 solUpd_[globalElemIdx] = simulator_.model().solution(/*timeIdx=*/0)[globalElemIdx];
581
582 // Ensure each element is zero
583 std::ranges::fill(solUpd_[globalElemIdx], 0.0);
584 }
585}
586
587template <class TypeTag>
588void
591{
592 const auto& elemMapper = simulator_.model().elementMapper();
593 const auto& gridView = simulator_.gridView();
594 for (const auto& elem : elements(gridView, Dune::Partitions::interior)) {
595 // Get cell vectors
596 unsigned globalElemIdx = elemMapper.index(elem);
597 auto& value = solUpd_[globalElemIdx];
598 const auto& update = dx[globalElemIdx];
599 assert(value.size() == update.size());
600
601 // Transfer update from dx to solution update container (SolutionVector type)
602 std::ranges::copy(update, value.begin());
603 }
604}
605
606template <class TypeTag>
609getMaxSolutionUpdate(const std::vector<unsigned>& ixCells)
610{
611 static constexpr bool enableSolvent = Indices::solventSaturationIdx >= 0;
612 static constexpr bool enableBrine = Indices::saltConcentrationIdx >= 0;
613
614 // Init output
615 Scalar dPMax = 0.0;
616 Scalar dSMax = 0.0;
617 Scalar dRsMax = 0.0;
618 Scalar dRvMax = 0.0;
619
620 // Loop over solution update, get the correct variables and calculate max.
621 for (const auto& ix : ixCells) {
622 const auto& value = solUpd_[ix];
623 for (int pvIdx = 0; pvIdx < static_cast<int>(value.size()); ++pvIdx) {
624 if (pvIdx == Indices::pressureSwitchIdx) {
625 dPMax = std::max(dPMax, std::abs(value[pvIdx]));
626 }
627 else if ( (pvIdx == Indices::waterSwitchIdx
628 && value.primaryVarsMeaningWater() == PrimaryVariables::WaterMeaning::Sw)
629 || (pvIdx == Indices::compositionSwitchIdx
630 && value.primaryVarsMeaningGas() == PrimaryVariables::GasMeaning::Sg)
631 || (enableSolvent && pvIdx == Indices::solventSaturationIdx
632 && value.primaryVarsMeaningSolvent() == PrimaryVariables::SolventMeaning::Ss)
633 || (enableBrine && enableSaltPrecipitation && pvIdx == Indices::saltConcentrationIdx
634 && value.primaryVarsMeaningBrine() == PrimaryVariables::BrineMeaning::Sp) ) {
635 dSMax = std::max(dSMax, std::abs(value[pvIdx]));
636 }
637 else if (pvIdx == Indices::compositionSwitchIdx
638 && value.primaryVarsMeaningGas() == PrimaryVariables::GasMeaning::Rs) {
639 dRsMax = std::max(dRsMax, std::abs(value[pvIdx]));
640 }
641 else if (pvIdx == Indices::compositionSwitchIdx
642 && value.primaryVarsMeaningGas() == PrimaryVariables::GasMeaning::Rv) {
643 dRvMax = std::max(dRvMax, std::abs(value[pvIdx]));
644 }
645 }
646 }
647
648 // Communicate max values
649 dPMax = this->grid_.comm().max(dPMax);
650 dSMax = this->grid_.comm().max(dSMax);
651 dRsMax = this->grid_.comm().max(dRsMax);
652 dRvMax = this->grid_.comm().max(dRvMax);
653
654 return { dPMax, dSMax, dRsMax, dRvMax };
655}
656
657template <class TypeTag>
658std::tuple<typename BlackoilModel<TypeTag>::Scalar,
662 const Scalar pvSumLocal,
663 const Scalar numAquiferPvSumLocal,
664 std::vector< Scalar >& R_sum,
665 std::vector< Scalar >& maxCoeff,
666 std::vector< Scalar >& B_avg)
667{
668 OPM_TIMEBLOCK(convergenceReduction);
669 // Compute total pore volume (use only owned entries)
670 Scalar pvSum = pvSumLocal;
671 Scalar numAquiferPvSum = numAquiferPvSumLocal;
672
673 if (comm.size() > 1) {
674 // global reduction
675 std::vector< Scalar > sumBuffer;
676 std::vector< Scalar > maxBuffer;
677 const int numComp = B_avg.size();
678 sumBuffer.reserve( 2*numComp + 2 ); // +2 for (numAquifer)pvSum
679 maxBuffer.reserve( numComp );
680 for (int compIdx = 0; compIdx < numComp; ++compIdx) {
681 sumBuffer.push_back(B_avg[compIdx]);
682 sumBuffer.push_back(R_sum[compIdx]);
683 maxBuffer.push_back(maxCoeff[ compIdx]);
684 }
685
686 // Compute total pore volume
687 sumBuffer.push_back( pvSum );
688 sumBuffer.push_back( numAquiferPvSum );
689
690 // compute global sum
691 comm.sum( sumBuffer.data(), sumBuffer.size() );
692
693 // compute global max
694 comm.max( maxBuffer.data(), maxBuffer.size() );
695
696 // restore values to local variables
697 for (int compIdx = 0, buffIdx = 0; compIdx < numComp; ++compIdx, ++buffIdx) {
698 B_avg[compIdx] = sumBuffer[buffIdx];
699 ++buffIdx;
700
701 R_sum[compIdx] = sumBuffer[buffIdx];
702 }
703
704 for (int compIdx = 0; compIdx < numComp; ++compIdx) {
705 maxCoeff[compIdx] = maxBuffer[compIdx];
706 }
707
708 // restore global pore volume
709 pvSum = sumBuffer[sumBuffer.size()-2];
710 numAquiferPvSum = sumBuffer.back();
711 }
712
713 // return global pore volume
714 return {pvSum, numAquiferPvSum};
715}
716
717template <class TypeTag>
718std::pair<typename BlackoilModel<TypeTag>::Scalar,
721localConvergenceData(std::vector<Scalar>& R_sum,
722 std::vector<Scalar>& maxCoeff,
723 std::vector<Scalar>& B_avg,
724 std::vector<int>& maxCoeffCell)
725{
726 OPM_TIMEBLOCK(localConvergenceData);
727 Scalar pvSumLocal = 0.0;
728 Scalar numAquiferPvSumLocal = 0.0;
729 const auto& model = simulator_.model();
730 const auto& problem = simulator_.problem();
731
732 const auto& residual = simulator_.model().linearizer().residual();
733
734 ElementContext elemCtx(simulator_);
735 const auto& gridView = simulator().gridView();
736 IsNumericalAquiferCell isNumericalAquiferCell(gridView.grid());
738 for (const auto& elem : elements(gridView, Dune::Partitions::interior)) {
739 elemCtx.updatePrimaryStencil(elem);
740 elemCtx.updatePrimaryIntensiveQuantities(/*timeIdx=*/0);
741
742 const unsigned cell_idx = elemCtx.globalSpaceIndex(/*spaceIdx=*/0, /*timeIdx=*/0);
743 const auto& intQuants = elemCtx.intensiveQuantities(/*spaceIdx=*/0, /*timeIdx=*/0);
744 const auto& fs = intQuants.fluidState();
745
746 const auto pvValue = problem.referencePorosity(cell_idx, /*timeIdx=*/0) *
747 model.dofTotalVolume(cell_idx);
748 pvSumLocal += pvValue;
749
750 if (isNumericalAquiferCell(elem)) {
751 numAquiferPvSumLocal += pvValue;
752 }
753
754 this->getMaxCoeff(cell_idx, intQuants, fs, residual, pvValue,
755 B_avg, R_sum, maxCoeff, maxCoeffCell);
756 }
757
758 OPM_END_PARALLEL_TRY_CATCH("BlackoilModel::localConvergenceData() failed: ", grid_.comm());
759
760 // compute local average in terms of global number of elements
761 const int bSize = B_avg.size();
762 for (int i = 0; i < bSize; ++i) {
763 B_avg[i] /= Scalar(global_nc_);
764 }
765
766 return {pvSumLocal, numAquiferPvSumLocal};
767}
768
769template <class TypeTag>
772characteriseCnvPvSplit(const std::vector<Scalar>& B_avg, const double dt)
773{
774 OPM_TIMEBLOCK(computeCnvErrorPv);
775
776 // 0: cnv <= tolerance_cnv
777 // 1: tolerance_cnv < cnv <= tolerance_cnv_relaxed
778 // 2: tolerance_cnv_relaxed < cnv
779 constexpr auto numPvGroups = std::vector<double>::size_type{3};
780
781 auto cnvPvSplit = std::pair<std::vector<double>, std::vector<int>> {
782 std::piecewise_construct,
783 std::forward_as_tuple(numPvGroups),
784 std::forward_as_tuple(numPvGroups)
785 };
786
787 auto maxCNV = [&B_avg, dt](const auto& residual, const double pvol)
788 {
789 return (dt / pvol) *
790 std::inner_product(residual.begin(), residual.end(),
791 B_avg.begin(), Scalar{0},
792 [](const Scalar m, const auto& x)
793 {
794 using std::abs;
795 return std::max(m, abs(x));
796 }, std::multiplies<>{});
797 };
798
799 auto& [splitPV, cellCntPV] = cnvPvSplit;
800
801 const auto& model = this->simulator().model();
802 const auto& problem = this->simulator().problem();
803 const auto& residual = model.linearizer().residual();
804 const auto& gridView = this->simulator().gridView();
805
806 const IsNumericalAquiferCell isNumericalAquiferCell(gridView.grid());
807
808 ElementContext elemCtx(this->simulator());
809
810 std::vector<unsigned> ixCells;
811
813 for (const auto& elem : elements(gridView, Dune::Partitions::interior)) {
814 // Skip cells of numerical Aquifer
815 if (isNumericalAquiferCell(elem)) {
816 continue;
817 }
818
819 elemCtx.updatePrimaryStencil(elem);
820
821 const unsigned cell_idx = elemCtx.globalSpaceIndex(/*spaceIdx=*/0, /*timeIdx=*/0);
822 const auto pvValue = problem.referencePorosity(cell_idx, /*timeIdx=*/0)
823 * model.dofTotalVolume(cell_idx);
824
825 const auto maxCnv = maxCNV(residual[cell_idx], pvValue);
826
827 const auto ix = (maxCnv > this->param_.tolerance_cnv_)
828 + (maxCnv > this->param_.tolerance_cnv_relaxed_);
829
830 splitPV[ix] += static_cast<double>(pvValue);
831 ++cellCntPV[ix];
832
833 // For dP and dS check, we need cell indices of [1] violations
834 if ( ix > 0 &&
835 (this->param_.tolerance_max_dp_ > 0.0 || this->param_.tolerance_max_ds_ > 0.0
836 || this->param_.tolerance_max_drs_ > 0.0 || this->param_.tolerance_max_drv_ > 0.0 ) ) {
837 ixCells.push_back(cell_idx);
838 }
839 }
840
841 OPM_END_PARALLEL_TRY_CATCH("BlackoilModel::characteriseCnvPvSplit() failed: ",
842 this->grid_.comm());
843
844 this->grid_.comm().sum(splitPV .data(), splitPV .size());
845 this->grid_.comm().sum(cellCntPV.data(), cellCntPV.size());
846
847 return { cnvPvSplit, ixCells };
848}
849
850template <class TypeTag>
851void
853updateTUNING(const Tuning& tuning)
854{
855 this->param_.tolerance_cnv_ = tuning.TRGCNV;
856 this->param_.tolerance_cnv_relaxed_ = tuning.XXXCNV;
857 this->param_.tolerance_mb_ = tuning.TRGMBE;
858 this->param_.tolerance_mb_relaxed_ = tuning.XXXMBE;
859 this->param_.newton_max_iter_ = tuning.NEWTMX;
860 this->param_.newton_min_iter_ = tuning.NEWTMN;
861}
862
863template <class TypeTag>
864void
866updateTUNINGDP(const TuningDp& tuning_dp)
867{
868 // NOTE: If TUNINGDP item is _not_ set it should be 0.0
869 this->param_.tolerance_max_dp_ = tuning_dp.TRGDDP;
870 this->param_.tolerance_max_ds_ = tuning_dp.TRGDDS;
871 this->param_.tolerance_max_drs_ = tuning_dp.TRGDDRS;
872 this->param_.tolerance_max_drv_ = tuning_dp.TRGDDRV;
873}
874
875template <class TypeTag>
878getReservoirConvergence(const double reportTime,
879 const double dt,
880 const int maxIter,
881 std::vector<Scalar>& B_avg,
882 std::vector<Scalar>& residual_norms)
883{
884 OPM_TIMEBLOCK(getReservoirConvergence);
885 using Vector = std::vector<Scalar>;
886
887 const auto& iterCtx = simulator_.problem().iterationContext();
888
889 ConvergenceReport report{reportTime};
890
891 const int numComp = numEq;
892
893 Vector R_sum(numComp, Scalar{0});
894 Vector maxCoeff(numComp, std::numeric_limits<Scalar>::lowest());
895 std::vector<int> maxCoeffCell(numComp, -1);
896
897 const auto [pvSumLocal, numAquiferPvSumLocal] =
898 this->localConvergenceData(R_sum, maxCoeff, B_avg, maxCoeffCell);
899
900 // compute global sum and max of quantities
901 const auto& [pvSum, numAquiferPvSum] =
902 this->convergenceReduction(this->grid_.comm(),
903 pvSumLocal,
904 numAquiferPvSumLocal,
905 R_sum, maxCoeff, B_avg);
906
907 auto cnvSplitData = this->characteriseCnvPvSplit(B_avg, dt);
908 report.setCnvPoreVolSplit(cnvSplitData.cnvPvSplit,
909 pvSum - numAquiferPvSum);
910
911 // For each iteration, we need to determine whether to use the
912 // relaxed tolerances. To disable the usage of relaxed
913 // tolerances, you can set the relaxed tolerances as the strict
914 // tolerances. If min_strict_mb_iter = -1 (default) we use a
915 // relaxed tolerance for the mass balance for the last
916 // iterations. For positive values we use the relaxed tolerance
917 // after the given number of iterations
918 const bool relax_final_iteration_mb =
919 this->param_.min_strict_mb_iter_ < 0 && iterCtx.iteration() == maxIter;
920
921 const bool relax_iter_mb = this->param_.min_strict_mb_iter_ >= 0 &&
922 iterCtx.shouldRelax(this->param_.min_strict_mb_iter_);
923
924 const bool use_relaxed_mb = relax_final_iteration_mb
925 || relax_iter_mb;
926
927 // If min_strict_cnv_iter = -1 we use a relaxed tolerance for
928 // the cnv for the last iterations. For positive values we use
929 // the relaxed tolerance after the given number of iterations.
930 // We also use relaxed tolerances for cells with total
931 // pore-volume less than relaxed_max_pv_fraction_. Default
932 // value of relaxed_max_pv_fraction_ is 0.03
933 const bool relax_final_iteration_cnv =
934 this->param_.min_strict_cnv_iter_ < 0 && iterCtx.iteration() == maxIter;
935
936 const bool relax_iter_cnv = this->param_.min_strict_cnv_iter_ >= 0 &&
937 iterCtx.shouldRelax(this->param_.min_strict_cnv_iter_);
938
939 // Note trailing parentheses here, just before the final
940 // semicolon. This is an immediately invoked function
941 // expression which calculates a single boolean value.
942 const auto relax_pv_fraction_cnv =
943 [&report, this, eligible = pvSum - numAquiferPvSum]()
944 {
945 const auto& cnvPvSplit = report.cnvPvSplit().first;
946
947 // [1]: tol < cnv <= relaxed
948 // [2]: relaxed < cnv
949 Scalar cnvPvSum = static_cast<Scalar>(cnvPvSplit[1] + cnvPvSplit[2]);
950 return cnvPvSum < this->param_.relaxed_max_pv_fraction_ * eligible &&
951 cnvPvSum > 0.0;
952 }();
953
954 // If tolerances for solution changes are met, we use the
955 // relaxed cnv tolerance. Note that all tolerances > 0.0
956 // must be met to enable use of relaxed cnv tolerance.
957 MaxSolutionUpdateData maxSolUpd;
958 const bool use_dp_tol = this->param_.tolerance_max_dp_ > 0.0;
959 const bool use_ds_tol = this->param_.tolerance_max_ds_ > 0.0;
960 const bool use_drs_tol = this->param_.tolerance_max_drs_ > 0.0;
961 const bool use_drv_tol = this->param_.tolerance_max_drv_ > 0.0;
962 const bool use_dsol_tol = use_dp_tol || use_ds_tol || use_drs_tol || use_drv_tol;
963 bool relax_dsol_cnv = false;
964 if (!iterCtx.isFirstGlobalIteration() && use_dsol_tol) {
965 maxSolUpd = getMaxSolutionUpdate(cnvSplitData.ixCells);
966 relax_dsol_cnv =
967 (!use_dp_tol || (maxSolUpd.dPMax > 0.0 && maxSolUpd.dPMax < this->param_.tolerance_max_dp_)) &&
968 (!use_ds_tol || (maxSolUpd.dSMax > 0.0 && maxSolUpd.dSMax < this->param_.tolerance_max_ds_)) &&
969 (!use_drs_tol || (maxSolUpd.dRsMax > 0.0 && maxSolUpd.dRsMax < this->param_.tolerance_max_drs_)) &&
970 (!use_drv_tol || (maxSolUpd.dRvMax > 0.0 && maxSolUpd.dRvMax < this->param_.tolerance_max_drv_));
971 }
972
973 // Determine if relaxed CNV tolerances should be used
974 const bool use_relaxed_cnv = relax_final_iteration_cnv
975 || relax_pv_fraction_cnv
976 || relax_iter_cnv
977 || relax_dsol_cnv;
978
979 // Ensure that CNV convergence criteria is met when max.
980 // solution change tolerances have been fulfilled
981 Scalar tolerance_cnv_relaxed = relax_dsol_cnv ? 1e20 : param_.tolerance_cnv_relaxed_;
982
983 const auto tol_cnv = use_relaxed_cnv ? tolerance_cnv_relaxed : param_.tolerance_cnv_;
984 const auto tol_mb = use_relaxed_mb ? param_.tolerance_mb_relaxed_ : param_.tolerance_mb_;
985 const auto tol_cnv_energy = use_relaxed_cnv ? param_.tolerance_cnv_energy_relaxed_ : param_.tolerance_cnv_energy_;
986 const auto tol_eb = use_relaxed_mb ? param_.tolerance_energy_balance_relaxed_ : param_.tolerance_energy_balance_;
987
988 // Finish computation
989 std::vector<Scalar> CNV(numComp);
990 std::vector<Scalar> mass_balance_residual(numComp);
991 for (int compIdx = 0; compIdx < numComp; ++compIdx)
992 {
993 CNV[compIdx] = B_avg[compIdx] * dt * maxCoeff[compIdx];
994 mass_balance_residual[compIdx] = std::abs(B_avg[compIdx]*R_sum[compIdx]) * dt / pvSum;
995 residual_norms.push_back(CNV[compIdx]);
996 }
997
998 using CR = ConvergenceReport;
999 for (int compIdx = 0; compIdx < numComp; ++compIdx) {
1000 const Scalar res[2] = {
1001 mass_balance_residual[compIdx], CNV[compIdx],
1002 };
1003
1004 const CR::ReservoirFailure::Type types[2] = {
1005 CR::ReservoirFailure::Type::MassBalance,
1006 CR::ReservoirFailure::Type::Cnv,
1007 };
1008
1009 Scalar tol[2] = { tol_mb, tol_cnv, };
1010 if (has_energy_ && compIdx == contiEnergyEqIdx) {
1011 tol[0] = tol_eb;
1012 tol[1] = tol_cnv_energy;
1013 }
1014
1015 for (int ii : {0, 1}) {
1016 if (std::isnan(res[ii])) {
1017 report.setReservoirFailed({types[ii], CR::Severity::NotANumber, compIdx});
1018 if (this->terminal_output_) {
1019 OpmLog::debug("NaN residual for " + this->compNames_.name(compIdx) + " equation.");
1020 }
1021 }
1022 else if (res[ii] > maxResidualAllowed()) {
1023 report.setReservoirFailed({types[ii], CR::Severity::TooLarge, compIdx});
1024 if (this->terminal_output_) {
1025 OpmLog::debug("Too large residual for " + this->compNames_.name(compIdx) + " equation.");
1026 }
1027 }
1028 else if (res[ii] < 0.0) {
1029 report.setReservoirFailed({types[ii], CR::Severity::Normal, compIdx});
1030 if (this->terminal_output_) {
1031 OpmLog::debug("Negative residual for " + this->compNames_.name(compIdx) + " equation.");
1032 }
1033 }
1034 else if (res[ii] > tol[ii]) {
1035 report.setReservoirFailed({types[ii], CR::Severity::Normal, compIdx});
1036 }
1037
1038 report.setReservoirConvergenceMetric(types[ii], compIdx, res[ii], tol[ii]);
1039 }
1040 }
1041
1042 // Compute the Newton convergence per cell.
1043 this->convergencePerCell(B_avg, dt, tol_cnv, tol_cnv_energy);
1044
1045 // Output of residuals.
1046 if (this->terminal_output_) {
1047 // Only rank 0 does print to std::cout
1048 if (iterCtx.isFirstGlobalIteration()) {
1049 std::string msg = "Iter";
1050 for (int compIdx = 0; compIdx < numComp; ++compIdx) {
1051 msg += " MB(";
1052 msg += this->compNames_.name(compIdx)[0];
1053 msg += ") ";
1054 }
1055
1056 for (int compIdx = 0; compIdx < numComp; ++compIdx) {
1057 msg += " CNV(";
1058 msg += this->compNames_.name(compIdx)[0];
1059 msg += ") ";
1060 }
1061
1062 if (use_dsol_tol) {
1063 msg += use_dp_tol ? " DP " : "";
1064 msg += use_ds_tol ? " DS " : "";
1065 msg += use_drs_tol ? " DRS " : "";
1066 msg += use_drv_tol ? " DRV " : "";
1067 }
1068
1069 msg += " MBFLAG";
1070 msg += " CNVFLAG";
1071
1072 OpmLog::debug(msg);
1073 }
1074
1075 std::ostringstream ss;
1076 const std::streamsize oprec = ss.precision(3);
1077 const std::ios::fmtflags oflags = ss.setf(std::ios::scientific);
1078
1079 ss << std::setw(4) << iterCtx.iteration();
1080 for (int compIdx = 0; compIdx < numComp; ++compIdx) {
1081 ss << std::setw(11) << mass_balance_residual[compIdx];
1082 }
1083
1084 for (int compIdx = 0; compIdx < numComp; ++compIdx) {
1085 ss << std::setw(11) << CNV[compIdx];
1086 }
1087
1088 if (use_dsol_tol) {
1089 auto print_dsol =
1090 [&] (bool use_tol, Scalar dsol) {
1091 if (!use_tol) {
1092 return;
1093 }
1094 if (iterCtx.isFirstGlobalIteration() || dsol <= 0.0) {
1095 ss << std::string(5, ' ') << "-" << std::string(5, ' ');
1096 }
1097 else {
1098 ss << std::setw(11) << dsol;
1099 }
1100 };
1101 print_dsol(use_dp_tol, maxSolUpd.dPMax);
1102 print_dsol(use_ds_tol, maxSolUpd.dSMax);
1103 print_dsol(use_drs_tol, maxSolUpd.dRsMax);
1104 print_dsol(use_drv_tol, maxSolUpd.dRvMax);
1105 }
1106
1107 const auto mb_flag = use_relaxed_mb
1108 ? DebugFlags::RELAXED
1109 : DebugFlags::STRICT;
1110
1111 const auto cnv_flag = relax_dsol_cnv ?
1112 DebugFlags::TUNINGDP
1113 : (use_relaxed_cnv
1114 ? DebugFlags::RELAXED
1115 : DebugFlags::STRICT);
1116
1117 ss << std::setw(9) << make_string<TypeTag>(mb_flag)
1118 << std::setw(9) << make_string<TypeTag>(cnv_flag);
1119
1120 ss.precision(oprec);
1121 ss.flags(oflags);
1122
1123 OpmLog::debug(ss.str());
1124 }
1125
1126 return report;
1127}
1128
1129template <class TypeTag>
1130void
1132convergencePerCell(const std::vector<Scalar>& B_avg,
1133 const double dt,
1134 const double tol_cnv,
1135 const double tol_cnv_energy)
1136{
1137 auto& rst_conv = simulator_.problem().eclWriter().mutableOutputModule().getConv();
1138 if (!rst_conv.hasConv()) {
1139 return;
1140 }
1141
1142 if (simulator_.problem().iterationContext().isFirstGlobalIteration()) {
1143 rst_conv.prepareConv();
1144 }
1145
1146 const auto& residual = simulator_.model().linearizer().residual();
1147 const auto& gridView = this->simulator().gridView();
1148 const IsNumericalAquiferCell isNumericalAquiferCell(gridView.grid());
1149 ElementContext elemCtx(this->simulator());
1150 std::vector<int> convNewt(residual.size(), 0);
1152 unsigned idx = 0;
1153 const int numComp = B_avg.size();
1154 for (const auto& elem : elements(gridView, Dune::Partitions::interior)) {
1155 elemCtx.updatePrimaryStencil(elem);
1156
1157 const unsigned cell_idx = elemCtx.globalSpaceIndex(/*spaceIdx=*/0, /*timeIdx=*/0);
1158 const auto pvValue = simulator_.problem().referencePorosity(cell_idx, /*timeIdx=*/0) *
1159 simulator_.model().dofTotalVolume(cell_idx);
1160 for (int compIdx = 0; compIdx < numComp; ++compIdx) {
1161 const auto tol = (has_energy_ && compIdx == contiEnergyEqIdx) ? tol_cnv_energy : tol_cnv;
1162 const Scalar cnv = std::abs(B_avg[compIdx] * residual[cell_idx][compIdx]) * dt / pvValue;
1163 if (std::isnan(cnv) || cnv > maxResidualAllowed() || cnv < 0.0 || cnv > tol) {
1164 convNewt[idx] = 1;
1165 break;
1166 }
1167 }
1168 ++idx;
1169 }
1170 OPM_END_PARALLEL_TRY_CATCH("BlackoilModel::convergencePerCell() failed: ",
1171 this->grid_.comm());
1172 rst_conv.updateNewton(convNewt);
1173}
1174
1175template <class TypeTag>
1179 const int maxIter,
1180 std::vector<Scalar>& residual_norms)
1181{
1182 OPM_TIMEBLOCK(getConvergence);
1183 // Get convergence reports for reservoir and wells.
1184 std::vector<Scalar> B_avg(numEq, 0.0);
1185 auto report = getReservoirConvergence(timer.simulationTimeElapsed(),
1186 timer.currentStepLength(),
1187 maxIter, B_avg, residual_norms);
1188 {
1189 OPM_TIMEBLOCK(getWellConvergence);
1190 report += wellModel().getWellConvergence(B_avg,
1191 /*checkWellGroupControlsAndNetwork*/report.converged());
1192 }
1193
1194 conv_monitor_.checkPenaltyCard(report, simulator_.problem().iterationContext().iteration());
1195
1196 return report;
1197}
1198
1199template <class TypeTag>
1200std::vector<std::vector<typename BlackoilModel<TypeTag>::Scalar> >
1202computeFluidInPlace(const std::vector<int>& /*fipnum*/) const
1203{
1204 OPM_TIMEBLOCK(computeFluidInPlace);
1205 // assert(true)
1206 // return an empty vector
1207 std::vector<std::vector<Scalar> > regionValues(0, std::vector<Scalar>(0,0.0));
1208 return regionValues;
1209}
1210
1211template <class TypeTag>
1212const SimulatorReport&
1215{
1216 if (!hasNlddSolver()) {
1217 OPM_THROW(std::runtime_error, "Cannot get local reports from a model without NLDD solver");
1218 }
1219 return nlddSolver_->localAccumulatedReports();
1220}
1221
1222template <class TypeTag>
1223const std::vector<SimulatorReport>&
1226{
1227 if (!nlddSolver_)
1228 OPM_THROW(std::runtime_error, "Cannot get domain reports from a model without NLDD solver");
1229 return nlddSolver_->domainAccumulatedReports();
1230}
1231
1232template <class TypeTag>
1233void
1235writeNonlinearIterationsPerCell(const std::filesystem::path& odir) const
1236{
1237 if (hasNlddSolver()) {
1238 nlddSolver_->writeNonlinearIterationsPerCell(odir);
1239 }
1240}
1241
1242template <class TypeTag>
1243void
1245writePartitions(const std::filesystem::path& odir) const
1246{
1247 if (hasNlddSolver()) {
1248 nlddSolver_->writePartitions(odir);
1249 return;
1250 }
1251
1252 const auto& elementMapper = this->simulator().model().elementMapper();
1253 const auto& cartMapper = this->simulator().vanguard().cartesianIndexMapper();
1254
1255 const auto& grid = this->simulator().vanguard().grid();
1256 const auto& comm = grid.comm();
1257 const auto nDigit = 1 + static_cast<int>(std::floor(std::log10(comm.size())));
1258
1259 std::ofstream pfile {odir / fmt::format("{1:0>{0}}", nDigit, comm.rank())};
1260
1261 for (const auto& cell : elements(grid.leafGridView(), Dune::Partitions::interior)) {
1262 pfile << comm.rank() << ' '
1263 << cartMapper.cartesianIndex(elementMapper.index(cell)) << ' '
1264 << comm.rank() << '\n';
1265 }
1266}
1267
1268template <class TypeTag>
1269template<class FluidState, class Residual>
1270void
1272getMaxCoeff(const unsigned cell_idx,
1273 const IntensiveQuantities& intQuants,
1274 const FluidState& fs,
1275 const Residual& modelResid,
1276 const Scalar pvValue,
1277 std::vector<Scalar>& B_avg,
1278 std::vector<Scalar>& R_sum,
1279 std::vector<Scalar>& maxCoeff,
1280 std::vector<int>& maxCoeffCell)
1281{
1282 for (unsigned phaseIdx = 0; phaseIdx < FluidSystem::numPhases; ++phaseIdx)
1283 {
1284 if (!FluidSystem::phaseIsActive(phaseIdx)) {
1285 continue;
1286 }
1287
1288 const unsigned sIdx = FluidSystem::solventComponentIndex(phaseIdx);
1289 const unsigned compIdx = FluidSystem::canonicalToActiveCompIdx(sIdx);
1290
1291 B_avg[compIdx] += 1.0 / fs.invB(phaseIdx).value();
1292 const auto R2 = modelResid[cell_idx][compIdx];
1293
1294 R_sum[compIdx] += R2;
1295 const Scalar Rval = std::abs(R2) / pvValue;
1296 if (Rval > maxCoeff[compIdx]) {
1297 maxCoeff[compIdx] = Rval;
1298 maxCoeffCell[compIdx] = cell_idx;
1299 }
1300 }
1301
1302 if constexpr (has_solvent_) {
1303 B_avg[contiSolventEqIdx] +=
1304 1.0 / intQuants.solventInverseFormationVolumeFactor().value();
1305 const auto R2 = modelResid[cell_idx][contiSolventEqIdx];
1306 R_sum[contiSolventEqIdx] += R2;
1307 maxCoeff[contiSolventEqIdx] = std::max(maxCoeff[contiSolventEqIdx],
1308 std::abs(R2) / pvValue);
1309 }
1310 if constexpr (has_extbo_) {
1311 B_avg[contiZfracEqIdx] += 1.0 / fs.invB(FluidSystem::gasPhaseIdx).value();
1312 const auto R2 = modelResid[cell_idx][contiZfracEqIdx];
1313 R_sum[ contiZfracEqIdx ] += R2;
1314 maxCoeff[contiZfracEqIdx] = std::max(maxCoeff[contiZfracEqIdx],
1315 std::abs(R2) / pvValue);
1316 }
1317 if constexpr (has_polymer_) {
1318 B_avg[contiPolymerEqIdx] += 1.0 / fs.invB(FluidSystem::waterPhaseIdx).value();
1319 const auto R2 = modelResid[cell_idx][contiPolymerEqIdx];
1320 R_sum[contiPolymerEqIdx] += R2;
1321 maxCoeff[contiPolymerEqIdx] = std::max(maxCoeff[contiPolymerEqIdx],
1322 std::abs(R2) / pvValue);
1323 }
1324 if constexpr (has_foam_) {
1325 B_avg[ contiFoamEqIdx ] += 1.0 / fs.invB(FluidSystem::gasPhaseIdx).value();
1326 const auto R2 = modelResid[cell_idx][contiFoamEqIdx];
1327 R_sum[contiFoamEqIdx] += R2;
1328 maxCoeff[contiFoamEqIdx] = std::max(maxCoeff[contiFoamEqIdx],
1329 std::abs(R2) / pvValue);
1330 }
1331 if constexpr (has_brine_) {
1332 B_avg[ contiBrineEqIdx ] += 1.0 / fs.invB(FluidSystem::waterPhaseIdx).value();
1333 const auto R2 = modelResid[cell_idx][contiBrineEqIdx];
1334 R_sum[contiBrineEqIdx] += R2;
1335 maxCoeff[contiBrineEqIdx] = std::max(maxCoeff[contiBrineEqIdx],
1336 std::abs(R2) / pvValue);
1337 }
1338
1339 if constexpr (has_polymermw_) {
1340 static_assert(has_polymer_);
1341
1342 B_avg[contiPolymerMWEqIdx] += 1.0 / fs.invB(FluidSystem::waterPhaseIdx).value();
1343 // the residual of the polymer molecular equation is scaled down by a 100, since molecular weight
1344 // can be much bigger than 1, and this equation shares the same tolerance with other mass balance equations
1345 // TODO: there should be a more general way to determine the scaling-down coefficient
1346 const auto R2 = modelResid[cell_idx][contiPolymerMWEqIdx] / 100.;
1347 R_sum[contiPolymerMWEqIdx] += R2;
1348 maxCoeff[contiPolymerMWEqIdx] = std::max(maxCoeff[contiPolymerMWEqIdx],
1349 std::abs(R2) / pvValue);
1350 }
1351
1352 if constexpr (has_energy_) {
1353 B_avg[contiEnergyEqIdx] += 1.0;
1354 const auto R2 = modelResid[cell_idx][contiEnergyEqIdx];
1355 R_sum[contiEnergyEqIdx] += R2;
1356 maxCoeff[contiEnergyEqIdx] = std::max(maxCoeff[contiEnergyEqIdx],
1357 std::abs(R2) / pvValue);
1358 }
1359
1360 if constexpr (has_bioeffects_) {
1361 B_avg[contiMicrobialEqIdx] += 1.0 / fs.invB(FluidSystem::waterPhaseIdx).value();
1362 const auto R1 = modelResid[cell_idx][contiMicrobialEqIdx];
1363 R_sum[contiMicrobialEqIdx] += R1;
1364 maxCoeff[contiMicrobialEqIdx] = std::max(maxCoeff[contiMicrobialEqIdx],
1365 std::abs(R1) / pvValue);
1366 B_avg[contiBiofilmEqIdx] += 1.0 / fs.invB(FluidSystem::waterPhaseIdx).value();
1367 const auto R2 = modelResid[cell_idx][contiBiofilmEqIdx];
1368 R_sum[contiBiofilmEqIdx] += R2;
1369 maxCoeff[contiBiofilmEqIdx] = std::max(maxCoeff[contiBiofilmEqIdx],
1370 std::abs(R2) / pvValue);
1371 if constexpr (has_micp_) {
1372 B_avg[contiOxygenEqIdx] += 1.0 / fs.invB(FluidSystem::waterPhaseIdx).value();
1373 const auto R3 = modelResid[cell_idx][contiOxygenEqIdx];
1374 R_sum[contiOxygenEqIdx] += R3;
1375 maxCoeff[contiOxygenEqIdx] = std::max(maxCoeff[contiOxygenEqIdx],
1376 std::abs(R3) / pvValue);
1377 B_avg[contiUreaEqIdx] += 1.0 / fs.invB(FluidSystem::waterPhaseIdx).value();
1378 const auto R4 = modelResid[cell_idx][contiUreaEqIdx];
1379 R_sum[contiUreaEqIdx] += R4;
1380 maxCoeff[contiUreaEqIdx] = std::max(maxCoeff[contiUreaEqIdx],
1381 std::abs(R4) / pvValue);
1382 B_avg[contiCalciteEqIdx] += 1.0 / fs.invB(FluidSystem::waterPhaseIdx).value();
1383 const auto R5 = modelResid[cell_idx][contiCalciteEqIdx];
1384 R_sum[contiCalciteEqIdx] += R5;
1385 maxCoeff[contiCalciteEqIdx] = std::max(maxCoeff[contiCalciteEqIdx],
1386 std::abs(R5) / pvValue);
1387 }
1388 }
1389}
1390
1391} // namespace Opm
1392
1393#endif // OPM_BLACKOILMODEL_IMPL_HEADER_INCLUDED
#define OPM_END_PARALLEL_TRY_CATCH(prefix, comm)
Catch exception and throw in a parallel try-catch clause.
Definition: DeferredLoggingErrorHelpers.hpp:192
#define OPM_BEGIN_PARALLEL_TRY_CATCH()
Macro to setup the try of a parallel try-catch.
Definition: DeferredLoggingErrorHelpers.hpp:158
BlackoilModel(Simulator &simulator, const ModelParameters &param, BlackoilWellModel< TypeTag > &well_model, const bool terminal_output)
Definition: BlackoilModel_impl.hpp:78
std::vector< std::vector< Scalar > > computeFluidInPlace(const T &, const std::vector< int > &fipnum) const
Wrapper required due to not following generic API.
Definition: BlackoilModel.hpp:264
ModelParameters param_
Definition: BlackoilModel.hpp:354
void writeNonlinearIterationsPerCell(const std::filesystem::path &odir) const
Write the number of nonlinear iterations per cell to a file in ResInsight compatible format.
Definition: BlackoilModel_impl.hpp:1235
GetPropType< TypeTag, Properties::IntensiveQuantities > IntensiveQuantities
Definition: BlackoilModel.hpp:67
SimulatorReportSingle nonlinearIterationNewton(const SimulatorTimerInterface &timer, NonlinearSolverType &nonlinear_solver)
Definition: BlackoilModel_impl.hpp:281
const std::vector< SimulatorReport > & domainAccumulatedReports() const
return the statistics of local solves accumulated for each domain on this rank
Definition: BlackoilModel_impl.hpp:1225
std::pair< Scalar, Scalar > localConvergenceData(std::vector< Scalar > &R_sum, std::vector< Scalar > &maxCoeff, std::vector< Scalar > &B_avg, std::vector< int > &maxCoeffCell)
Get reservoir quantities on this process needed for convergence calculations.
Definition: BlackoilModel_impl.hpp:721
GetPropType< TypeTag, Properties::ElementContext > ElementContext
Definition: BlackoilModel.hpp:66
std::vector< StepReport > convergence_reports_
Definition: BlackoilModel.hpp:371
const SimulatorReport & localAccumulatedReports() const
return the statistics of local solves accumulated for this rank
Definition: BlackoilModel_impl.hpp:1214
std::tuple< Scalar, Scalar > convergenceReduction(Parallel::Communication comm, const Scalar pvSumLocal, const Scalar numAquiferPvSumLocal, std::vector< Scalar > &R_sum, std::vector< Scalar > &maxCoeff, std::vector< Scalar > &B_avg)
Definition: BlackoilModel_impl.hpp:661
const Grid & grid_
Definition: BlackoilModel.hpp:343
GetPropType< TypeTag, Properties::SolutionVector > SolutionVector
Definition: BlackoilModel.hpp:69
void updateSolution(const BVector &dx)
Apply an update to the primary variables.
Definition: BlackoilModel_impl.hpp:534
long int global_nc_
The number of cells of the global grid.
Definition: BlackoilModel.hpp:363
DebugFlags
Definition: BlackoilModel.hpp:127
void updateTUNING(const Tuning &tuning)
Definition: BlackoilModel_impl.hpp:853
void getMaxCoeff(const unsigned cell_idx, const IntensiveQuantities &intQuants, const FluidState &fs, const Residual &modelResid, const Scalar pvValue, std::vector< Scalar > &B_avg, std::vector< Scalar > &R_sum, std::vector< Scalar > &maxCoeff, std::vector< int > &maxCoeffCell)
Definition: BlackoilModel_impl.hpp:1272
std::unique_ptr< BlackoilModelNldd< TypeTag > > nlddSolver_
Non-linear DD solver.
Definition: BlackoilModel.hpp:374
void writePartitions(const std::filesystem::path &odir) const
Definition: BlackoilModel_impl.hpp:1245
ConvergenceReport getConvergence(const SimulatorTimerInterface &timer, const int maxIter, std::vector< Scalar > &residual_norms)
Definition: BlackoilModel_impl.hpp:1178
GetPropType< TypeTag, Properties::Simulator > Simulator
Definition: BlackoilModel.hpp:64
void updateTUNINGDP(const TuningDp &tuning_dp)
Definition: BlackoilModel_impl.hpp:866
void initialLinearization(SimulatorReportSingle &report, const int minIter, const int maxIter, const SimulatorTimerInterface &timer)
Definition: BlackoilModel_impl.hpp:184
CnvPvSplitData characteriseCnvPvSplit(const std::vector< Scalar > &B_avg, const double dt)
Compute pore-volume/cell count split among "converged", "relaxed converged", "unconverged" cells base...
Definition: BlackoilModel_impl.hpp:772
SimulatorReportSingle assembleReservoir(const SimulatorTimerInterface &timer)
Assemble the residual and Jacobian of the nonlinear system.
Definition: BlackoilModel_impl.hpp:373
void convergencePerCell(const std::vector< Scalar > &B_avg, const double dt, const double tol_cnv, const double tol_cnv_energy)
Compute the number of Newtons required by each cell in order to satisfy the solution change convergen...
Definition: BlackoilModel_impl.hpp:1132
SimulatorReportSingle nonlinearIteration(const SimulatorTimerInterface &timer, NonlinearSolverType &nonlinear_solver)
Definition: BlackoilModel_impl.hpp:246
void storeSolutionUpdate(const BVector &dx)
Definition: BlackoilModel_impl.hpp:590
void solveJacobianSystem(BVector &x)
Definition: BlackoilModel_impl.hpp:472
MaxSolutionUpdateData getMaxSolutionUpdate(const std::vector< unsigned > &ixCells)
Definition: BlackoilModel_impl.hpp:609
SimulatorReportSingle prepareStep(const SimulatorTimerInterface &timer)
Definition: BlackoilModel_impl.hpp:113
ConvergenceReport getReservoirConvergence(const double reportTime, const double dt, const int maxIter, std::vector< Scalar > &B_avg, std::vector< Scalar > &residual_norms)
Definition: BlackoilModel_impl.hpp:878
Dune::BlockVector< VectorBlockType > BVector
Definition: BlackoilModel.hpp:109
void prepareStoringSolutionUpdate()
Get solution update vector as a PrimaryVariable.
Definition: BlackoilModel_impl.hpp:569
GetPropType< TypeTag, Properties::Scalar > Scalar
Definition: BlackoilModel.hpp:75
Scalar relativeChange() const
Definition: BlackoilModel_impl.hpp:385
Class for handling the blackoil well model.
Definition: BlackoilWellModel.hpp:98
Definition: ConvergenceReport.hpp:38
Severity
Definition: ConvergenceReport.hpp:49
void setReservoirFailed(const ReservoirFailure &rf)
Definition: ConvergenceReport.hpp:266
Interface class for SimulatorTimer objects, to be improved.
Definition: SimulatorTimerInterface.hpp:34
virtual int reportStepNum() const
Current report step number. This might differ from currentStepNum in case of sub stepping.
Definition: SimulatorTimerInterface.hpp:109
virtual bool lastStepFailed() const =0
Return true if last time step failed.
virtual double currentStepLength() const =0
virtual double simulationTimeElapsed() const =0
virtual int currentStepNum() const =0
Dune::Communication< MPIComm > Communication
Definition: ParallelCommunication.hpp:30
std::size_t countGlobalCells(const Grid &grid)
Get the number of cells of a global grid.
Definition: countGlobalCells.hpp:80
Definition: blackoilbioeffectsmodules.hh:45
std::string to_string(const ConvergenceReport::ReservoirFailure::Type t)
Definition: BlackoilModel.hpp:114
Definition: BlackoilModel.hpp:119
Scalar dRvMax
Definition: BlackoilModel.hpp:123
Scalar dPMax
Definition: BlackoilModel.hpp:120
Scalar dSMax
Definition: BlackoilModel.hpp:121
Scalar dRsMax
Definition: BlackoilModel.hpp:122
Solver parameters for the BlackoilModel.
Definition: BlackoilModelParameters.hpp:196
std::string nonlinear_solver_
Nonlinear solver type: newton or nldd.
Definition: BlackoilModelParameters.hpp:363
Definition: AquiferGridUtils.hpp:35
Definition: SimulatorReport.hpp:122
A struct for returning timing data from a simulator to its caller.
Definition: SimulatorReport.hpp:34
double linear_solve_time
Definition: SimulatorReport.hpp:43
double assemble_time
Definition: SimulatorReport.hpp:39
bool converged
Definition: SimulatorReport.hpp:55
double pre_post_time
Definition: SimulatorReport.hpp:40
double linear_solve_setup_time
Definition: SimulatorReport.hpp:42
unsigned int total_newton_iterations
Definition: SimulatorReport.hpp:50
double update_time
Definition: SimulatorReport.hpp:45
unsigned int total_linearizations
Definition: SimulatorReport.hpp:49
unsigned int total_linear_iterations
Definition: SimulatorReport.hpp:51