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