IncompFlowSolverHybrid.hpp
Go to the documentation of this file.
1// -*- tab-width: 8; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set ts=8 sw=4 et sts=4:
3//===========================================================================
4//
5// File: IncompFlowSolverHybrid.hpp
6//
7// Created: Tue Jun 30 10:25:40 2009
8//
9// Author(s): B�rd Skaflestad <bard.skaflestad@sintef.no>
10// Atgeirr F Rasmussen <atgeirr@sintef.no>
11//
12// $Date$
13//
14// $Revision$
15//
16//===========================================================================
17
18/*
19 Copyright 2009, 2010 SINTEF ICT, Applied Mathematics.
20 Copyright 2009, 2010 Statoil ASA.
21
22 This file is part of The Open Reservoir Simulator Project (OpenRS).
23
24 OpenRS is free software: you can redistribute it and/or modify
25 it under the terms of the GNU General Public License as published by
26 the Free Software Foundation, either version 3 of the License, or
27 (at your option) any later version.
28
29 OpenRS is distributed in the hope that it will be useful,
30 but WITHOUT ANY WARRANTY; without even the implied warranty of
31 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
32 GNU General Public License for more details.
33
34 You should have received a copy of the GNU General Public License
35 along with OpenRS. If not, see <http://www.gnu.org/licenses/>.
36*/
37
38#ifndef OPENRS_INCOMPFLOWSOLVERHYBRID_HEADER
39#define OPENRS_INCOMPFLOWSOLVERHYBRID_HEADER
40
41#include <opm/common/ErrorMacros.hpp>
42#include <opm/grid/utility/SparseTable.hpp>
45
46#include <opm/common/utility/platform_dependent/disable_warnings.h>
47
48#include <unordered_map>
49
50#include <dune/common/fvector.hh>
51#include <dune/common/fmatrix.hh>
52
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>
67
68#include <opm/common/utility/platform_dependent/reenable_warnings.h>
69
70
71#include <algorithm>
72#include <functional>
73#include <map>
74#include <memory>
75#include <numeric>
76#include <ostream>
77#include <stdexcept>
78#include <utility>
79#include <vector>
80#include <iostream>
81
82namespace Opm {
83 namespace {
107 template<class GI>
108 bool topologyIsSane(const GI& g)
109 {
110 typedef typename GI::CellIterator CI;
111 typedef typename CI::FaceIterator FI;
112
113 bool sane = g.numberOfCells() >= 0;
114
115 for (CI c = g.cellbegin(); sane && c != g.cellend(); ++c) {
116 std::vector<int> n;
117
118 for (FI f = c->facebegin(); f != c->faceend(); ++f) {
119 if (!f->boundary()) {
120 n.push_back(f->neighbourCellIndex());
121 }
122 }
123 std::ranges::sort(n);
124
125 sane = std::unique(n.begin(), n.end()) == n.end();
126 }
127
128 return sane;
129 }
130
131
157 template<typename T>
158 class axpby {
159 public:
165 axpby(const T& a, const T& b) : a_(a), b_(b) {}
166
179 T operator()(const T& x, const T& y)
180 {
181 return a_*x + b_*y;
182 }
183 private:
184 T a_, b_;
185 };
186 }
187
188
375 template <class GridInterface,
376 class RockInterface,
377 class BCInterface,
378 template <class GridIF, class RockIF> class InnerProduct>
385 typedef typename GridInterface::Scalar Scalar;
386
391 enum FaceType { Internal, Dirichlet, Neumann, Periodic };
392
395 class FlowSolution {
396 public:
402 typedef typename GridInterface::Scalar Scalar;
403
407 typedef typename GridInterface::CellIterator CI;
408
411 typedef typename CI ::FaceIterator FI;
412
413 friend class IncompFlowSolverHybrid;
414
424 Scalar pressure(const CI& c) const
425 {
426 return pressure_[cellno_[c->index()]];
427 }
428
439 Scalar outflux (const FI& f) const
440 {
441 return outflux_[cellno_[f->cellIndex()]][f->localIndex()];
442 }
443 Scalar outflux (int hf) const
444 {
445 return outflux_.data(hf);
446 }
447 private:
448 std::vector< int > cellno_;
449 Opm::SparseTable< int > cellFaces_;
450 std::vector<Scalar> pressure_;
451 Opm::SparseTable<Scalar> outflux_;
452
453 void clear() {
454 std::vector<int>().swap(cellno_);
455 cellFaces_.clear();
456
457 std::vector<Scalar>().swap(pressure_);
458 outflux_.clear();
459 }
460 };
461
462 public:
489 template<class Point>
490 void init(const GridInterface& g,
491 const RockInterface& r,
492 const Point& grav,
493 const BCInterface& bc)
494 {
495 clear();
496
497 if (g.numberOfCells() > 0) {
498 initSystemStructure(g, bc);
499 computeInnerProducts(r, grav);
500 }
501 }
502
503
511 void clear()
512 {
513 pgrid_ = 0;
514 max_ncf_ = -1;
515 num_internal_faces_ = 0;
516 total_num_faces_ = 0;
517 matrix_structure_valid_ = false;
518 do_regularization_ = true; // Assume pure Neumann by default.
519
520 bdry_id_map_.clear();
521
522 std::vector<Scalar>().swap(L_);
523 std::vector<Scalar>().swap(g_);
524 F_.clear();
525
526 flowSolution_.clear();
527
528 cleared_state_ = true;
529 }
530
531
547 void initSystemStructure(const GridInterface& g,
548 const BCInterface& bc)
549 {
550 // You must call clear() prior to initSystemStructure()
551 assert (cleared_state_);
552
553 assert (topologyIsSane(g));
554
555 enumerateDof(g, bc);
556 allocateConnections(bc);
557 setConnections(bc);
558 }
559
560
578 template<class Point>
579 void computeInnerProducts(const RockInterface& r,
580 const Point& grav)
581 {
582 // You must call connectionsCompleted() prior to computeInnerProducts()
583 assert (matrix_structure_valid_);
584
585 typedef typename GridInterface::CellIterator CI;
586 const Opm::SparseTable<int>& cf = flowSolution_.cellFaces_;
587
588 int i = 0;
589 for (CI c = pgrid_->cellbegin(); c != pgrid_->cellend(); ++c, ++i) {
590 ip_.buildStaticContrib(c, r, grav, cf.rowSize(i));
591 }
592 }
593
594
671 template<class FluidInterface>
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)
683 {
684 assembleDynamic(r, sat, bc, src);
685// static int count = 0;
686// ++count;
687// printSystem(std::string("linsys_mimetic-") + std::to_string(count));
688 switch (linsolver_type) {
689 case 0: // ILU0 preconditioned CG
690 solveLinearSystem(residual_tolerance, linsolver_verbosity, linsolver_maxit);
691 break;
692 case 1: // AMG preconditioned CG
693 solveLinearSystemAMG(residual_tolerance, linsolver_verbosity,
694 linsolver_maxit, prolongate_factor, same_matrix, smooth_steps);
695 break;
696
697 case 2: // KAMG preconditioned CG
698 solveLinearSystemKAMG(residual_tolerance, linsolver_verbosity,
699 linsolver_maxit, prolongate_factor, same_matrix,smooth_steps);
700 break;
701 case 3: // CG preconditioned with AMG that uses less memory badwidth
702 solveLinearSystemFastAMG(residual_tolerance, linsolver_verbosity,
703 linsolver_maxit, prolongate_factor, same_matrix,smooth_steps);
704 break;
705 default:
706 std::cerr << "Unknown linsolver_type: " << linsolver_type << '\n';
707 throw std::runtime_error("Unknown linsolver_type");
708 }
709 computePressureAndFluxes(r, sat);
710 }
711
712 private:
714 class FaceFluxes
715 {
716 public:
717 explicit FaceFluxes(int sz)
718 : fluxes_(sz, 0.0), visited_(sz, 0), max_modification_(0.0)
719 {
720 }
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;
725 ++visited_[f_ix];
726 }
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);
733 flux = new_flux;
734 ++visited_[f_ix];
735 }
736 void resetVisited()
737 {
738 std::ranges::fill(visited_, 0);
739 }
740
741 double maxMod() const
742 {
743 return max_modification_;
744 }
745 private:
746 std::vector<double> fluxes_;
747 std::vector<int> visited_;
748 double max_modification_;
749
750 };
751
752 public:
762 {
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_;
768
769 FaceFluxes face_fluxes(pgrid_->numberOfFaces());
770 // First pass: compute projected fluxes.
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()];
776 if (f->boundary()) {
777 if (ppartner_dof_.empty()) {
778 continue;
779 }
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);
784 }
785 } else {
786 face_fluxes.put(flux, f_ix);
787 }
788 }
789 }
790 face_fluxes.resetVisited();
791 // Second pass: set all fluxes to the projected ones.
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()];
797 if (f->boundary()) {
798 if (ppartner_dof_.empty()) {
799 continue;
800 }
801 int partner_f_ix = ppartner_dof_[f_ix];
802 if (partner_f_ix != -1) {
803 face_fluxes.get(flux, f_ix);
804 double dummy = flux;
805 face_fluxes.get(dummy, partner_f_ix);
806 assert(dummy == flux);
807 }
808 } else {
809 face_fluxes.get(flux, f_ix);
810 }
811 }
812 }
813 return face_fluxes.maxMod();
814 }
815
816
825 typedef const FlowSolution& SolutionType;
826
836 {
837 return flowSolution_;
838 }
839
840
855 template<typename charT, class traits>
856 void printStats(std::basic_ostream<charT,traits>& os)
857 {
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';
862
863 const std::vector<int>& cell = flowSolution_.cellno_;
864 os << "cell index map = [";
865 std::ranges::copy(cell, std::ostream_iterator<int>(os, " "));
866 os << "\b]\n";
867
868 const Opm::SparseTable<int>& cf = flowSolution_.cellFaces_;
869 os << "cell faces =\n";
870 for (int i = 0; i < cf.size(); ++i)
871 {
872 os << "\t[" << i << "] -> [";
873 std::ranges::copy(cf[i], std::ostream_iterator<int>(os, ","));
874 os << "\b]\n";
875 }
876 }
877
878
900 void printSystem(const std::string& prefix)
901 {
902 writeMatrixToMatlab(S_, prefix + "-mat.dat");
903
904 std::string rhsfile(prefix + "-rhs.dat");
905 std::ofstream rhs(rhsfile.c_str());
906 rhs.precision(15);
907 rhs.setf(std::ios::scientific | std::ios::showpos);
908 std::ranges::copy(rhs_, std::ostream_iterator<VectorBlockType>(rhs, "\n"));
909 }
910
911 private:
912 typedef std::pair<int,int> DofID;
913 typedef std::unordered_map<int,DofID> BdryIdMapType;
914 typedef BdryIdMapType::const_iterator BdryIdMapIterator;
915
916 const GridInterface* pgrid_;
917 BdryIdMapType bdry_id_map_;
918 std::vector<int> ppartner_dof_;
919
920 InnerProduct<GridInterface, RockInterface> ip_;
921
922 // ----------------------------------------------------------------
923 bool cleared_state_;
924 int max_ncf_;
925 int num_internal_faces_;
926 int total_num_faces_;
927
928 // ----------------------------------------------------------------
929 std::vector<Scalar> L_, g_;
930 Opm::SparseTable<Scalar> F_ ;
931
932 // ----------------------------------------------------------------
933 // Actual, assembled system of linear equations
934 typedef Dune::FieldVector<Scalar, 1 > VectorBlockType;
935 typedef Dune::FieldMatrix<Scalar, 1, 1> MatrixBlockType;
936
937 Dune::BCRSMatrix <MatrixBlockType> S_; // System matrix
938 Dune::BlockVector<VectorBlockType> rhs_; // System RHS
939 Dune::BlockVector<VectorBlockType> soln_; // System solution (contact pressure)
940 bool matrix_structure_valid_;
941 bool do_regularization_;
942
943 // ----------------------------------------------------------------
944 // Physical quantities (derived)
945 FlowSolution flowSolution_;
946
947
948 // ----------------------------------------------------------------
949 void enumerateDof(const GridInterface& g, const BCInterface& bc)
950 // ----------------------------------------------------------------
951 {
952 enumerateGridDof(g);
953 enumerateBCDof(g, bc);
954
955 pgrid_ = &g;
956 cleared_state_ = false;
957 }
958
959 // ----------------------------------------------------------------
960 void enumerateGridDof(const GridInterface& g)
961 // ----------------------------------------------------------------
962 {
963 typedef typename GridInterface::CellIterator CI;
964 typedef typename CI ::FaceIterator FI;
965
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 ;
970
971 // Allocate cell structures.
972 std::vector<int>(nc, -1).swap(flowSolution_.cellno_);
973
974 std::vector<int>& cell = flowSolution_.cellno_;
975
976 // First pass: enumerate internal faces.
977 int cellno = 0; fpos.push_back(0);
978 int tot_ncf = 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));
982
983 cell[c0] = cellno;
984
985 int& ncf = num_cf[c0];
986
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));
991
992 if (cell[c1] == -1) {
993 // Previously undiscovered internal face.
994 faces.push_back(c1);
995 }
996 }
997 ++ncf;
998 }
999
1000 fpos.push_back(int(faces.size()));
1001 max_ncf_ = std::max(max_ncf_, ncf);
1002 tot_ncf += ncf;
1003 }
1004 assert (cellno == nc);
1005
1006 total_num_faces_ = num_internal_faces_ = int(faces.size());
1007
1008 ip_.init(max_ncf_); ip_.reserveMatrices(num_cf);
1009 F_ .reserve(nc, tot_ncf);
1010
1011 flowSolution_.cellFaces_.reserve(nc, tot_ncf);
1012 flowSolution_.outflux_ .reserve(nc, tot_ncf);
1013
1014 Opm::SparseTable<int>& cf = flowSolution_.cellFaces_;
1015
1016 // Avoid (most) allocation(s) inside 'c' loop.
1017 std::vector<int> l2g; l2g .reserve(max_ncf_);
1018 std::vector<Scalar> F_alloc; F_alloc .reserve(max_ncf_);
1019
1020 // Second pass: build cell-to-face mapping, including boundary.
1021 typedef std::vector<int>::iterator VII;
1022 for (CI c = g.cellbegin(); c != g.cellend(); ++c) {
1023 const int c0 = c->index();
1024
1025 assert ((0 <= c0 ) && ( c0 < nc) &&
1026 (0 <= cell[c0]) && (cell[c0] < nc));
1027
1028 const int ncf = num_cf[cell[c0]];
1029 l2g .resize(ncf , 0 );
1030 F_alloc .resize(ncf , Scalar(0.0));
1031
1032 for (FI f = c->facebegin(); f != c->faceend(); ++f) {
1033 if (f->boundary()) {
1034 // External, not counted before. Add new face...
1035 l2g[f->localIndex()] = total_num_faces_++;
1036 } else {
1037 // Internal face. Need to determine during
1038 // traversal of which cell we discovered this
1039 // face first, and extract the face number
1040 // from the 'faces' table range of that cell.
1041
1042 // Note: std::find() below is potentially
1043 // *VERY* expensive (e.g., large number of
1044 // seeks in moderately sized data in case of
1045 // faulted cells).
1046
1047 const int c1 = f->neighbourCellIndex();
1048 assert ((0 <= c1 ) && ( c1 < nc) &&
1049 (0 <= cell[c1]) && (cell[c1] < nc));
1050
1051 int t = c0, seek = c1;
1052 if (cell[seek] < cell[t])
1053 std::swap(t, seek);
1054
1055 int s = fpos[cell[t]], e = fpos[cell[t] + 1];
1056
1057 VII p = std::find(faces.begin() + s, faces.begin() + e, seek);
1058 assert(p != faces.begin() + e);
1059
1060 l2g[f->localIndex()] = s + (p - (faces.begin() + s));
1061 }
1062 }
1063
1064 cf.appendRow (l2g .begin(), l2g .end());
1065 F_.appendRow (F_alloc.begin(), F_alloc.end());
1066
1067 flowSolution_.outflux_
1068 .appendRow(F_alloc.begin(), F_alloc.end());
1069 }
1070 }
1071
1072
1073 // ----------------------------------------------------------------
1074 void enumerateBCDof(const GridInterface& g, const BCInterface& bc)
1075 // ----------------------------------------------------------------
1076 {
1077 typedef typename GridInterface::CellIterator CI;
1078 typedef typename CI ::FaceIterator FI;
1079
1080 const std::vector<int>& cell = flowSolution_.cellno_;
1081 const Opm::SparseTable<int>& cf = flowSolution_.cellFaces_;
1082
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));
1090 }
1091 }
1092 }
1093
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()];
1101
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];
1106
1107 ppartner_dof_[dof1] = dof2;
1108 ppartner_dof_[dof2] = dof1;
1109 }
1110 }
1111 }
1112 }
1113 }
1114
1115
1116
1117 // ----------------------------------------------------------------
1118 void allocateConnections(const BCInterface& bc)
1119 // ----------------------------------------------------------------
1120 {
1121 // You must call enumerateDof() prior to allocateConnections()
1122 assert(!cleared_state_);
1123
1124 assert (!matrix_structure_valid_);
1125
1126 // Clear any residual data, prepare for assembling structure.
1127 S_.setSize(total_num_faces_, total_num_faces_);
1128 S_.setBuildMode(Dune::BCRSMatrix<MatrixBlockType>::random);
1129
1130 // Compute row sizes
1131 for (int f = 0; f < total_num_faces_; ++f) {
1132 S_.setrowsize(f, 1);
1133 }
1134
1135 allocateGridConnections();
1136 allocateBCConnections(bc);
1137
1138 S_.endrowsizes();
1139
1140 rhs_ .resize(total_num_faces_);
1141 soln_.resize(total_num_faces_);
1142 }
1143
1144
1145 // ----------------------------------------------------------------
1146 void allocateGridConnections()
1147 // ----------------------------------------------------------------
1148 {
1149 const Opm::SparseTable<int>& cf = flowSolution_.cellFaces_;
1150
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);
1155 }
1156 }
1157 }
1158
1159
1160 // ----------------------------------------------------------------
1161 void allocateBCConnections(const BCInterface& bc)
1162 // ----------------------------------------------------------------
1163 {
1164 // Include system connections for periodic boundary
1165 // conditions, if any. We make an arbitrary choice in
1166 // that the face/degree-of-freedom with the minimum
1167 // numerical id of the two periodic partners represents
1168 // the coupling. Suppose <i_p> is this minimum dof-id.
1169 // We then need to introduce a *symmetric* coupling to
1170 // <i_p> to each of the dof's of the cell *NOT* containing
1171 // <i_p>. This choice is implemented in the following
1172 // loop by introducing couplings only when dof1 (self) is
1173 // less than dof2 (other).
1174 //
1175 // See also: setBCConnections() and addCellContrib().
1176 //
1177 typedef typename GridInterface::CellIterator CI;
1178 typedef typename CI ::FaceIterator FI;
1179
1180 const std::vector<int>& cell = flowSolution_.cellno_;
1181 const Opm::SparseTable<int>& cf = flowSolution_.cellFaces_;
1182
1183 if (!bdry_id_map_.empty()) {
1184 // At least one periodic BC. Allocate corresponding
1185 // connections.
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()) {
1189 // dof-id of self
1190 const int dof1 = cf[cell[c->index()]][f->localIndex()];
1191
1192 // dof-id of other
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];
1198
1199 if (dof1 < dof2) {
1200 // Allocate storage for the actual
1201 // couplings.
1202 //
1203 const int ndof = cf.rowSize(c2);
1204 S_.incrementrowsize(dof1, ndof); // self->other
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);
1210 }
1211 S_.incrementrowsize(ii, 1); // other->self
1212 }
1213 }
1214 }
1215 }
1216 }
1217 }
1218 }
1219
1220
1221
1222 // ----------------------------------------------------------------
1223 void setConnections(const BCInterface& bc)
1224 // ----------------------------------------------------------------
1225 {
1226 setGridConnections();
1227 setBCConnections(bc);
1228
1229 S_.endindices();
1230
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_);
1235
1236 matrix_structure_valid_ = true;
1237 }
1238
1239
1240 // ----------------------------------------------------------------
1241 void setGridConnections()
1242 // ----------------------------------------------------------------
1243 {
1244 const Opm::SparseTable<int>& cf = flowSolution_.cellFaces_;
1245
1246 // Compute actual connections (the non-zero structure).
1247 for (int c = 0; c < cf.size(); ++c) {
1248 auto fb = cf[c].begin(), fe = cf[c].end();
1249
1250 for (auto i = fb; i != fe; ++i) {
1251 for (auto j = fb; j != fe; ++j) {
1252 S_.addindex(*i, *j);
1253 }
1254 }
1255 }
1256 }
1257
1258
1259 // ----------------------------------------------------------------
1260 void setBCConnections(const BCInterface& bc)
1261 // ----------------------------------------------------------------
1262 {
1263 // Include system connections for periodic boundary
1264 // conditions, if any. We make an arbitrary choice in
1265 // that the face/degree-of-freedom with the minimum
1266 // numerical id of the two periodic partners represents
1267 // the coupling. Suppose <i_p> is this minimum dof-id.
1268 // We then need to introduce a *symmetric* coupling to
1269 // <i_p> to each of the dof's of the cell *NOT* containing
1270 // <i_p>. This choice is implemented in the following
1271 // loop by introducing couplings only when dof1 (self) is
1272 // less than dof2 (other).
1273 //
1274 // See also: allocateBCConnections() and addCellContrib().
1275 //
1276 typedef typename GridInterface::CellIterator CI;
1277 typedef typename CI ::FaceIterator FI;
1278
1279 const std::vector<int>& cell = flowSolution_.cellno_;
1280 const Opm::SparseTable<int>& cf = flowSolution_.cellFaces_;
1281
1282 if (!bdry_id_map_.empty()) {
1283 // At least one periodic BC. Assign periodic
1284 // connections.
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()) {
1288 // dof-id of self
1289 const int dof1 = cf[cell[c->index()]][f->localIndex()];
1290
1291 // dof-id of other
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];
1297
1298 if (dof1 < dof2) {
1299 // Assign actual couplings.
1300 //
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)) {
1306 ii = pp;
1307 }
1308 S_.addindex(dof1, ii); // self->other
1309 S_.addindex(ii, dof1); // other->self
1310 S_.addindex(dof2, ii);
1311 S_.addindex(ii, dof2);
1312 }
1313 }
1314 }
1315 }
1316 }
1317 }
1318 }
1319
1320
1321
1322 // ----------------------------------------------------------------
1323 template<class FluidInterface>
1324 void assembleDynamic(const FluidInterface& fl ,
1325 const std::vector<double>& sat,
1326 const BCInterface& bc ,
1327 const std::vector<double>& src)
1328 // ----------------------------------------------------------------
1329 {
1330 typedef typename GridInterface::CellIterator CI;
1331
1332 const std::vector<int>& cell = flowSolution_.cellno_;
1333 const Opm::SparseTable<int>& cf = flowSolution_.cellFaces_;
1334
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_);
1339
1340 std::vector<FaceType> facetype (max_ncf_);
1341 std::vector<Scalar> condval (max_ncf_);
1342 std::vector<int> ppartner (max_ncf_);
1343
1344 // Clear residual data
1345 S_ = 0.0;
1346 rhs_ = 0.0;
1347
1348 std::ranges::fill(g_, Scalar(0.0));
1349 std::ranges::fill(L_, Scalar(0.0));
1350 std::ranges::fill(e, Scalar(1.0));
1351
1352 // We will have to regularize resulting system if there
1353 // are no prescribed pressures (i.e., Dirichlet BC's).
1354 do_regularization_ = true;
1355
1356 // Assemble dynamic contributions for each cell
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();
1361
1362 rhs .resize(nf);
1363 gflux.resize(nf);
1364
1365 setExternalContrib(c, c0, bc, src[ci], rhs,
1366 facetype, condval, ppartner);
1367
1368 ip_.computeDynamicParams(c, fl, sat);
1369
1370 SharedFortranMatrix S(nf, nf, &data_store[0]);
1371 ip_.getInverseMatrix(c, S);
1372
1373 std::ranges::fill(gflux, Scalar(0.0));
1374 ip_.gravityFlux(c, gflux);
1375
1376 ImmutableFortranMatrix one(nf, 1, &e[0]);
1377 buildCellContrib(c0, one, gflux, S, rhs);
1378
1379 addCellContrib(S, rhs, facetype, condval, ppartner, cf[c0]);
1380 }
1381 }
1382
1383
1384
1385 // ----------------------------------------------------------------
1386 void solveLinearSystem(double residual_tolerance, int verbosity_level, int maxit)
1387 // ----------------------------------------------------------------
1388 {
1389 // Adapted from DuMux...
1390 Scalar residTol = residual_tolerance;
1391
1392 typedef Dune::BCRSMatrix <MatrixBlockType> MatrixT;
1393 typedef Dune::BlockVector<VectorBlockType> VectorT;
1394 typedef Dune::MatrixAdapter<Matrix,Vector,Vector> Adapter;
1395
1396 // Regularize the matrix (only for pure Neumann problems...)
1397 if (do_regularization_) {
1398 S_[0][0] *= 2;
1399 }
1400 Adapter opS(S_);
1401
1402 // Construct preconditioner.
1403 Dune::SeqILU<MatrixT,VectorT,VectorT> precond(S_, 1.0);
1404
1405 // Construct solver for system of linear equations.
1406 Dune::CGSolver<VectorT> linsolve(opS, precond, residTol,
1407 (maxit>0)?maxit:S_.N(), verbosity_level);
1408
1409 Dune::InverseOperatorResult result;
1410 soln_ = 0.0;
1411
1412 // Solve system of linear equations to recover
1413 // face/contact pressure values (soln_).
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');
1421 }
1422 }
1423
1424
1425
1426 // ------------------ AMG typedefs --------------------
1427
1428 // Representation types for linear system.
1429 typedef Dune::BCRSMatrix <MatrixBlockType> Matrix;
1430 typedef Dune::BlockVector<VectorBlockType> Vector;
1431 typedef Dune::MatrixAdapter<Matrix,Vector,Vector> Operator;
1432
1433 // AMG specific types.
1434 // Old: FIRST_DIAGONAL 1, SYMMETRIC 1, SMOOTHER_ILU 1, ANISOTROPIC_3D 0
1435 // SPE10: FIRST_DIAGONAL 0, SYMMETRIC 1, SMOOTHER_ILU 0, ANISOTROPIC_3D 1
1436#ifndef FIRST_DIAGONAL
1437#define FIRST_DIAGONAL 1
1438#endif
1439#ifndef SYMMETRIC
1440#define SYMMETRIC 1
1441#endif
1442#ifndef SMOOTHER_ILU
1443#define SMOOTHER_ILU 1
1444#endif
1445#ifndef SMOOTHER_BGS
1446#define SMOOTHER_BGS 0
1447#endif
1448#ifndef ANISOTROPIC_3D
1449#define ANISOTROPIC_3D 0
1450#endif
1451
1452#if FIRST_DIAGONAL
1453 typedef Dune::Amg::FirstDiagonal CouplingMetric;
1454#else
1455 typedef Dune::Amg::RowSum CouplingMetric;
1456#endif
1457
1458#if SYMMETRIC
1459 typedef Dune::Amg::SymmetricCriterion<Matrix,CouplingMetric> CriterionBase;
1460#else
1461 typedef Dune::Amg::UnSymmetricCriterion<Matrix,CouplingMetric> CriterionBase;
1462#endif
1463
1464#if SMOOTHER_BGS
1465 typedef Dune::SeqOverlappingSchwarz<Matrix,Vector,Dune::MultiplicativeSchwarzMode> Smoother;
1466#else
1467#if SMOOTHER_ILU
1468 typedef Dune::SeqILU<Matrix,Vector,Vector> Smoother;
1469#else
1470 typedef Dune::SeqSSOR<Matrix,Vector,Vector> Smoother;
1471#endif
1472#endif
1473 typedef Dune::Amg::CoarsenCriterion<CriterionBase> Criterion;
1474
1475
1476 // --------- storing the AMG operator and preconditioner --------
1477 std::unique_ptr<Operator> opS_;
1478 typedef Dune::Preconditioner<Vector,Vector> PrecondBase;
1479 std::unique_ptr<PrecondBase> precond_;
1480
1481
1482 // ----------------------------------------------------------------
1483 void solveLinearSystemAMG(double residual_tolerance, int verbosity_level,
1484 int maxit, double prolong_factor, bool same_matrix, int smooth_steps)
1485 // ----------------------------------------------------------------
1486 {
1487 typedef Dune::Amg::AMG<Operator,Vector,Smoother,Dune::Amg::SequentialInformation>
1488 Precond;
1489
1490 // Adapted from upscaling.cc by Arne Rekdal, 2009
1491 Scalar residTol = residual_tolerance;
1492
1493 if (!same_matrix) {
1494 // Regularize the matrix (only for pure Neumann problems...)
1495 if (do_regularization_) {
1496 S_[0][0] *= 2;
1497 }
1498 opS_.reset(new Operator(S_));
1499
1500 // Construct preconditioner.
1501 double relax = 1;
1502 typename Precond::SmootherArgs smootherArgs;
1503 smootherArgs.relaxationFactor = relax;
1504#if SMOOTHER_BGS
1505 smootherArgs.overlap = Precond::SmootherArgs::none;
1506 smootherArgs.onthefly = false;
1507#endif
1508 Criterion criterion;
1509 criterion.setDebugLevel(verbosity_level);
1510#if ANISOTROPIC_3D
1511 criterion.setDefaultValuesAnisotropic(3, 2);
1512#endif
1513 criterion.setProlongationDampingFactor(prolong_factor);
1514 criterion.setBeta(1e-10);
1515 criterion.setNoPreSmoothSteps(smooth_steps);
1516 criterion.setNoPostSmoothSteps(smooth_steps);
1517 criterion.setGamma(1); // V-cycle; this is the default
1518 precond_.reset(new Precond(*opS_, criterion, smootherArgs));
1519 }
1520 // Construct solver for system of linear equations.
1521 Dune::CGSolver<Vector> linsolve(*opS_, dynamic_cast<Precond&>(*precond_), residTol, (maxit>0)?maxit:S_.N(), verbosity_level);
1522
1523 Dune::InverseOperatorResult result;
1524 soln_ = 0.0;
1525 // Adapt initial guess such Dirichlet boundary conditions are
1526 // represented, i.e. soln_i=A_{ii}^-1 rhs_i
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)
1533 isDirichlet=false;
1534 if(isDirichlet)
1535 soln_[ri.index()]=rhs_[ri.index()]/S_[ri.index()][ri.index()];
1536 }
1537 // Solve system of linear equations to recover
1538 // face/contact pressure values (soln_).
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');
1546 }
1547
1548 }
1549
1550
1551 // ----------------------------------------------------------------
1552 void solveLinearSystemFastAMG(double residual_tolerance, int verbosity_level,
1553 int maxit, double prolong_factor, bool same_matrix, int smooth_steps)
1554 // ----------------------------------------------------------------
1555 {
1556 typedef Dune::Amg::FastAMG<Operator,Vector> Precond;
1557
1558 // Adapted from upscaling.cc by Arne Rekdal, 2009
1559 Scalar residTol = residual_tolerance;
1560
1561 if (!same_matrix) {
1562 // Regularize the matrix (only for pure Neumann problems...)
1563 if (do_regularization_) {
1564 S_[0][0] *= 2;
1565 }
1566 opS_.reset(new Operator(S_));
1567
1568 // Construct preconditioner.
1569 typedef Dune::Amg::AggregationCriterion<Dune::Amg::SymmetricMatrixDependency<Matrix,CouplingMetric> > CritBase;
1570
1571 typedef Dune::Amg::CoarsenCriterion<CritBase> Crit;
1572 Crit criterion;
1573 criterion.setDebugLevel(verbosity_level);
1574#if ANISOTROPIC_3D
1575 criterion.setDefaultValuesAnisotropic(3, 2);
1576#endif
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));
1584 }
1585 // Construct solver for system of linear equations.
1586 Dune::GeneralizedPCGSolver<Vector> linsolve(*opS_, dynamic_cast<Precond&>(*precond_), residTol, (maxit>0)?maxit:S_.N(), verbosity_level);
1587
1588 Dune::InverseOperatorResult result;
1589 soln_ = 0.0;
1590
1591 // Adapt initial guess such Dirichlet boundary conditions are
1592 // represented, i.e. soln_i=A_{ii}^-1 rhs_i
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)
1599 isDirichlet=false;
1600 if(isDirichlet)
1601 soln_[ri.index()]=rhs_[ri.index()]/S_[ri.index()][ri.index()];
1602 }
1603 // Solve system of linear equations to recover
1604 // face/contact pressure values (soln_).
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');
1612 }
1613
1614 }
1615
1616 // ----------------------------------------------------------------
1617 void solveLinearSystemKAMG(double residual_tolerance, int verbosity_level,
1618 int maxit, double prolong_factor, bool same_matrix, int smooth_steps)
1619 // ----------------------------------------------------------------
1620 {
1621
1622 typedef Dune::Amg::KAMG<Operator,Vector,Smoother,Dune::Amg::SequentialInformation> Precond;
1623 // Adapted from upscaling.cc by Arne Rekdal, 2009
1624 Scalar residTol = residual_tolerance;
1625 if (!same_matrix) {
1626 // Regularize the matrix (only for pure Neumann problems...)
1627 if (do_regularization_) {
1628 S_[0][0] *= 2;
1629 }
1630 opS_.reset(new Operator(S_));
1631
1632 // Construct preconditioner.
1633 double relax = 1;
1634 typename Precond::SmootherArgs smootherArgs;
1635 smootherArgs.relaxationFactor = relax;
1636#if SMOOTHER_BGS
1637 smootherArgs.overlap = Precond::SmootherArgs::none;
1638 smootherArgs.onthefly = false;
1639#endif
1640 Criterion criterion;
1641 criterion.setDebugLevel(verbosity_level);
1642#if ANISOTROPIC_3D
1643 criterion.setDefaultValuesAnisotropic(3, 2);
1644#endif
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));
1651 }
1652 // Construct solver for system of linear equations.
1653 Dune::CGSolver<Vector> linsolve(*opS_, dynamic_cast<Precond&>(*precond_), residTol, (maxit>0)?maxit:S_.N(), verbosity_level);
1654
1655 Dune::InverseOperatorResult result;
1656 soln_ = 0.0;
1657 // Adapt initial guess such Dirichlet boundary conditions are
1658 // represented, i.e. soln_i=A_{ii}^-1 rhs_i
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)
1665 isDirichlet=false;
1666 if(isDirichlet)
1667 soln_[ri.index()]=rhs_[ri.index()]/S_[ri.index()][ri.index()];
1668 }
1669 // Solve system of linear equations to recover
1670 // face/contact pressure values (soln_).
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');
1678 }
1679
1680 }
1681
1682
1683
1684 // ----------------------------------------------------------------
1685 template<class FluidInterface>
1686 void computePressureAndFluxes(const FluidInterface& r ,
1687 const std::vector<double>& sat)
1688 // ----------------------------------------------------------------
1689 {
1690 typedef typename GridInterface::CellIterator CI;
1691
1692 const std::vector<int>& cell = flowSolution_.cellno_;
1693 const Opm::SparseTable<int>& cf = flowSolution_.cellFaces_;
1694
1695 std::vector<Scalar>& p = flowSolution_.pressure_;
1696 Opm::SparseTable<Scalar>& v = flowSolution_.outflux_;
1697
1698 //std::vector<double> mob(FluidInterface::NumberOfPhases);
1699 std::vector<double> pi (max_ncf_);
1700 std::vector<double> gflux(max_ncf_);
1701 std::vector<double> Binv_storage(max_ncf_ * max_ncf_);
1702
1703 // Assemble dynamic contributions for each cell
1704 for (CI c = pgrid_->cellbegin(); c != pgrid_->cellend(); ++c) {
1705 const int c0 = cell[c->index()];
1706 const int nf = cf.rowSize(c0);
1707
1708 pi .resize(nf);
1709 gflux.resize(nf);
1710
1711 // Extract contact pressures for cell 'c'.
1712 for (int i = 0; i < nf; ++i) {
1713 pi[i] = soln_[cf[c0][i]];
1714 }
1715
1716 // Compute cell pressure in cell 'c'.
1717 p[c0] = (g_[c0] +
1718 std::inner_product(F_[c0].begin(), F_[c0].end(),
1719 pi.begin(), 0.0)) / L_[c0];
1720
1721 std::ranges::transform(pi, pi.begin(),
1722 [&p, c0](const double& input) { return p[c0] - input; });
1723
1724 // Recover fluxes from local system
1725 // Bv = Bv_g + Cp - D\pi
1726 //
1727 // 1) Solve system Bv = Cp - D\pi
1728 //
1729 ip_.computeDynamicParams(c, r, sat);
1730
1731 SharedFortranMatrix Binv(nf, nf, &Binv_storage[0]);
1732 ip_.getInverseMatrix(c, Binv);
1733 vecMulAdd_N(Scalar(1.0), Binv, &pi[0], Scalar(0.0), &v[c0][0]);
1734
1735 // 2) Add gravity flux contributions (v <- v + v_g)
1736 //
1737 ip_.gravityFlux(c, gflux);
1738 std::ranges::transform(gflux, v[c0], v[c0].begin(), std::plus<Scalar>());
1739 }
1740 }
1741
1742
1743
1744
1745 // ----------------------------------------------------------------
1746 void setExternalContrib(const typename GridInterface::CellIterator c,
1747 const int c0, const BCInterface& bc,
1748 const double src,
1749 std::vector<Scalar>& rhs,
1750 std::vector<FaceType>& facetype,
1751 std::vector<double>& condval,
1752 std::vector<int>& ppartner)
1753 // ----------------------------------------------------------------
1754 {
1755 typedef typename GridInterface::CellIterator::FaceIterator FI;
1756
1757 const Opm::SparseTable<int>& cf = flowSolution_.cellFaces_;
1758
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);
1763
1764 g_[c0] = src;
1765
1766 int k = 0;
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());
1778
1779 facetype[k] = Periodic;
1780 condval[k] = bcond.pressureDifference();
1781 ppartner[k] = cf[j->second.first][j->second.second];
1782 } else {
1783 assert (bcond.isNeumann());
1784 facetype[k] = Neumann;
1785 rhs[k] = bcond.outflux();
1786 }
1787 }
1788 }
1789 }
1790
1791
1792
1793
1794 // ----------------------------------------------------------------
1795 void buildCellContrib(const int c ,
1796 const ImmutableFortranMatrix& one ,
1797 const std::vector<Scalar>& gflux,
1799 std::vector<Scalar>& rhs)
1800 // ----------------------------------------------------------------
1801 {
1802 // Ft <- B^{-t} * ones([size(S,2),1])
1803 SharedFortranMatrix Ft(S.numRows(), 1, &F_[c][0]);
1804 matMulAdd_TN(Scalar(1.0), S, one, Scalar(0.0), Ft);
1805
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));
1808
1809 // rhs <- v_g - rhs (== v_g - h)
1810 std::ranges::transform(gflux, rhs, rhs.begin(), std::minus<Scalar>());
1811
1812 // rhs <- rhs + g_[c]/L_[c]*F
1813 std::transform(rhs.begin(), rhs.end(), Ft.data(),
1814 rhs.begin(),
1815 axpby<Scalar>(Scalar(1.0), Scalar(g_[c] / L_[c])));
1816
1817 // S <- S - F'*F/L_c
1818 symmetricUpdate(-Scalar(1.0)/L_[c], Ft, Scalar(1.0), S);
1819 }
1820
1821
1822
1823 // ----------------------------------------------------------------
1825 template<class L2G>
1826 void addCellContrib(const SharedFortranMatrix& S ,
1827 const std::vector<Scalar>& rhs ,
1828 const std::vector<FaceType>& facetype,
1829 const std::vector<Scalar>& condval ,
1830 const std::vector<int>& ppartner,
1831 const L2G& l2g)
1832 // ----------------------------------------------------------------
1833 {
1834 int r = 0;
1835 for (auto i = l2g.begin(); i != l2g.end(); ++i, ++r) {
1836 // Indirection for periodic BC handling.
1837 int ii = *i;
1838
1839 switch (facetype[r]) {
1840 case Dirichlet:
1841 // Pressure is already known. Assemble trivial
1842 // equation of the form: a*x = a*p where 'p' is
1843 // the known pressure value (i.e., condval[r]).
1844 //
1845 S_ [ii][ii] = S(r,r);
1846 rhs_[ii] = S(r,r) * condval[r];
1847 continue;
1848 case Periodic:
1849 // Periodic boundary condition. Contact pressures
1850 // linked by relations of the form
1851 //
1852 // a*(x0 - x1) = a*pd
1853 //
1854 // where 'pd' is the known pressure difference
1855 // x0-x1 (condval[r]). Preserve matrix symmetry
1856 // by assembling both of the equations
1857 //
1858 // a*(x0-x1) = a* pd, and (1)
1859 // a*(x1-x0) = a*(-pd) (2)
1860 //
1861 assert ((0 <= ppartner[r]) && (ppartner[r] < int(rhs_.size())));
1862 assert (ii != ppartner[r]);
1863
1864 {
1865 const double a = S(r,r), b = a * condval[r];
1866
1867 // Equation (1)
1868 S_ [ ii][ ii] += a;
1869 S_ [ ii][ppartner[r]] -= a;
1870 rhs_[ ii] += b;
1871
1872 // Equation (2)
1873 S_ [ppartner[r]][ ii] -= a;
1874 S_ [ppartner[r]][ppartner[r]] += a;
1875 rhs_[ppartner[r]] -= b;
1876 }
1877
1878 ii = std::min(ii, ppartner[r]);
1879
1880 // INTENTIONAL FALL-THROUGH!
1881 // IOW: Don't insert <break;> here!
1882 //
1883 default:
1884 int c = 0;
1885 for (auto j = l2g.begin(); j != l2g.end(); ++j, ++c) {
1886 // Indirection for periodic BC handling.
1887 int jj = *j;
1888
1889 if (facetype[c] == Dirichlet) {
1890 rhs_[ii] -= S(r,c) * condval[c];
1891 continue;
1892 }
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];
1898 jj = ppartner[c];
1899 }
1900 }
1901 S_[ii][jj] += S(r,c);
1902 }
1903 break;
1904 }
1905 rhs_[ii] += rhs[r];
1906 }
1907 }
1908 };
1909} // namespace Opm
1910
1911#endif // OPENRS_INCOMPFLOWSOLVERHYBRID_HEADER
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