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