24#ifndef OPM_BLACKOILMODEL_NLDD_HEADER_INCLUDED
25#define OPM_BLACKOILMODEL_NLDD_HEADER_INCLUDED
27#include <dune/common/timer.hh>
28#include <dune/istl/istlexception.hh>
30#include <opm/common/Exceptions.hpp>
32#include <opm/grid/common/SubGridPart.hpp>
61#include <fmt/format.h>
80template<
class TypeTag>
class BlackoilModel;
83template <
class TypeTag>
100 static constexpr int numEq = Indices::numEq;
108 , wellModel_(model.wellModel())
109 , rank_(model_.simulator().vanguard().grid().comm().rank())
112 const auto& [partition_vector_initial, num_domains_initial] = this->partitionCells();
114 int num_domains = num_domains_initial;
115 std::vector<int> partition_vector = partition_vector_initial;
121 bool isolated_cells =
false;
122 for (
auto& domainId : partition_vector) {
124 domainId = num_domains;
125 isolated_cells =
true;
128 if (isolated_cells) {
133 model.
wellModel().setNlddAdapter(&wellModel_);
136 std::vector<int> sizes(num_domains, 0);
137 for (
const auto& p : partition_vector) {
142 using EntitySeed =
typename Grid::template Codim<0>::EntitySeed;
143 std::vector<std::vector<EntitySeed>> seeds(num_domains);
144 std::vector<std::vector<int>> partitions(num_domains);
145 for (
int domain = 0; domain < num_domains; ++domain) {
146 seeds[domain].resize(sizes[domain]);
147 partitions[domain].resize(sizes[domain]);
152 const auto& grid = model_.simulator().vanguard().grid();
154 std::vector<int> count(num_domains, 0);
155 const auto& gridView = grid.leafGridView();
156 const auto beg = gridView.template begin<0, Dune::Interior_Partition>();
157 const auto end = gridView.template end<0, Dune::Interior_Partition>();
159 for (
auto it = beg; it != end; ++it, ++cell) {
160 const int p = partition_vector[cell];
161 seeds[p][count[p]] = it->seed();
162 partitions[p][count[p]] = cell;
165 assert(count == sizes);
168 for (
int index = 0; index < num_domains; ++index) {
169 std::vector<bool> interior(partition_vector.size(),
false);
170 for (
int ix : partitions[index]) {
174 Dune::SubGridPart<Grid> view{grid, std::move(seeds[index])};
177 const bool skip = isolated_cells && (index == num_domains - 1);
178 this->domains_.emplace_back(index,
179 std::move(partitions[index]),
186 const auto numCells = grid.size(0);
187 previousMobilities_.resize(numCells * FluidSystem::numActivePhases(), 0.0);
188 for (
const auto& domain : domains_) {
189 updateMobilities(domain);
193 domain_needs_solving_.resize(num_domains,
true);
196 domain_matrices_.resize(num_domains);
199 for (
int index = 0; index < num_domains; ++index) {
203 const auto& eclState = model_.simulator().vanguard().eclState();
206 loc_param.
init(eclState.getSimulationConfig().useCPR());
208 if (domains_[index].cells.size() < 200) {
212 const bool force_serial =
true;
213 domain_linsolvers_.emplace_back(model_.simulator(), loc_param, force_serial);
214 domain_linsolvers_.back().setDomainIndex(index);
217 assert(
int(domains_.size()) == num_domains);
219 domain_reports_accumulated_.resize(num_domains);
225 local_reports_accumulated_,
226 domain_reports_accumulated_,
228 wellModel_.numLocalWellsEnd());
235 wellModel_.setupDomains(domains_);
239 template <
class NonlinearSolverType>
242 NonlinearSolverType& nonlinear_solver)
247 if (iteration < model_.param().nldd_num_initial_newton_iter_) {
248 report = model_.nonlinearIterationNewton(iteration,
254 model_.initialLinearization(report, iteration, model_.param().newton_min_iter_,
255 model_.param().newton_max_iter_, timer);
262 Dune::Timer localSolveTimer;
263 Dune::Timer detailTimer;
264 localSolveTimer.start();
266 auto& solution = model_.simulator().model().solution(0);
267 auto initial_solution = solution;
268 auto locally_solved = initial_solution;
271 const auto domain_order = this->getSubdomainOrder();
276 std::vector<SimulatorReportSingle> domain_reports(domains_.size());
279 for (
const int domain_index : domain_order) {
280 const auto& domain = domains_[domain_index];
285 domain_needs_solving_[domain_index] = checkIfSubdomainNeedsSolving(domain, iteration);
287 updateMobilities(domain);
289 if (domain.skip || !domain_needs_solving_[domain_index]) {
292 domain_reports[domain.index] = local_report;
295 switch (model_.param().local_solve_approach_) {
297 solveDomainJacobi(solution, locally_solved, local_report, logger,
298 iteration, timer, domain);
302 solveDomainGaussSeidel(solution, locally_solved, local_report, logger,
303 iteration, timer, domain);
311 logger.
debug(fmt::format(
"Convergence failure in domain {} on rank {}." , domain.index, rank_));
314 domain_reports[domain.index] = local_report;
322 global_logger.logMessages();
327 std::array<int, 5> counts{ 0, 0, 0,
static_cast<int>(domain_reports.size()), 0 };
328 int& num_converged = counts[0];
329 int& num_converged_already = counts[1];
330 int& num_local_newtons = counts[2];
331 int& num_domains = counts[3];
332 int& num_skipped = counts[4];
334 auto step_newtons = 0;
335 const auto dr_size = domain_reports.size();
336 for (
auto i = 0*dr_size; i < dr_size; ++i) {
338 domain_needs_solving_[i] =
false;
339 const auto& dr = domain_reports[i];
342 if (dr.total_newton_iterations == 0) {
343 ++num_converged_already;
347 domain_needs_solving_[i] =
true;
350 if (dr.skipped_domains) {
353 step_newtons += dr.total_newton_iterations;
355 domain_reports_accumulated_[i] += dr;
357 local_reports_accumulated_ += dr;
359 num_local_newtons = step_newtons;
363 solution = locally_solved;
364 model_.simulator().model().invalidateAndUpdateIntensiveQuantities(0);
375 const auto& comm = model_.simulator().vanguard().grid().comm();
376 if (comm.size() > 1) {
377 const auto* ccomm = model_.simulator().model().newtonMethod().linearSolver().comm();
380 ccomm->copyOwnerToAll(solution, solution);
383 const std::size_t num = solution.size();
384 Dune::BlockVector<std::size_t> allmeanings(num);
385 for (std::size_t ii = 0; ii < num; ++ii) {
388 ccomm->copyOwnerToAll(allmeanings, allmeanings);
389 for (std::size_t ii = 0; ii < num; ++ii) {
394 model_.simulator().model().invalidateAndUpdateIntensiveQuantitiesOverlap(0);
397 comm.sum(counts.data(), counts.size());
401 const bool is_iorank = this->rank_ == 0;
403 OpmLog::debug(fmt::format(
"Local solves finished. Converged for {}/{} domains. {} domains were skipped. {} domains did no work. {} total local Newton iterations.\n",
404 num_converged, num_domains, num_skipped, num_converged_already, num_local_newtons));
406 auto total_local_solve_time = localSolveTimer.stop();
416 auto rep = model_.nonlinearIterationNewton(iteration + 100, timer, nonlinear_solver);
427 return local_reports_accumulated_;
435 const auto dr_size = domain_reports_accumulated_.size();
437 for (
auto i = 0*dr_size; i < dr_size; ++i) {
438 domain_reports_accumulated_[i].success.num_wells = 0;
441 for (
const auto& [wname, domain] : wellModel_.well_domain()) {
442 domain_reports_accumulated_[domain].success.num_wells++;
444 return domain_reports_accumulated_;
451 const auto& grid = this->model_.simulator().vanguard().grid();
452 const auto& elementMapper = this->model_.simulator().model().elementMapper();
453 const auto& cartMapper = this->model_.simulator().vanguard().cartesianIndexMapper();
466 const auto& grid = this->model_.simulator().vanguard().grid();
467 const auto& elementMapper = this->model_.simulator().model().elementMapper();
468 const auto& cartMapper = this->model_.simulator().vanguard().cartesianIndexMapper();
473 domain_reports_accumulated_,
482 solveDomain(
const Domain& domain,
486 [[maybe_unused]]
const int global_iteration,
487 const bool initial_assembly_required)
489 auto& modelSimulator = model_.simulator();
490 Dune::Timer detailTimer;
492 modelSimulator.model().newtonMethod().setIterationIndex(0);
499 if (initial_assembly_required) {
501 modelSimulator.model().newtonMethod().setIterationIndex(iter);
504 wellModel_.assemble(modelSimulator.model().newtonMethod().numIterations(),
505 modelSimulator.timeStepSize(),
507 const double tt0 = detailTimer.stop();
513 this->assembleReservoirDomain(domain);
519 std::vector<Scalar> resnorms;
520 auto convreport = this->getDomainConvergence(domain, timer, 0, logger, resnorms);
522 if (convreport.converged()) {
532 model_.wellModel().linearizeDomain(domain,
533 modelSimulator.model().linearizer().jacobian(),
534 modelSimulator.model().linearizer().residual());
535 const double tt1 = detailTimer.stop();
540 const int max_iter = model_.param().max_local_solve_iterations_;
541 const auto& grid = modelSimulator.vanguard().grid();
542 double damping_factor = 1.0;
543 std::vector<std::vector<Scalar>> convergence_history;
544 convergence_history.reserve(20);
545 convergence_history.push_back(resnorms);
549 const int nc = grid.size(0);
553 double setup_time = 0.0;
555 this->solveJacobianSystemDomain(domain, x, setup_time);
557 catch (
const NumericalProblem& e) {
559 logger.
debug(fmt::format(
560 "Local linear solver failed in domain {} on rank {}: {}",
561 domain.index, rank_, e.what()));
566 modelSimulator.problem().endIteration();
572 model_.wellModel().postSolveDomain(x, domain);
573 if (damping_factor != 1.0) {
583 this->updateDomainSolution(domain, x);
590 modelSimulator.model().newtonMethod().setIterationIndex(iter);
594 wellModel_.assemble(modelSimulator.model().newtonMethod().numIterations(),
595 modelSimulator.timeStepSize(),
597 const double tt3 = detailTimer.stop();
602 this->assembleReservoirDomain(domain);
609 convreport = this->getDomainConvergence(domain, timer, iter, logger, resnorms);
610 convergence_history.push_back(resnorms);
617 model_.wellModel().linearizeDomain(domain,
618 modelSimulator.model().linearizer().jacobian(),
619 modelSimulator.model().linearizer().residual());
620 const double tt2 = detailTimer.stop();
625 if (!convreport.converged() && !convreport.wellFailed()) {
626 bool oscillate =
false;
627 bool stagnate =
false;
628 const auto num_residuals = convergence_history.front().size();
630 Scalar{0.2}, 1, oscillate, stagnate);
632 damping_factor *= 0.85;
633 logger.
debug(fmt::format(
"| Damping factor is now {}", damping_factor));
636 }
while (!convreport.converged() && iter <= max_iter);
638 modelSimulator.problem().endIteration();
640 local_report.
converged = convreport.converged();
648 void assembleReservoirDomain(
const Domain& domain)
650 OPM_TIMEBLOCK(assembleReservoirDomain);
652 model_.simulator().model().linearizer().linearizeDomain(domain,
true);
656 void solveJacobianSystemDomain(
const Domain& domain,
BVector& global_x,
double& setup_time)
658 const auto& modelSimulator = model_.simulator();
660 Dune::Timer perfTimer;
663 const Mat& main_matrix = modelSimulator.model().linearizer().jacobian().istlMatrix();
664 if (domain_matrices_[domain.index]) {
667 domain_matrices_[domain.index] = std::make_unique<Mat>(
Details::extractMatrix(main_matrix, domain.cells));
669 auto& jac = *domain_matrices_[domain.index];
678 auto& linsolver = domain_linsolvers_[domain.index];
680 linsolver.prepare(jac, res);
681 setup_time = perfTimer.stop();
682 linsolver.setResidual(res);
689 void updateDomainSolution(
const Domain& domain,
const BVector& dx)
691 OPM_TIMEBLOCK(updateDomainSolution);
692 auto& simulator = model_.simulator();
693 auto& newtonMethod = simulator.model().newtonMethod();
696 newtonMethod.update_(solution,
705 simulator.model().invalidateAndUpdateIntensiveQuantities(0, domain);
709 std::pair<Scalar, Scalar> localDomainConvergenceData(
const Domain& domain,
710 std::vector<Scalar>& R_sum,
711 std::vector<Scalar>& maxCoeff,
712 std::vector<Scalar>& B_avg,
713 std::vector<int>& maxCoeffCell)
715 const auto& modelSimulator = model_.simulator();
718 Scalar numAquiferPvSumLocal = 0.0;
719 const auto& model = modelSimulator.model();
720 const auto& problem = modelSimulator.problem();
722 const auto& modelResid = modelSimulator.model().linearizer().residual();
725 const auto& gridView = domain.view;
726 const auto& elemEndIt = gridView.template end<0>();
727 IsNumericalAquiferCell isNumericalAquiferCell(gridView.grid());
729 for (
auto elemIt = gridView.template begin</*codim=*/0>();
733 if (elemIt->partitionType() != Dune::InteriorEntity) {
736 const auto& elem = *elemIt;
737 elemCtx.updatePrimaryStencil(elem);
738 elemCtx.updatePrimaryIntensiveQuantities(0);
740 const unsigned cell_idx = elemCtx.globalSpaceIndex(0, 0);
741 const auto& intQuants = elemCtx.intensiveQuantities(0, 0);
742 const auto& fs = intQuants.fluidState();
744 const auto pvValue = problem.referencePorosity(cell_idx, 0) *
745 model.dofTotalVolume(cell_idx);
746 pvSumLocal += pvValue;
748 if (isNumericalAquiferCell(elem))
750 numAquiferPvSumLocal += pvValue;
753 model_.getMaxCoeff(cell_idx, intQuants, fs, modelResid, pvValue,
754 B_avg, R_sum, maxCoeff, maxCoeffCell);
758 const int bSize = B_avg.size();
759 for (
int i = 0; i<bSize; ++i )
761 B_avg[ i ] /=
Scalar(domain.cells.size());
764 return {pvSumLocal, numAquiferPvSumLocal};
767 ConvergenceReport getDomainReservoirConvergence(
const double reportTime,
771 DeferredLogger& logger,
772 std::vector<Scalar>& B_avg,
773 std::vector<Scalar>& residual_norms)
775 using Vector = std::vector<Scalar>;
777 const int numComp =
numEq;
778 Vector R_sum(numComp, 0.0 );
779 Vector maxCoeff(numComp, std::numeric_limits<Scalar>::lowest() );
780 std::vector<int> maxCoeffCell(numComp, -1);
781 const auto [ pvSum, numAquiferPvSum]
782 = this->localDomainConvergenceData(domain, R_sum, maxCoeff, B_avg, maxCoeffCell);
784 auto cnvErrorPvFraction = computeCnvErrorPvLocal(domain, B_avg, dt);
785 cnvErrorPvFraction /= (pvSum - numAquiferPvSum);
790 const bool use_relaxed_cnv = cnvErrorPvFraction < model_.param().relaxed_max_pv_fraction_ &&
791 iteration >= model_.param().min_strict_cnv_iter_;
794 const Scalar tol_cnv = model_.param().local_tolerance_scaling_cnv_
795 * (use_relaxed_cnv ? model_.param().tolerance_cnv_relaxed_
796 : model_.param().tolerance_cnv_);
798 const bool use_relaxed_mb = iteration >= model_.param().min_strict_mb_iter_;
799 const Scalar tol_mb = model_.param().local_tolerance_scaling_mb_
800 * (use_relaxed_mb ? model_.param().tolerance_mb_relaxed_ : model_.param().tolerance_mb_);
803 std::vector<Scalar> CNV(numComp);
804 std::vector<Scalar> mass_balance_residual(numComp);
805 for (
int compIdx = 0; compIdx < numComp; ++compIdx )
807 CNV[compIdx] = B_avg[compIdx] * dt * maxCoeff[compIdx];
808 mass_balance_residual[compIdx] = std::abs(B_avg[compIdx]*R_sum[compIdx]) * dt / pvSum;
809 residual_norms.push_back(CNV[compIdx]);
813 ConvergenceReport report{reportTime};
814 using CR = ConvergenceReport;
815 for (
int compIdx = 0; compIdx < numComp; ++compIdx) {
816 Scalar res[2] = { mass_balance_residual[compIdx], CNV[compIdx] };
817 CR::ReservoirFailure::Type types[2] = { CR::ReservoirFailure::Type::MassBalance,
818 CR::ReservoirFailure::Type::Cnv };
819 Scalar tol[2] = { tol_mb, tol_cnv };
820 for (
int ii : {0, 1}) {
821 if (std::isnan(res[ii])) {
822 report.setReservoirFailed({types[ii], CR::Severity::NotANumber, compIdx});
823 logger.debug(
"NaN residual for " + model_.compNames().name(compIdx) +
" equation.");
824 }
else if (res[ii] > model_.param().max_residual_allowed_) {
825 report.setReservoirFailed({types[ii], CR::Severity::TooLarge, compIdx});
826 logger.debug(
"Too large residual for " + model_.compNames().name(compIdx) +
" equation.");
827 }
else if (res[ii] < 0.0) {
828 report.setReservoirFailed({types[ii], CR::Severity::Normal, compIdx});
829 logger.debug(
"Negative residual for " + model_.compNames().name(compIdx) +
" equation.");
830 }
else if (res[ii] > tol[ii]) {
831 report.setReservoirFailed({types[ii], CR::Severity::Normal, compIdx});
834 report.setReservoirConvergenceMetric(types[ii], compIdx, res[ii], tol[ii]);
839 const bool converged_at_initial_state = (report.converged() && iteration == 0);
840 if (!converged_at_initial_state) {
841 if (iteration == 0) {
843 std::string msg = fmt::format(
"Domain {} on rank {}, size {}, containing cell {}\n| Iter",
844 domain.index, this->rank_, domain.cells.size(), domain.cells[0]);
845 for (
int compIdx = 0; compIdx < numComp; ++compIdx) {
847 msg += model_.compNames().name(compIdx)[0];
850 for (
int compIdx = 0; compIdx < numComp; ++compIdx) {
852 msg += model_.compNames().name(compIdx)[0];
858 std::ostringstream ss;
860 const std::streamsize oprec = ss.precision(3);
861 const std::ios::fmtflags oflags = ss.setf(std::ios::scientific);
862 ss << std::setw(4) << iteration;
863 for (
int compIdx = 0; compIdx < numComp; ++compIdx) {
864 ss << std::setw(11) << mass_balance_residual[compIdx];
866 for (
int compIdx = 0; compIdx < numComp; ++compIdx) {
867 ss << std::setw(11) << CNV[compIdx];
871 logger.debug(ss.str());
877 ConvergenceReport getDomainConvergence(
const Domain& domain,
878 const SimulatorTimerInterface& timer,
880 DeferredLogger& logger,
881 std::vector<Scalar>& residual_norms)
883 OPM_TIMEBLOCK(getDomainConvergence);
884 std::vector<Scalar> B_avg(
numEq, 0.0);
885 auto report = this->getDomainReservoirConvergence(timer.simulationTimeElapsed(),
886 timer.currentStepLength(),
892 report += wellModel_.getWellConvergence(domain, B_avg, logger);
897 std::vector<int> getSubdomainOrder()
899 const auto& modelSimulator = model_.simulator();
900 const auto& solution = modelSimulator.model().solution(0);
902 std::vector<int> domain_order(domains_.size());
903 std::iota(domain_order.begin(), domain_order.end(), 0);
910 std::vector<Scalar> measure_per_domain(domains_.size());
911 switch (model_.param().local_domains_ordering_) {
914 for (
const auto& domain : domains_) {
916 std::accumulate(domain.cells.begin(), domain.cells.end(),
Scalar{0},
917 [&solution](
const auto acc,
const auto c)
918 { return acc + solution[c][Indices::pressureSwitchIdx]; });
919 const Scalar avgpress = press_sum / domain.cells.size();
920 measure_per_domain[domain.index] = avgpress;
926 for (
const auto& domain : domains_) {
927 measure_per_domain[domain.index] =
928 std::accumulate(domain.cells.begin(), domain.cells.end(),
Scalar{0},
929 [&solution](
const auto acc,
const auto c)
930 { return std::max(acc, solution[c][Indices::pressureSwitchIdx]); });
936 const auto& residual = modelSimulator.model().linearizer().residual();
937 const int num_vars = residual[0].size();
938 for (
const auto& domain : domains_) {
940 for (
const int c : domain.cells) {
941 for (
int ii = 0; ii < num_vars; ++ii) {
942 maxres = std::max(maxres, std::fabs(residual[c][ii]));
945 measure_per_domain[domain.index] = maxres;
952 const auto& m = measure_per_domain;
953 std::stable_sort(domain_order.begin(), domain_order.end(),
954 [&m](
const int i1,
const int i2){ return m[i1] > m[i2]; });
957 throw std::logic_error(
"Domain solve approach must be Jacobi or Gauss-Seidel");
961 template<
class GlobalEqVector>
962 void solveDomainJacobi(GlobalEqVector& solution,
963 GlobalEqVector& locally_solved,
964 SimulatorReportSingle& local_report,
965 DeferredLogger& logger,
967 const SimulatorTimerInterface& timer,
970 auto initial_local_well_primary_vars = wellModel_.getPrimaryVarsDomain(domain.index);
972 auto convrep = solveDomain(domain, timer, local_report, logger, iteration,
false);
973 if (local_report.converged) {
977 model_.simulator().model().invalidateAndUpdateIntensiveQuantities(0, domain);
979 wellModel_.setPrimaryVarsDomain(domain.index, initial_local_well_primary_vars);
981 model_.simulator().model().invalidateAndUpdateIntensiveQuantities(0, domain);
985 template<
class GlobalEqVector>
986 void solveDomainGaussSeidel(GlobalEqVector& solution,
987 GlobalEqVector& locally_solved,
988 SimulatorReportSingle& local_report,
989 DeferredLogger& logger,
991 const SimulatorTimerInterface& timer,
994 auto initial_local_well_primary_vars = wellModel_.getPrimaryVarsDomain(domain.index);
996 auto convrep = solveDomain(domain, timer, local_report, logger, iteration,
true);
997 if (!local_report.converged) {
1001 if (!convrep.wellFailed()) {
1005 for (
const auto& rc : convrep.reservoirConvergence()) {
1007 mb_sum += rc.value();
1009 cnv_sum += rc.value();
1013 const Scalar acceptable_local_mb_sum = 1e-3;
1014 const Scalar acceptable_local_cnv_sum = 1.0;
1015 if (mb_sum < acceptable_local_mb_sum && cnv_sum < acceptable_local_cnv_sum) {
1016 local_report.converged =
true;
1017 local_report.accepted_unconverged_domains += 1;
1018 logger.debug(fmt::format(
"Accepting solution in unconverged domain {} on rank {}.", domain.index, rank_));
1019 logger.debug(fmt::format(
"Value of mb_sum: {} cnv_sum: {}", mb_sum, cnv_sum));
1021 logger.debug(
"Unconverged local solution.");
1024 logger.debug(
"Unconverged local solution with well convergence failures:");
1025 for (
const auto& wf : convrep.wellFailures()) {
1030 if (local_report.converged) {
1031 local_report.converged_domains += 1;
1035 local_report.unconverged_domains += 1;
1036 wellModel_.setPrimaryVarsDomain(domain.index, initial_local_well_primary_vars);
1038 model_.simulator().model().invalidateAndUpdateIntensiveQuantities(0, domain);
1043 const std::vector<Scalar>& B_avg,
double dt)
const
1046 const auto& simulator = model_.simulator();
1047 const auto& model = simulator.model();
1048 const auto& problem = simulator.problem();
1049 const auto& residual = simulator.model().linearizer().residual();
1051 for (
const int cell_idx : domain.cells) {
1052 const Scalar pvValue = problem.referencePorosity(cell_idx, 0) *
1053 model.dofTotalVolume(cell_idx);
1054 const auto& cellResidual = residual[cell_idx];
1055 bool cnvViolated =
false;
1057 for (
unsigned eqIdx = 0; eqIdx < cellResidual.size(); ++eqIdx) {
1059 Scalar CNV = cellResidual[eqIdx] * dt * B_avg[eqIdx] / pvValue;
1060 cnvViolated = cnvViolated || (fabs(CNV) > model_.param().tolerance_cnv_);
1070 decltype(
auto) partitionCells()
const
1072 const auto& grid = this->model_.simulator().vanguard().grid();
1074 using GridView = std::remove_cv_t<std::remove_reference_t<
1075 decltype(grid.leafGridView())>>;
1077 using Element = std::remove_cv_t<std::remove_reference_t<
1078 typename GridView::template Codim<0>::Entity>>;
1080 const auto& param = this->model_.param();
1082 auto zoltan_ctrl = ZoltanPartitioningControl<Element>{};
1084 zoltan_ctrl.domain_imbalance = param.local_domains_partition_imbalance_;
1087 [elementMapper = &this->model_.simulator().model().elementMapper()]
1088 (
const Element& element)
1090 return elementMapper->index(element);
1093 zoltan_ctrl.local_to_global =
1094 [cartMapper = &this->model_.simulator().vanguard().cartesianIndexMapper()]
1097 return cartMapper->cartesianIndex(elemIdx);
1101 const auto need_wells = param.local_domains_partition_method_ ==
"zoltan";
1103 const auto wells = need_wells
1104 ? this->model_.simulator().vanguard().schedule().getWellsatEnd()
1105 : std::vector<Well>{};
1107 const auto& possibleFutureConnectionSet = need_wells
1108 ? this->model_.simulator().vanguard().schedule().getPossibleFutureConnections()
1109 : std::unordered_map<std::string, std::set<int>> {};
1112 constexpr int default_cells_per_domain = 1000;
1113 const int num_domains = (param.num_local_domains_ > 0)
1114 ? param.num_local_domains_
1117 num_domains, grid.leafGridView(), wells,
1118 possibleFutureConnectionSet, zoltan_ctrl,
1119 param.local_domains_partition_well_neighbor_levels_);
1122 void updateMobilities(
const Domain& domain)
1124 if (domain.skip || model_.param().nldd_relative_mobility_change_tol_ == 0.0) {
1127 const auto numActivePhases = FluidSystem::numActivePhases();
1128 for (
const auto globalDofIdx : domain.cells) {
1129 const auto& intQuants = model_.simulator().model().intensiveQuantities(globalDofIdx, 0);
1131 for (
unsigned activePhaseIdx = 0; activePhaseIdx < numActivePhases; ++activePhaseIdx) {
1132 const auto phaseIdx = FluidSystem::activeToCanonicalPhaseIdx(activePhaseIdx);
1133 const auto mobIdx = globalDofIdx * numActivePhases + activePhaseIdx;
1134 previousMobilities_[mobIdx] = getValue(intQuants.mobility(phaseIdx));
1139 bool checkIfSubdomainNeedsSolving(
const Domain& domain,
const int iteration)
1147 if (domain_needs_solving_[domain.index]) {
1152 if (model_.param().nldd_relative_mobility_change_tol_ == 0.0) {
1157 if (iteration == 0) {
1161 return checkSubdomainChangeRelative(domain);
1164 bool checkSubdomainChangeRelative(
const Domain& domain)
1166 const auto numActivePhases = FluidSystem::numActivePhases();
1169 for (
const auto globalDofIdx : domain.cells) {
1170 const auto& intQuants = model_.simulator().model().intensiveQuantities(globalDofIdx, 0);
1174 for (
unsigned activePhaseIdx = 0; activePhaseIdx < numActivePhases; ++activePhaseIdx) {
1175 const auto mobIdx = globalDofIdx * numActivePhases + activePhaseIdx;
1176 cellMob += previousMobilities_[mobIdx] / numActivePhases;
1180 for (
unsigned activePhaseIdx = 0; activePhaseIdx < numActivePhases; ++activePhaseIdx) {
1181 const auto phaseIdx = FluidSystem::activeToCanonicalPhaseIdx(activePhaseIdx);
1182 const auto mobIdx = globalDofIdx * numActivePhases + activePhaseIdx;
1183 const auto mobility = getValue(intQuants.mobility(phaseIdx));
1184 const auto relDiff = std::abs(mobility - previousMobilities_[mobIdx]) / cellMob;
1185 if (relDiff > model_.param().nldd_relative_mobility_change_tol_) {
1193 BlackoilModel<TypeTag>& model_;
1194 BlackoilWellModelNldd<TypeTag> wellModel_;
1195 std::vector<Domain> domains_;
1196 std::vector<std::unique_ptr<Mat>> domain_matrices_;
1197 std::vector<ISTLSolverType> domain_linsolvers_;
1198 SimulatorReport local_reports_accumulated_;
1200 mutable std::vector<SimulatorReport> domain_reports_accumulated_;
1203 std::vector<Scalar> previousMobilities_;
1205 std::vector<bool> domain_needs_solving_;
#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
A NLDD implementation for three-phase black oil.
Definition: BlackoilModelNldd.hpp:85
GetPropType< TypeTag, Properties::SolutionVector > SolutionVector
Definition: BlackoilModelNldd.hpp:93
GetPropType< TypeTag, Properties::Indices > Indices
Definition: BlackoilModelNldd.hpp:90
const std::vector< SimulatorReport > & domainAccumulatedReports() const
return the statistics of local solves accumulated for each domain on this rank
Definition: BlackoilModelNldd.hpp:431
void writePartitions(const std::filesystem::path &odir) const
Definition: BlackoilModelNldd.hpp:449
BlackoilModelNldd(BlackoilModel< TypeTag > &model)
The constructor sets up the subdomains.
Definition: BlackoilModelNldd.hpp:106
GetPropType< TypeTag, Properties::Scalar > Scalar
Definition: BlackoilModelNldd.hpp:91
typename BlackoilModel< TypeTag >::Mat Mat
Definition: BlackoilModelNldd.hpp:98
const SimulatorReport & localAccumulatedReports() const
return the statistics of local solves accumulated for this rank
Definition: BlackoilModelNldd.hpp:425
SimulatorReportSingle nonlinearIterationNldd(const int iteration, const SimulatorTimerInterface &timer, NonlinearSolverType &nonlinear_solver)
Do one non-linear NLDD iteration.
Definition: BlackoilModelNldd.hpp:240
SubDomain< Grid > Domain
Definition: BlackoilModelNldd.hpp:96
GetPropType< TypeTag, Properties::ElementContext > ElementContext
Definition: BlackoilModelNldd.hpp:87
GetPropType< TypeTag, Properties::FluidSystem > FluidSystem
Definition: BlackoilModelNldd.hpp:88
void prepareStep()
Called before starting a time step.
Definition: BlackoilModelNldd.hpp:232
void writeNonlinearIterationsPerCell(const std::filesystem::path &odir) const
Write the number of nonlinear iterations per cell to a file in ResInsight compatible format.
Definition: BlackoilModelNldd.hpp:464
static constexpr int numEq
Definition: BlackoilModelNldd.hpp:100
GetPropType< TypeTag, Properties::Grid > Grid
Definition: BlackoilModelNldd.hpp:89
typename BlackoilModel< TypeTag >::BVector BVector
Definition: BlackoilModelNldd.hpp:95
Definition: BlackoilModel.hpp:61
BlackoilWellModel< TypeTag > & wellModel()
return the StandardWells object
Definition: BlackoilModel.hpp:305
typename SparseMatrixAdapter::IstlMatrix Mat
Definition: BlackoilModel.hpp:108
Dune::BlockVector< VectorBlockType > BVector
Definition: BlackoilModel.hpp:109
Definition: ConvergenceReport.hpp:38
Definition: DeferredLogger.hpp:57
void debug(const std::string &tag, const std::string &message)
Definition: ISTLSolver.hpp:149
Interface class for SimulatorTimer objects, to be improved.
Definition: SimulatorTimerInterface.hpp:34
Vector extractVector(const Vector &x, const std::vector< int > &indices)
Definition: extractMatrix.hpp:104
void copySubMatrix(const Matrix &A, const std::vector< int > &indices, Matrix &B)
Definition: extractMatrix.hpp:35
void setGlobal(const Vector &x, const std::vector< int > &indices, Vector &global_x)
Definition: extractMatrix.hpp:115
Matrix extractMatrix(const Matrix &m, const std::vector< int > &indices)
Definition: extractMatrix.hpp:47
void unPack(PV &privar, const std::size_t meanings)
Definition: priVarsPacking.hpp:41
std::size_t pack(const PV &privar)
Definition: priVarsPacking.hpp:31
std::size_t countGlobalCells(const Grid &grid)
Get the number of cells of a global grid.
Definition: countGlobalCells.hpp:80
void detectOscillations(const std::vector< std::vector< Scalar > > &residualHistory, const int it, const int numPhases, const Scalar relaxRelTol, const int minimumOscillatingPhases, bool &oscillate, bool &stagnate)
Detect oscillation or stagnation in a given residual history.
Definition: blackoilbioeffectsmodules.hh:43
Opm::DeferredLogger gatherDeferredLogger(const Opm::DeferredLogger &local_deferredlogger, Parallel::Communication communicator)
Create a global log combining local logs.
std::pair< std::vector< int >, int > partitionCells(const std::string &method, const int num_local_domains, const GridView &grid_view, const std::vector< Well > &wells, const std::unordered_map< std::string, std::set< int > > &possibleFutureConnections, const ZoltanPartitioningControl< Element > &zoltan_ctrl, const int num_neighbor_levels)
void printDomainDistributionSummary(const std::vector< int > &partition_vector, const std::vector< Domain > &domains, SimulatorReport &local_reports_accumulated, std::vector< SimulatorReport > &domain_reports_accumulated, const Grid &grid, int num_wells)
Definition: NlddReporting.hpp:219
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property (equivalent to old macro GET_PROP_TYPE(....
Definition: propertysystem.hh:233
void writePartitions(const std::filesystem::path &odir, const std::vector< Domain > &domains, const Grid &grid, const ElementMapper &elementMapper, const CartMapper &cartMapper)
Definition: NlddReporting.hpp:164
std::string to_string(const ConvergenceReport::ReservoirFailure::Type t)
void writeNonlinearIterationsPerCell(const std::filesystem::path &odir, const std::vector< Domain > &domains, const std::vector< SimulatorReport > &domain_reports, const Grid &grid, const ElementMapper &elementMapper, const CartMapper &cartMapper)
Definition: NlddReporting.hpp:104
Solver parameters for the BlackoilModel.
Definition: BlackoilModelParameters.hpp:194
This class carries all parameters for the NewtonIterationBlackoilInterleaved class.
Definition: FlowLinearSolverParameters.hpp:98
bool is_nldd_local_solver_
Definition: FlowLinearSolverParameters.hpp:112
void init(bool cprRequestedInDataFile)
bool linear_solver_print_json_definition_
Definition: FlowLinearSolverParameters.hpp:114
std::string linsolver_
Definition: FlowLinearSolverParameters.hpp:113
Definition: SimulatorReport.hpp:122
SimulatorReportSingle success
Definition: SimulatorReport.hpp:123
A struct for returning timing data from a simulator to its caller.
Definition: SimulatorReport.hpp:34
double linear_solve_time
Definition: SimulatorReport.hpp:43
int skipped_domains
Definition: SimulatorReport.hpp:71
double assemble_time
Definition: SimulatorReport.hpp:39
double assemble_time_well
Definition: SimulatorReport.hpp:41
double solver_time
Definition: SimulatorReport.hpp:38
bool converged
Definition: SimulatorReport.hpp:55
double pre_post_time
Definition: SimulatorReport.hpp:40
double total_time
Definition: SimulatorReport.hpp:37
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
double local_solve_time
Definition: SimulatorReport.hpp:44
unsigned int total_linearizations
Definition: SimulatorReport.hpp:49
unsigned int total_linear_iterations
Definition: SimulatorReport.hpp:51
Definition: SubDomain.hpp:85