25#ifndef TPSA_LINEARIZER_HPP
26#define TPSA_LINEARIZER_HPP
28#include <dune/common/fvector.hh>
30#include <opm/common/TimingMacros.hpp>
32#include <opm/grid/utility/SparseTable.hpp>
34#include <opm/input/eclipse/Schedule/BCProp.hpp>
36#include <opm/material/materialstates/MaterialStateTPSA.hpp>
51template<
class TypeTag>
67 using MaterialState = MaterialStateTPSA<Evaluation>;
69 enum { numEq = getPropValue<TypeTag, Properties::NumEqTPSA>() };
71 using ADVectorBlock = Dune::FieldVector<Evaluation, numEq>;
72 using MatrixBlock =
typename SparseMatrixAdapter::MatrixBlock;
73 using VectorBlock = Dune::FieldVector<Scalar, numEq>;
84 simulatorPtr_ =
nullptr;
95 void init(Simulator& simulator)
97 simulatorPtr_ = &simulator;
117 jacobian_->finalize();
125 template <
class SubDomainType>
129 initFirstIteration_();
131 for (
int globI : domain.cells) {
132 residual_[globI] = 0.0;
133 jacobian_->clearRow(globI, 0.0);
159 catch (
const std::exception& e) {
160 std::cout <<
"rank " << simulator_().gridView().comm().rank()
161 <<
" caught an exception while linearizing TPSA system:" << e.what()
162 <<
"\n" << std::flush;
166 std::cout <<
"rank " << simulator_().gridView().comm().rank()
167 <<
" caught an exception while linearizing TPSA system"
168 <<
"\n" << std::flush;
171 succeeded = simulator_().gridView().comm().min(succeeded);
174 throw NumericalProblem(
"A process did not succeed in linearizing the TPSA system");
186 template <
class SubDomainType>
194 initFirstIteration_();
197 if (domain.cells.size() == flowModel_().numTotalDof()) {
221 void setResAndJacobi(VectorBlock& res, MatrixBlock& bMat,
const ADVectorBlock& resid)
const
224 for (
unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
225 res[eqIdx] = resid[eqIdx].value();
230 for (
unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
231 for (
unsigned pvIdx = 0; pvIdx < numEq; ++pvIdx) {
232 bMat[eqIdx][pvIdx] = resid[eqIdx].derivative(pvIdx);
242 for (
auto& bdyInfo : boundaryInfo_) {
244 const auto [type, displacementAD] = problem_().mechBoundaryCondition(bdyInfo.cell, bdyInfo.dir);
247 std::vector<double> displacement(3, 0.0);
248 for (std::size_t ii = 0; ii < displacement.size(); ++ii) {
249 displacement[ii] = displacementAD[ii].value();
253 bdyInfo.bcdata.type = type;
254 bdyInfo.bcdata.displacement = displacement;
267 {
return *jacobian_; }
275 {
return *jacobian_; }
283 {
return residual_; }
291 {
return residual_; }
302 {
return linearizationType_; }
320 { linearizationType_ = linearizationType; }
334 OPM_TIMEBLOCK(createMatrixTPSA);
337 if (!neighborInfo_.empty()) {
342 const auto& flowModel = flowModel_();
343 Stencil stencil(gridView_(), flowModel.dofMapper());
346 std::vector<std::set<unsigned>> sparsityPattern(flowModel.numTotalDof());
347 unsigned numCells = flowModel.numTotalDof();
348 neighborInfo_.reserve(numCells, 6 * numCells);
349 std::vector<NeighborInfo> loc_nbinfo;
350 for (
const auto& elem : elements(gridView_())) {
352 stencil.update(elem);
354 for (
unsigned primaryDofIdx = 0; primaryDofIdx < stencil.numPrimaryDof(); ++primaryDofIdx) {
356 const unsigned myIdx = stencil.globalSpaceIndex(primaryDofIdx);
357 loc_nbinfo.resize(stencil.numDof() - 1);
359 for (
unsigned dofIdx = 0; dofIdx < stencil.numDof(); ++dofIdx) {
362 const unsigned neighborIdx = stencil.globalSpaceIndex(dofIdx);
363 sparsityPattern[myIdx].insert(neighborIdx);
365 const auto scvfIdx = dofIdx - 1;
366 const auto& scvf = stencil.interiorFace(scvfIdx);
367 const Scalar area = scvf.area();
368 loc_nbinfo[dofIdx - 1] = NeighborInfo{ neighborIdx, area,
nullptr };
373 neighborInfo_.appendRow(loc_nbinfo.begin(), loc_nbinfo.end());
376 unsigned bfIndex = 0;
377 for (
const auto& intersection : intersections(gridView_(), elem)) {
378 if (intersection.boundary()) {
380 const auto& bf = stencil.boundaryFace(bfIndex);
381 const int dir_id = bf.dirId();
389 const auto [type, displacementAD] = problem_().mechBoundaryCondition(myIdx, dir_id);
392 std::vector<double> displacement(3, 0.0);
393 for (std::size_t ii = 0; ii < displacement.size(); ++ii) {
394 displacement[ii] = displacementAD[ii].value();
398 BoundaryConditionData bcdata { type, displacement, bfIndex, bf.area() };
399 boundaryInfo_.push_back( { myIdx, dir_id, bfIndex, bcdata } );
403 if (!intersection.neighbor()) {
412 jacobian_ = std::make_unique<SparseMatrixAdapter>(simulator_());
413 diagMatAddress_.resize(numCells);
414 jacobian_->reserve(sparsityPattern);
415 for (
unsigned globI = 0; globI < numCells; globI++) {
416 const auto& nbInfos = neighborInfo_[globI];
417 diagMatAddress_[globI] = jacobian_->blockAddress(globI, globI);
418 for (
auto& nbInfo : nbInfos) {
419 nbInfo.matBlockAddress = jacobian_->blockAddress(nbInfo.neighbor, globI);
424 fullDomain_.cells.resize(numCells);
425 std::iota(fullDomain_.cells.begin(), fullDomain_.cells.end(), 0);
433 template <
class SubDomainType>
434 void linearize_(
const SubDomainType& domain)
436 OPM_TIMEBLOCK(linearizeTPSA);
439 const auto& flowModel = flowModel_();
440 const auto& geoMechModel = geoMechModel_();
441 auto& problem = problem_();
442 const unsigned int numCells = domain.cells.size();
445#pragma omp parallel for
448 for (
unsigned ii = 0; ii < numCells; ++ii) {
449 OPM_TIMEBLOCK_LOCAL(linearizationForEachCellTPSA, Subsystem::Assembly);
451 const unsigned globI = domain.cells[ii];
452 const auto& nbInfos = neighborInfo_[globI];
453 VectorBlock res(0.0);
454 MatrixBlock bMat(0.0);
455 ADVectorBlock adres(0.0);
456 const MaterialState& materialStateIn = geoMechModel.materialState(globI, 0);
462 OPM_TIMEBLOCK_LOCAL(faceCalculationForEachCellTPSA, Subsystem::Assembly);
465 for (
const auto& nbInfo : nbInfos) {
466 OPM_TIMEBLOCK_LOCAL(calculationForEachFaceTPSA, Subsystem::Assembly);
474 const unsigned globJ = nbInfo.neighbor;
475 assert(globJ != globI);
476 const MaterialState& materialStateEx = geoMechModel.materialState(globJ, 0);
479 LocalResidual::computeFaceTerm(adres,
485 adres *= nbInfo.faceArea;
491 residual_[globI] += res;
495 *diagMatAddress_[globI] += bMat;
503 *nbInfo.matBlockAddress += bMat;
512 OPM_TIMEBLOCK_LOCAL(computeVolumeTerm, Subsystem::Assembly);
515 LocalResidual::computeVolumeTerm(adres,
520 const double volume = flowModel.dofTotalVolume(globI);
527 residual_[globI] += res;
531 *diagMatAddress_[globI] += bMat;
541 LocalResidual::computeSourceTerm(adres,
551 residual_[globI] += res;
555 *diagMatAddress_[globI] += bMat;
561 for (
const auto& bdyInfo : boundaryInfo_) {
562 VectorBlock res(0.0);
563 MatrixBlock bMat(0.0);
564 ADVectorBlock adres(0.0);
565 const unsigned globI = bdyInfo.cell;
566 const MaterialState& materialStateIn = geoMechModel.materialState(globI, 0);
569 LocalResidual::computeBoundaryTerm(adres,
574 adres *= bdyInfo.bcdata.faceArea;
580 residual_[globI] += res;
584 *diagMatAddress_[globI] += bMat;
591 void initFirstIteration_()
597 residual_.resize(flowModel_().numTotalDof());
621 Simulator& simulator_()
622 {
return *simulatorPtr_; }
629 const Simulator& simulator_()
const
630 {
return *simulatorPtr_; }
638 {
return simulator_().problem(); }
645 const Problem& problem_()
const
646 {
return simulator_().problem(); }
653 FlowModel& flowModel_()
654 {
return simulator_().model(); }
661 const FlowModel& flowModel_()
const
662 {
return simulator_().model(); }
669 GeomechModel& geoMechModel_()
670 {
return problem_().geoMechModel(); }
677 const GeomechModel& geoMechModel_()
const
678 {
return problem_().geoMechModel(); }
685 const GridView& gridView_()
const
686 {
return problem_().gridView(); }
688 Simulator* simulatorPtr_{};
689 LinearizationType linearizationType_{};
691 std::vector<MatrixBlock*> diagMatAddress_{};
692 std::unique_ptr<SparseMatrixAdapter> jacobian_{};
693 GlobalEqVector residual_;
703 unsigned int neighbor;
705 MatrixBlock* matBlockAddress;
707 SparseTable<NeighborInfo> neighborInfo_{};
712 struct BoundaryConditionData
715 std::vector<double> displacement;
716 unsigned boundaryFaceIndex;
727 unsigned int bfIndex;
728 BoundaryConditionData bcdata;
730 std::vector<BoundaryInfo> boundaryInfo_;
737 std::vector<int> cells;
738 std::vector<bool> interior;
740 FullDomain fullDomain_;
Linearizes TPSA equations and generates system matrix and residual for linear solver.
Definition: tpsalinearizer.hpp:53
GlobalEqVector & residual()
Get residual vector.
Definition: tpsalinearizer.hpp:290
const SparseMatrixAdapter & jacobian() const
Get Jacobian matrix.
Definition: tpsalinearizer.hpp:266
const GlobalEqVector & residual() const
Get residual vector.
Definition: tpsalinearizer.hpp:282
void eraseMatrix()
Causes the Jacobian matrix to be recreated from scratch before the next iteration.
Definition: tpsalinearizer.hpp:107
void resetSystem_(const SubDomainType &domain)
Initializing and/or reset residual and Jacobian.
Definition: tpsalinearizer.hpp:126
std::map< unsigned, Constraints > constraintsMap() const
Get constraints map.
Definition: tpsalinearizer.hpp:311
void updateBoundaryConditionData()
Update boundary condition information.
Definition: tpsalinearizer.hpp:240
SparseMatrixAdapter & jacobian()
Get Jacobian matrix.
Definition: tpsalinearizer.hpp:274
TpsaLinearizer()
Constructor.
Definition: tpsalinearizer.hpp:82
void linearizeDomain(const SubDomainType &domain)
Linearize the non-linear system for the spatial domain.
Definition: tpsalinearizer.hpp:187
void linearizeDomain()
Linearize the non-linear system for the spatial domain.
Definition: tpsalinearizer.hpp:152
void finalize()
Finalize creation of Jacobian matrix and make ready for linear solver.
Definition: tpsalinearizer.hpp:115
void setLinearizationType(LinearizationType linearizationType)
Set linearization type.
Definition: tpsalinearizer.hpp:319
const LinearizationType & getLinearizationType() const
Get linearization type.
Definition: tpsalinearizer.hpp:301
void linearizeAuxiliaryEquations()
Linearize auxillary equation.
Definition: tpsalinearizer.hpp:211
void init(Simulator &simulator)
Initialize the linearizer.
Definition: tpsalinearizer.hpp:95
void linearize()
Linearize the non-linear system.
Definition: tpsalinearizer.hpp:140
void setResAndJacobi(VectorBlock &res, MatrixBlock &bMat, const ADVectorBlock &resid) const
Extract local residuals and sub-block Jacobians from locally computed AD residual.
Definition: tpsalinearizer.hpp:221
Definition: blackoilbioeffectsmodules.hh:45
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
Definition: linearizationtype.hh:34