38#ifndef OPENRS_INCOMPFLOWSOLVERHYBRID_HEADER
39#define OPENRS_INCOMPFLOWSOLVERHYBRID_HEADER
41#include <opm/common/ErrorMacros.hpp>
42#include <opm/grid/utility/SparseTable.hpp>
46#include <opm/common/utility/platform_dependent/disable_warnings.h>
48#include <unordered_map>
50#include <dune/common/fvector.hh>
51#include <dune/common/fmatrix.hh>
53#include <dune/istl/bvector.hh>
54#include <dune/istl/bcrsmatrix.hh>
55#include <dune/istl/operators.hh>
56#include <dune/istl/io.hh>
57#include <dune/istl/overlappingschwarz.hh>
58#include <dune/istl/schwarz.hh>
59#include <dune/istl/preconditioners.hh>
60#include <dune/istl/solvers.hh>
61#include <dune/istl/owneroverlapcopy.hh>
62#include <dune/istl/paamg/amg.hh>
63#include <dune/common/version.hh>
64#include <dune/istl/paamg/fastamg.hh>
65#include <dune/istl/paamg/kamg.hh>
66#include <dune/istl/paamg/pinfo.hh>
68#include <opm/common/utility/platform_dependent/reenable_warnings.h>
108 bool topologyIsSane(
const GI& g)
110 typedef typename GI::CellIterator CI;
111 typedef typename CI::FaceIterator FI;
113 bool sane = g.numberOfCells() >= 0;
115 for (CI c = g.cellbegin(); sane && c != g.cellend(); ++c) {
118 for (FI f = c->facebegin(); f != c->faceend(); ++f) {
119 if (!f->boundary()) {
120 n.push_back(f->neighbourCellIndex());
123 std::ranges::sort(
n);
125 sane = std::unique(
n.begin(),
n.end()) ==
n.end();
165 axpby(
const T&
a,
const T& b) : a_(
a), b_(b) {}
179 T operator()(
const T&
x,
const T&
y)
375 template <
class GridInterface,
378 template <
class Gr
idIF,
class RockIF>
class InnerProduct>
385 typedef typename GridInterface::Scalar Scalar;
391 enum FaceType { Internal, Dirichlet, Neumann, Periodic };
402 typedef typename GridInterface::Scalar Scalar;
407 typedef typename GridInterface::CellIterator CI;
411 typedef typename CI ::FaceIterator FI;
424 Scalar pressure(
const CI& c)
const
426 return pressure_[cellno_[c->index()]];
439 Scalar outflux (
const FI& f)
const
441 return outflux_[cellno_[f->cellIndex()]][f->localIndex()];
443 Scalar outflux (
int hf)
const
445 return outflux_.data(hf);
448 std::vector< int > cellno_;
449 Opm::SparseTable< int > cellFaces_;
450 std::vector<Scalar> pressure_;
451 Opm::SparseTable<Scalar> outflux_;
454 std::vector<int>().swap(cellno_);
457 std::vector<Scalar>().swap(pressure_);
489 template<
class Po
int>
490 void init(
const GridInterface& g,
491 const RockInterface& r,
493 const BCInterface& bc)
497 if (g.numberOfCells() > 0) {
515 num_internal_faces_ = 0;
516 total_num_faces_ = 0;
517 matrix_structure_valid_ =
false;
518 do_regularization_ =
true;
520 bdry_id_map_.clear();
522 std::vector<Scalar>().swap(L_);
523 std::vector<Scalar>().swap(g_);
526 flowSolution_.clear();
528 cleared_state_ =
true;
548 const BCInterface& bc)
551 assert (cleared_state_);
553 assert (topologyIsSane(g));
556 allocateConnections(bc);
578 template<
class Po
int>
583 assert (matrix_structure_valid_);
585 typedef typename GridInterface::CellIterator CI;
586 const Opm::SparseTable<int>& cf = flowSolution_.cellFaces_;
589 for (CI c = pgrid_->cellbegin(); c != pgrid_->cellend(); ++c, ++i) {
590 ip_.buildStaticContrib(c, r, grav, cf.rowSize(i));
671 template<
class Flu
idInterface>
672 void solve(
const FluidInterface& r ,
673 const std::vector<double>& sat,
674 const BCInterface& bc ,
675 const std::vector<double>& src,
676 double residual_tolerance = 1e-8,
677 int linsolver_verbosity = 1,
678 int linsolver_type = 1,
679 bool same_matrix =
false,
680 int linsolver_maxit = 0,
681 double prolongate_factor = 1.6,
682 int smooth_steps = 1)
684 assembleDynamic(r, sat, bc, src);
688 switch (linsolver_type) {
690 solveLinearSystem(residual_tolerance, linsolver_verbosity, linsolver_maxit);
693 solveLinearSystemAMG(residual_tolerance, linsolver_verbosity,
694 linsolver_maxit, prolongate_factor, same_matrix, smooth_steps);
698 solveLinearSystemKAMG(residual_tolerance, linsolver_verbosity,
699 linsolver_maxit, prolongate_factor, same_matrix,smooth_steps);
702 solveLinearSystemFastAMG(residual_tolerance, linsolver_verbosity,
703 linsolver_maxit, prolongate_factor, same_matrix,smooth_steps);
706 std::cerr <<
"Unknown linsolver_type: " << linsolver_type <<
'\n';
707 throw std::runtime_error(
"Unknown linsolver_type");
709 computePressureAndFluxes(r, sat);
717 explicit FaceFluxes(
int sz)
718 : fluxes_(sz, 0.0), visited_(sz, 0), max_modification_(0.0)
721 void put(
double flux,
int f_ix) {
722 assert(visited_[f_ix] == 0 || visited_[f_ix] == 1);
723 double sign = visited_[f_ix] ? -1.0 : 1.0;
724 fluxes_[f_ix] += sign*flux;
727 void get(
double& flux,
int f_ix) {
728 assert(visited_[f_ix] == 0 || visited_[f_ix] == 1);
729 double sign = visited_[f_ix] ? -1.0 : 1.0;
730 double new_flux = 0.5*sign*fluxes_[f_ix];
731 double diff = std::fabs(flux - new_flux);
732 max_modification_ = std::max(max_modification_, diff);
738 std::ranges::fill(visited_, 0);
741 double maxMod()
const
743 return max_modification_;
746 std::vector<double> fluxes_;
747 std::vector<int> visited_;
748 double max_modification_;
763 typedef typename GridInterface::CellIterator CI;
764 typedef typename CI ::FaceIterator FI;
765 const std::vector<int>& cell = flowSolution_.cellno_;
766 const Opm::SparseTable<int>& cf = flowSolution_.cellFaces_;
767 Opm::SparseTable<double>& cflux = flowSolution_.outflux_;
769 FaceFluxes face_fluxes(pgrid_->numberOfFaces());
771 for (CI c = pgrid_->cellbegin(); c != pgrid_->cellend(); ++c) {
772 const int cell_index = cell[c->index()];
773 for (FI f = c->facebegin(); f != c->faceend(); ++f) {
774 int f_ix = cf[cell_index][f->localIndex()];
775 double flux = cflux[cell_index][f->localIndex()];
777 if (ppartner_dof_.empty()) {
780 int partner_f_ix = ppartner_dof_[f_ix];
781 if (partner_f_ix != -1) {
782 face_fluxes.put(flux, f_ix);
783 face_fluxes.put(flux, partner_f_ix);
786 face_fluxes.put(flux, f_ix);
790 face_fluxes.resetVisited();
792 for (CI c = pgrid_->cellbegin(); c != pgrid_->cellend(); ++c) {
793 const int cell_index = cell[c->index()];
794 for (FI f = c->facebegin(); f != c->faceend(); ++f) {
795 int f_ix = cf[cell_index][f->localIndex()];
796 double& flux = cflux[cell_index][f->localIndex()];
798 if (ppartner_dof_.empty()) {
801 int partner_f_ix = ppartner_dof_[f_ix];
802 if (partner_f_ix != -1) {
803 face_fluxes.get(flux, f_ix);
805 face_fluxes.get(dummy, partner_f_ix);
806 assert(dummy == flux);
809 face_fluxes.get(flux, f_ix);
813 return face_fluxes.maxMod();
837 return flowSolution_;
855 template<
typename charT,
class traits>
858 os <<
"IncompFlowSolverHybrid<>:\n"
859 <<
"\tMaximum number of cell faces = " << max_ncf_ <<
'\n'
860 <<
"\tNumber of internal faces = " << num_internal_faces_ <<
'\n'
861 <<
"\tTotal number of faces = " << total_num_faces_ <<
'\n';
863 const std::vector<int>& cell = flowSolution_.cellno_;
864 os <<
"cell index map = [";
865 std::ranges::copy(cell, std::ostream_iterator<int>(os,
" "));
868 const Opm::SparseTable<int>& cf = flowSolution_.cellFaces_;
869 os <<
"cell faces =\n";
870 for (
int i = 0; i < cf.size(); ++i)
872 os <<
"\t[" << i <<
"] -> [";
873 std::ranges::copy(cf[i], std::ostream_iterator<int>(os,
","));
902 writeMatrixToMatlab(S_, prefix +
"-mat.dat");
904 std::string rhsfile(prefix +
"-rhs.dat");
905 std::ofstream rhs(rhsfile.c_str());
907 rhs.setf(std::ios::scientific | std::ios::showpos);
908 std::ranges::copy(rhs_, std::ostream_iterator<VectorBlockType>(rhs,
"\n"));
912 typedef std::pair<int,int> DofID;
913 typedef std::unordered_map<int,DofID> BdryIdMapType;
914 typedef BdryIdMapType::const_iterator BdryIdMapIterator;
916 const GridInterface* pgrid_;
917 BdryIdMapType bdry_id_map_;
918 std::vector<int> ppartner_dof_;
920 InnerProduct<GridInterface, RockInterface> ip_;
925 int num_internal_faces_;
926 int total_num_faces_;
929 std::vector<Scalar> L_, g_;
930 Opm::SparseTable<Scalar> F_ ;
934 typedef Dune::FieldVector<Scalar, 1 > VectorBlockType;
935 typedef Dune::FieldMatrix<Scalar, 1, 1> MatrixBlockType;
937 Dune::BCRSMatrix <MatrixBlockType> S_;
938 Dune::BlockVector<VectorBlockType> rhs_;
939 Dune::BlockVector<VectorBlockType> soln_;
940 bool matrix_structure_valid_;
941 bool do_regularization_;
945 FlowSolution flowSolution_;
949 void enumerateDof(
const GridInterface& g,
const BCInterface& bc)
953 enumerateBCDof(g, bc);
956 cleared_state_ =
false;
960 void enumerateGridDof(
const GridInterface& g)
963 typedef typename GridInterface::CellIterator CI;
964 typedef typename CI ::FaceIterator FI;
966 const int nc = g.numberOfCells();
967 std::vector<int> fpos ; fpos.reserve(nc + 1);
968 std::vector<int> num_cf(nc) ;
969 std::vector<int> faces ;
972 std::vector<int>(nc, -1).swap(flowSolution_.cellno_);
974 std::vector<int>& cell = flowSolution_.cellno_;
977 int cellno = 0; fpos.push_back(0);
979 for (CI c = g.cellbegin(); c != g.cellend(); ++c, ++cellno) {
980 const int c0 = c->index();
981 assert((0 <= c0) && (c0 < nc) && (cell[c0] == -1));
985 int& ncf = num_cf[c0];
987 for (FI f = c->facebegin(); f != c-> faceend(); ++f) {
988 if (!f->boundary()) {
989 const int c1 = f->neighbourCellIndex();
990 assert((0 <= c1) && (c1 < nc) && (c1 != c0));
992 if (cell[c1] == -1) {
1000 fpos.push_back(
int(faces.size()));
1001 max_ncf_ = std::max(max_ncf_, ncf);
1004 assert (cellno == nc);
1006 total_num_faces_ = num_internal_faces_ = int(faces.size());
1008 ip_.init(max_ncf_); ip_.reserveMatrices(num_cf);
1009 F_ .reserve(nc, tot_ncf);
1011 flowSolution_.cellFaces_.reserve(nc, tot_ncf);
1012 flowSolution_.outflux_ .reserve(nc, tot_ncf);
1014 Opm::SparseTable<int>& cf = flowSolution_.cellFaces_;
1017 std::vector<int> l2g; l2g .reserve(max_ncf_);
1018 std::vector<Scalar> F_alloc; F_alloc .reserve(max_ncf_);
1021 typedef std::vector<int>::iterator VII;
1022 for (CI c = g.cellbegin(); c != g.cellend(); ++c) {
1023 const int c0 = c->index();
1025 assert ((0 <= c0 ) && ( c0 < nc) &&
1026 (0 <= cell[c0]) && (cell[c0] < nc));
1028 const int ncf = num_cf[cell[c0]];
1029 l2g .resize(ncf , 0 );
1030 F_alloc .resize(ncf , Scalar(0.0));
1032 for (FI f = c->facebegin(); f != c->faceend(); ++f) {
1033 if (f->boundary()) {
1035 l2g[f->localIndex()] = total_num_faces_++;
1047 const int c1 = f->neighbourCellIndex();
1048 assert ((0 <= c1 ) && ( c1 < nc) &&
1049 (0 <= cell[c1]) && (cell[c1] < nc));
1051 int t = c0, seek = c1;
1052 if (cell[seek] < cell[t])
1055 int s = fpos[cell[t]], e = fpos[cell[t] + 1];
1057 VII p = std::find(faces.begin() + s, faces.begin() + e, seek);
1058 assert(p != faces.begin() + e);
1060 l2g[f->localIndex()] = s + (p - (faces.begin() + s));
1064 cf.appendRow (l2g .begin(), l2g .end());
1065 F_.appendRow (F_alloc.begin(), F_alloc.end());
1067 flowSolution_.outflux_
1068 .appendRow(F_alloc.begin(), F_alloc.end());
1074 void enumerateBCDof(
const GridInterface& g,
const BCInterface& bc)
1077 typedef typename GridInterface::CellIterator CI;
1078 typedef typename CI ::FaceIterator FI;
1080 const std::vector<int>& cell = flowSolution_.cellno_;
1081 const Opm::SparseTable<int>& cf = flowSolution_.cellFaces_;
1083 bdry_id_map_.clear();
1084 for (CI c = g.cellbegin(); c != g.cellend(); ++c) {
1085 for (FI f = c->facebegin(); f != c->faceend(); ++f) {
1086 if (f->boundary() && bc.flowCond(*f).isPeriodic()) {
1087 const int bid = f->boundaryId();
1088 DofID dof(cell[c->index()], f->localIndex());
1089 bdry_id_map_.insert(std::make_pair(bid, dof));
1094 ppartner_dof_.clear();
1095 if (!bdry_id_map_.empty()) {
1096 ppartner_dof_.assign(total_num_faces_, -1);
1097 for (CI c = g.cellbegin(); c != g.cellend(); ++c) {
1098 for (FI f = c->facebegin(); f != c->faceend(); ++f) {
1099 if (f->boundary() && bc.flowCond(*f).isPeriodic()) {
1100 const int dof1 = cf[cell[c->index()]][f->localIndex()];
1102 BdryIdMapIterator
j =
1103 bdry_id_map_.find(bc.getPeriodicPartner(f->boundaryId()));
1104 assert (
j != bdry_id_map_.end());
1105 const int dof2 = cf[
j->second.first][
j->second.second];
1107 ppartner_dof_[dof1] = dof2;
1108 ppartner_dof_[dof2] = dof1;
1118 void allocateConnections(
const BCInterface& bc)
1122 assert(!cleared_state_);
1124 assert (!matrix_structure_valid_);
1127 S_.setSize(total_num_faces_, total_num_faces_);
1128 S_.setBuildMode(Dune::BCRSMatrix<MatrixBlockType>::random);
1131 for (
int f = 0; f < total_num_faces_; ++f) {
1132 S_.setrowsize(f, 1);
1135 allocateGridConnections();
1136 allocateBCConnections(bc);
1140 rhs_ .resize(total_num_faces_);
1141 soln_.resize(total_num_faces_);
1146 void allocateGridConnections()
1149 const Opm::SparseTable<int>& cf = flowSolution_.cellFaces_;
1151 for (
int c = 0; c < cf.size(); ++c) {
1152 const int nf = cf[c].size();
1153 for (
auto& f : cf[c]) {
1154 S_.incrementrowsize(f, nf - 1);
1161 void allocateBCConnections(
const BCInterface& bc)
1177 typedef typename GridInterface::CellIterator CI;
1178 typedef typename CI ::FaceIterator FI;
1180 const std::vector<int>& cell = flowSolution_.cellno_;
1181 const Opm::SparseTable<int>& cf = flowSolution_.cellFaces_;
1183 if (!bdry_id_map_.empty()) {
1186 for (CI c = pgrid_->cellbegin(); c != pgrid_->cellend(); ++c) {
1187 for (FI f = c->facebegin(); f != c->faceend(); ++f) {
1188 if (f->boundary() && bc.flowCond(*f).isPeriodic()) {
1190 const int dof1 = cf[cell[c->index()]][f->localIndex()];
1193 BdryIdMapIterator
j =
1194 bdry_id_map_.find(bc.getPeriodicPartner(f->boundaryId()));
1195 assert (
j != bdry_id_map_.end());
1196 const int c2 =
j->second.first;
1197 const int dof2 = cf[c2][
j->second.second];
1203 const int ndof = cf.rowSize(c2);
1204 S_.incrementrowsize(dof1, ndof);
1205 for (
int dof = 0; dof < ndof; ++dof) {
1206 int ii = cf[c2][dof];
1207 int pp = ppartner_dof_[ii];
1208 if ((pp != -1) && (pp != dof1) && (pp < ii)) {
1209 S_.incrementrowsize(pp, 1);
1211 S_.incrementrowsize(ii, 1);
1223 void setConnections(
const BCInterface& bc)
1226 setGridConnections();
1227 setBCConnections(bc);
1231 const int nc = pgrid_->numberOfCells();
1232 std::vector<Scalar>(nc).swap(flowSolution_.pressure_);
1233 std::vector<Scalar>(nc).swap(g_);
1234 std::vector<Scalar>(nc).swap(L_);
1236 matrix_structure_valid_ =
true;
1241 void setGridConnections()
1244 const Opm::SparseTable<int>& cf = flowSolution_.cellFaces_;
1247 for (
int c = 0; c < cf.size(); ++c) {
1248 auto fb = cf[c].begin(), fe = cf[c].end();
1250 for (
auto i = fb; i != fe; ++i) {
1251 for (
auto j = fb;
j != fe; ++
j) {
1252 S_.addindex(*i, *
j);
1260 void setBCConnections(
const BCInterface& bc)
1276 typedef typename GridInterface::CellIterator CI;
1277 typedef typename CI ::FaceIterator FI;
1279 const std::vector<int>& cell = flowSolution_.cellno_;
1280 const Opm::SparseTable<int>& cf = flowSolution_.cellFaces_;
1282 if (!bdry_id_map_.empty()) {
1285 for (CI c = pgrid_->cellbegin(); c != pgrid_->cellend(); ++c) {
1286 for (FI f = c->facebegin(); f != c->faceend(); ++f) {
1287 if (f->boundary() && bc.flowCond(*f).isPeriodic()) {
1289 const int dof1 = cf[cell[c->index()]][f->localIndex()];
1292 BdryIdMapIterator
j =
1293 bdry_id_map_.find(bc.getPeriodicPartner(f->boundaryId()));
1294 assert (
j != bdry_id_map_.end());
1295 const int c2 =
j->second.first;
1296 const int dof2 = cf[c2][
j->second.second];
1301 const int ndof = cf.rowSize(c2);
1302 for (
int dof = 0; dof < ndof; ++dof) {
1303 int ii = cf[c2][dof];
1304 int pp = ppartner_dof_[ii];
1305 if ((pp != -1) && (pp != dof1) && (pp < ii)) {
1308 S_.addindex(dof1, ii);
1309 S_.addindex(ii, dof1);
1310 S_.addindex(dof2, ii);
1311 S_.addindex(ii, dof2);
1323 template<
class Flu
idInterface>
1324 void assembleDynamic(
const FluidInterface& fl ,
1325 const std::vector<double>& sat,
1326 const BCInterface& bc ,
1327 const std::vector<double>& src)
1330 typedef typename GridInterface::CellIterator CI;
1332 const std::vector<int>& cell = flowSolution_.cellno_;
1333 const Opm::SparseTable<int>& cf = flowSolution_.cellFaces_;
1335 std::vector<Scalar> data_store(max_ncf_ * max_ncf_);
1336 std::vector<Scalar> e (max_ncf_);
1337 std::vector<Scalar> rhs (max_ncf_);
1338 std::vector<Scalar> gflux (max_ncf_);
1340 std::vector<FaceType> facetype (max_ncf_);
1341 std::vector<Scalar> condval (max_ncf_);
1342 std::vector<int> ppartner (max_ncf_);
1348 std::ranges::fill(g_, Scalar(0.0));
1349 std::ranges::fill(L_, Scalar(0.0));
1350 std::ranges::fill(e, Scalar(1.0));
1354 do_regularization_ =
true;
1357 for (CI c = pgrid_->cellbegin(); c != pgrid_->cellend(); ++c) {
1358 const int ci = c->index();
1359 const int c0 = cell[ci]; assert (c0 < cf.size());
1360 const int nf = cf[c0].size();
1365 setExternalContrib(c, c0, bc, src[ci], rhs,
1366 facetype, condval, ppartner);
1368 ip_.computeDynamicParams(c, fl, sat);
1371 ip_.getInverseMatrix(c, S);
1373 std::ranges::fill(gflux, Scalar(0.0));
1374 ip_.gravityFlux(c, gflux);
1377 buildCellContrib(c0, one, gflux, S, rhs);
1379 addCellContrib(S, rhs, facetype, condval, ppartner, cf[c0]);
1386 void solveLinearSystem(
double residual_tolerance,
int verbosity_level,
int maxit)
1390 Scalar residTol = residual_tolerance;
1392 typedef Dune::BCRSMatrix <MatrixBlockType> MatrixT;
1393 typedef Dune::BlockVector<VectorBlockType> VectorT;
1394 typedef Dune::MatrixAdapter<Matrix,Vector,Vector> Adapter;
1397 if (do_regularization_) {
1403 Dune::SeqILU<MatrixT,VectorT,VectorT> precond(S_, 1.0);
1406 Dune::CGSolver<VectorT> linsolve(opS, precond, residTol,
1407 (maxit>0)?maxit:S_.N(), verbosity_level);
1409 Dune::InverseOperatorResult result;
1414 linsolve.apply(soln_, rhs_, result);
1415 if (!result.converged) {
1416 OPM_THROW(std::runtime_error,
1417 "Linear solver failed to converge in " +
1418 std::to_string(result.iterations) +
" iterations.\n"
1419 "Residual reduction achieved is " +
1420 std::to_string(result.reduction) +
'\n');
1429 typedef Dune::BCRSMatrix <MatrixBlockType> Matrix;
1430 typedef Dune::BlockVector<VectorBlockType> Vector;
1431 typedef Dune::MatrixAdapter<Matrix,Vector,Vector> Operator;
1436#ifndef FIRST_DIAGONAL
1437#define FIRST_DIAGONAL 1
1443#define SMOOTHER_ILU 1
1446#define SMOOTHER_BGS 0
1448#ifndef ANISOTROPIC_3D
1449#define ANISOTROPIC_3D 0
1453 typedef Dune::Amg::FirstDiagonal CouplingMetric;
1455 typedef Dune::Amg::RowSum CouplingMetric;
1459 typedef Dune::Amg::SymmetricCriterion<Matrix,CouplingMetric> CriterionBase;
1461 typedef Dune::Amg::UnSymmetricCriterion<Matrix,CouplingMetric> CriterionBase;
1465 typedef Dune::SeqOverlappingSchwarz<Matrix,Vector,Dune::MultiplicativeSchwarzMode> Smoother;
1468 typedef Dune::SeqILU<Matrix,Vector,Vector> Smoother;
1470 typedef Dune::SeqSSOR<Matrix,Vector,Vector> Smoother;
1473 typedef Dune::Amg::CoarsenCriterion<CriterionBase> Criterion;
1477 std::unique_ptr<Operator> opS_;
1478 typedef Dune::Preconditioner<Vector,Vector> PrecondBase;
1479 std::unique_ptr<PrecondBase> precond_;
1483 void solveLinearSystemAMG(
double residual_tolerance,
int verbosity_level,
1484 int maxit,
double prolong_factor,
bool same_matrix,
int smooth_steps)
1487 typedef Dune::Amg::AMG<Operator,Vector,Smoother,Dune::Amg::SequentialInformation>
1491 Scalar residTol = residual_tolerance;
1495 if (do_regularization_) {
1498 opS_.reset(
new Operator(S_));
1502 typename Precond::SmootherArgs smootherArgs;
1503 smootherArgs.relaxationFactor = relax;
1505 smootherArgs.overlap = Precond::SmootherArgs::none;
1506 smootherArgs.onthefly =
false;
1508 Criterion criterion;
1509 criterion.setDebugLevel(verbosity_level);
1511 criterion.setDefaultValuesAnisotropic(3, 2);
1513 criterion.setProlongationDampingFactor(prolong_factor);
1514 criterion.setBeta(1e-10);
1515 criterion.setNoPreSmoothSteps(smooth_steps);
1516 criterion.setNoPostSmoothSteps(smooth_steps);
1517 criterion.setGamma(1);
1518 precond_.reset(
new Precond(*opS_, criterion, smootherArgs));
1521 Dune::CGSolver<Vector> linsolve(*opS_,
dynamic_cast<Precond&
>(*precond_), residTol, (maxit>0)?maxit:S_.N(), verbosity_level);
1523 Dune::InverseOperatorResult result;
1527 typedef typename Dune::BCRSMatrix <MatrixBlockType>::ConstRowIterator RowIter;
1528 typedef typename Dune::BCRSMatrix <MatrixBlockType>::ConstColIterator ColIter;
1529 for(RowIter ri=S_.begin(); ri!=S_.end(); ++ri){
1530 bool isDirichlet=
true;
1531 for(ColIter ci=ri->begin(); ci!=ri->end(); ++ci)
1532 if(ci.index()!=ri.index() && *ci!=0.0)
1535 soln_[ri.index()]=rhs_[ri.index()]/S_[ri.index()][ri.index()];
1539 linsolve.apply(soln_, rhs_, result);
1540 if (!result.converged) {
1541 OPM_THROW(std::runtime_error,
1542 "Linear solver failed to converge in " +
1543 std::to_string(result.iterations) +
" iterations.\n"
1544 "Residual reduction achieved is " +
1545 std::to_string(result.reduction) +
'\n');
1552 void solveLinearSystemFastAMG(
double residual_tolerance,
int verbosity_level,
1553 int maxit,
double prolong_factor,
bool same_matrix,
int smooth_steps)
1556 typedef Dune::Amg::FastAMG<Operator,Vector> Precond;
1559 Scalar residTol = residual_tolerance;
1563 if (do_regularization_) {
1569 typedef Dune::Amg::AggregationCriterion<Dune::Amg::SymmetricMatrixDependency<Matrix,CouplingMetric> > CritBase;
1571 typedef Dune::Amg::CoarsenCriterion<CritBase> Crit;
1573 criterion.setDebugLevel(verbosity_level);
1575 criterion.setDefaultValuesAnisotropic(3, 2);
1577 criterion.setProlongationDampingFactor(prolong_factor);
1578 criterion.setBeta(1e-10);
1579 Dune::Amg::Parameters parms;
1580 parms.setDebugLevel(verbosity_level);
1581 parms.setNoPreSmoothSteps(smooth_steps);
1582 parms.setNoPostSmoothSteps(smooth_steps);
1583 precond_.reset(
new Precond(*opS_, criterion, parms));
1586 Dune::GeneralizedPCGSolver<Vector> linsolve(*opS_,
dynamic_cast<Precond&
>(*precond_), residTol, (maxit>0)?maxit:S_.N(), verbosity_level);
1588 Dune::InverseOperatorResult result;
1593 typedef typename Dune::BCRSMatrix <MatrixBlockType>::ConstRowIterator RowIter;
1594 typedef typename Dune::BCRSMatrix <MatrixBlockType>::ConstColIterator ColIter;
1595 for(RowIter ri=S_.begin(); ri!=S_.end(); ++ri){
1596 bool isDirichlet=
true;
1597 for(ColIter ci=ri->begin(); ci!=ri->end(); ++ci)
1598 if(ci.index()!=ri.index() && *ci!=0.0)
1601 soln_[ri.index()]=rhs_[ri.index()]/S_[ri.index()][ri.index()];
1605 linsolve.apply(soln_, rhs_, result);
1606 if (!result.converged) {
1607 OPM_THROW(std::runtime_error,
1608 "Linear solver failed to converge in " +
1609 std::to_string(result.iterations) +
" iterations.\n"
1610 "Residual reduction achieved is " +
1611 std::to_string(result.reduction) +
'\n');
1617 void solveLinearSystemKAMG(
double residual_tolerance,
int verbosity_level,
1618 int maxit,
double prolong_factor,
bool same_matrix,
int smooth_steps)
1622 typedef Dune::Amg::KAMG<Operator,Vector,Smoother,Dune::Amg::SequentialInformation> Precond;
1624 Scalar residTol = residual_tolerance;
1627 if (do_regularization_) {
1634 typename Precond::SmootherArgs smootherArgs;
1635 smootherArgs.relaxationFactor = relax;
1637 smootherArgs.overlap = Precond::SmootherArgs::none;
1638 smootherArgs.onthefly =
false;
1640 Criterion criterion;
1641 criterion.setDebugLevel(verbosity_level);
1643 criterion.setDefaultValuesAnisotropic(3, 2);
1645 criterion.setProlongationDampingFactor(prolong_factor);
1646 criterion.setBeta(1e-10);
1647 criterion.setNoPreSmoothSteps(smooth_steps);
1648 criterion.setNoPostSmoothSteps(smooth_steps);
1649 criterion.setGamma(2);
1650 precond_.reset(
new Precond(*opS_, criterion, smootherArgs));
1653 Dune::CGSolver<Vector> linsolve(*opS_,
dynamic_cast<Precond&
>(*precond_), residTol, (maxit>0)?maxit:S_.N(), verbosity_level);
1655 Dune::InverseOperatorResult result;
1659 typedef typename Dune::BCRSMatrix <MatrixBlockType>::ConstRowIterator RowIter;
1660 typedef typename Dune::BCRSMatrix <MatrixBlockType>::ConstColIterator ColIter;
1661 for(RowIter ri=S_.begin(); ri!=S_.end(); ++ri){
1662 bool isDirichlet=
true;
1663 for(ColIter ci=ri->begin(); ci!=ri->end(); ++ci)
1664 if(ci.index()!=ri.index() && *ci!=0.0)
1667 soln_[ri.index()]=rhs_[ri.index()]/S_[ri.index()][ri.index()];
1671 linsolve.apply(soln_, rhs_, result);
1672 if (!result.converged) {
1673 OPM_THROW(std::runtime_error,
1674 "Linear solver failed to converge in " +
1675 std::to_string(result.iterations) +
" iterations.\n"
1676 "Residual reduction achieved is " +
1677 std::to_string(result.reduction) +
'\n');
1685 template<
class Flu
idInterface>
1686 void computePressureAndFluxes(
const FluidInterface& r ,
1687 const std::vector<double>& sat)
1690 typedef typename GridInterface::CellIterator CI;
1692 const std::vector<int>& cell = flowSolution_.cellno_;
1693 const Opm::SparseTable<int>& cf = flowSolution_.cellFaces_;
1695 std::vector<Scalar>& p = flowSolution_.pressure_;
1696 Opm::SparseTable<Scalar>& v = flowSolution_.outflux_;
1699 std::vector<double> pi (max_ncf_);
1700 std::vector<double> gflux(max_ncf_);
1701 std::vector<double> Binv_storage(max_ncf_ * max_ncf_);
1704 for (CI c = pgrid_->cellbegin(); c != pgrid_->cellend(); ++c) {
1705 const int c0 = cell[c->index()];
1706 const int nf = cf.rowSize(c0);
1712 for (
int i = 0; i < nf; ++i) {
1713 pi[i] = soln_[cf[c0][i]];
1718 std::inner_product(F_[c0].begin(), F_[c0].end(),
1719 pi.begin(), 0.0)) / L_[c0];
1721 std::ranges::transform(pi, pi.begin(),
1722 [&p, c0](
const double& input) { return p[c0] - input; });
1729 ip_.computeDynamicParams(c, r, sat);
1732 ip_.getInverseMatrix(c, Binv);
1733 vecMulAdd_N(Scalar(1.0), Binv, &pi[0], Scalar(0.0), &v[c0][0]);
1737 ip_.gravityFlux(c, gflux);
1738 std::ranges::transform(gflux, v[c0], v[c0].begin(), std::plus<Scalar>());
1746 void setExternalContrib(
const typename GridInterface::CellIterator c,
1747 const int c0,
const BCInterface& bc,
1749 std::vector<Scalar>& rhs,
1750 std::vector<FaceType>& facetype,
1751 std::vector<double>& condval,
1752 std::vector<int>& ppartner)
1755 typedef typename GridInterface::CellIterator::FaceIterator FI;
1757 const Opm::SparseTable<int>& cf = flowSolution_.cellFaces_;
1759 std::ranges::fill(rhs, Scalar(0.0));
1760 std::ranges::fill(facetype, Internal);
1761 std::ranges::fill(condval, Scalar(0.0));
1762 std::ranges::fill(ppartner, -1);
1767 for (FI f = c->facebegin(); f != c->faceend(); ++f, ++
k) {
1768 if (f->boundary()) {
1769 const FlowBC& bcond = bc.flowCond(*f);
1770 if (bcond.isDirichlet()) {
1771 facetype[
k] = Dirichlet;
1772 condval[
k] = bcond.pressure();
1773 do_regularization_ =
false;
1774 }
else if (bcond.isPeriodic()) {
1775 BdryIdMapIterator
j =
1776 bdry_id_map_.find(bc.getPeriodicPartner(f->boundaryId()));
1777 assert (
j != bdry_id_map_.end());
1779 facetype[
k] = Periodic;
1780 condval[
k] = bcond.pressureDifference();
1781 ppartner[
k] = cf[
j->second.first][
j->second.second];
1783 assert (bcond.isNeumann());
1784 facetype[
k] = Neumann;
1785 rhs[
k] = bcond.outflux();
1795 void buildCellContrib(
const int c ,
1797 const std::vector<Scalar>& gflux,
1799 std::vector<Scalar>& rhs)
1806 L_[c] = std::accumulate(Ft.data(), Ft.data() + Ft.numRows(), 0.0);
1807 g_[c] -= std::accumulate(gflux.begin(), gflux.end(), Scalar(0.0));
1810 std::ranges::transform(gflux, rhs, rhs.begin(), std::minus<Scalar>());
1813 std::transform(rhs.begin(), rhs.end(), Ft.data(),
1815 axpby<Scalar>(Scalar(1.0), Scalar(g_[c] / L_[c])));
1827 const std::vector<Scalar>& rhs ,
1828 const std::vector<FaceType>& facetype,
1829 const std::vector<Scalar>& condval ,
1830 const std::vector<int>& ppartner,
1835 for (
auto i = l2g.begin(); i != l2g.end(); ++i, ++r) {
1839 switch (facetype[r]) {
1845 S_ [ii][ii] = S(r,r);
1846 rhs_[ii] = S(r,r) * condval[r];
1861 assert ((0 <= ppartner[r]) && (ppartner[r] <
int(rhs_.size())));
1862 assert (ii != ppartner[r]);
1865 const double a = S(r,r), b =
a * condval[r];
1869 S_ [ ii][ppartner[r]] -=
a;
1873 S_ [ppartner[r]][ ii] -=
a;
1874 S_ [ppartner[r]][ppartner[r]] +=
a;
1875 rhs_[ppartner[r]] -= b;
1885 for (
auto j = l2g.begin();
j != l2g.end(); ++
j, ++c) {
1889 if (facetype[c] == Dirichlet) {
1890 rhs_[ii] -= S(r,c) * condval[c];
1893 if (facetype[c] == Periodic) {
1894 assert ((0 <= ppartner[c]) && (ppartner[c] <
int(rhs_.size())));
1895 assert (jj != ppartner[c]);
1896 if (ppartner[c] < jj) {
1897 rhs_[ii] -= S(r,c) * condval[c];
1901 S_[ii][jj] += S(r,c);
void const int const int const double * a
Definition: blas_lapack.hpp:79
void const int const int const double const double const int const double const int const double double * y
Definition: blas_lapack.hpp:59
void const int const int const int * k
Definition: blas_lapack.hpp:64
void const int const int * n
Definition: blas_lapack.hpp:56
void const int const int const double const double const int const double * x
Definition: blas_lapack.hpp:58
Solve mixed formulation of incompressible flow modelled by Darcy's law.
Definition: IncompFlowSolverHybrid.hpp:379
double postProcessFluxes()
Postprocess the solution fluxes. This method modifies the solution object so that out-fluxes of twin ...
Definition: IncompFlowSolverHybrid.hpp:761
SolutionType getSolution()
Recover the solution to the problem defined by the parameters to method.
Definition: IncompFlowSolverHybrid.hpp:835
void solve(const FluidInterface &r, const std::vector< double > &sat, const BCInterface &bc, const std::vector< double > &src, double residual_tolerance=1e-8, int linsolver_verbosity=1, int linsolver_type=1, bool same_matrix=false, int linsolver_maxit=0, double prolongate_factor=1.6, int smooth_steps=1)
Construct and solve system of linear equations for the pressure values on each interface/contact betw...
Definition: IncompFlowSolverHybrid.hpp:672
const FlowSolution & SolutionType
Type representing the solution to the problem defined by the parameters to.
Definition: IncompFlowSolverHybrid.hpp:825
void printSystem(const std::string &prefix)
Output current system of linear equations to permanent storage in files. One file for the coefficient...
Definition: IncompFlowSolverHybrid.hpp:900
void computeInnerProducts(const RockInterface &r, const Point &grav)
Compute static (i.e., independent of saturation) parts of the spatially varying inner products for e...
Definition: IncompFlowSolverHybrid.hpp:579
void clear()
Clear all topologic, geometric and rock-dependent information currently held in internal data structu...
Definition: IncompFlowSolverHybrid.hpp:511
void printStats(std::basic_ostream< charT, traits > &os)
Print statistics about the connections in the current model. This is mostly for debugging purposes an...
Definition: IncompFlowSolverHybrid.hpp:856
void init(const GridInterface &g, const RockInterface &r, const Point &grav, const BCInterface &bc)
All-in-one initialization routine. Enumerates all grid connections, allocates sufficient space,...
Definition: IncompFlowSolverHybrid.hpp:490
void initSystemStructure(const GridInterface &g, const BCInterface &bc)
Compute structure of coefficient matrix in final system of linear equations for this flow problem....
Definition: IncompFlowSolverHybrid.hpp:547
Dune::MatrixAdapter< Matrix, Vector, Vector > Operator
A linear operator.
Definition: elasticity_preconditioners.hpp:49
min[0]
Definition: elasticity_upscale_impl.hpp:146
int j
Definition: elasticity_upscale_impl.hpp:658
Definition: ImplicitAssembly.hpp:43
void symmetricUpdate(const T &a1, const FullMatrix< T, StoragePolicy, FortranOrdering > &A, const T &a2, FullMatrix< T, StoragePolicy, FortranOrdering > &C)
Symmetric, rank update of symmetric matrix. Specifically, .
Definition: Matrix.hpp:829
const FullMatrix< double, ImmutableSharedData, FortranOrdering > ImmutableFortranMatrix
Definition: Matrix.hpp:590
void vecMulAdd_N(const T &a1, const FullMatrix< T, SP, FortranOrdering > &A, const std::vector< T > &x, const T &a2, std::vector< T > &y)
GEneral Matrix-Vector product (GAXPY operation). Specifically, .
Definition: Matrix.hpp:913
FullMatrix< double, SharedData, FortranOrdering > SharedFortranMatrix
Definition: Matrix.hpp:589
void matMulAdd_TN(const T &a1, const FullMatrix< T, SP1, FortranOrdering > &A, const FullMatrix< T, SP2, FortranOrdering > &B, const T &a2, FullMatrix< T, SP3, FortranOrdering > &C)
GEneral Matrix-Matrix product update of other matrix. Specificlly, .
Definition: Matrix.hpp:1252