tpsalinearizer.hpp
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3/*
4 Copyright 2025 NORCE AS
5
6 This file is part of the Open Porous Media project (OPM).
7
8 OPM is free software: you can redistribute it and/or modify
9 it under the terms of the GNU General Public License as published by
10 the Free Software Foundation, either version 2 of the License, or
11 (at your option) any later version.
12
13 OPM is distributed in the hope that it will be useful,
14 but WITHOUT ANY WARRANTY; without even the implied warranty of
15 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
16 GNU General Public License for more details.
17
18 You should have received a copy of the GNU General Public License
19 along with OPM. If not, see <http://www.gnu.org/licenses/>.
20
21 Consult the COPYING file in the top-level source directory of this
22 module for the precise wording of the license and the list of
23 copyright holders.
24*/
25#ifndef TPSA_LINEARIZER_HPP
26#define TPSA_LINEARIZER_HPP
27
28#include <dune/common/fvector.hh>
29
30#include <opm/common/TimingMacros.hpp>
31
32#include <opm/grid/utility/SparseTable.hpp>
33
34#include <opm/input/eclipse/Schedule/BCProp.hpp>
35
36#include <opm/material/materialstates/MaterialStateTPSA.hpp>
37
40
41#include <cassert>
42#include <set>
43#include <vector>
44
45
46namespace Opm {
47
51template<class TypeTag>
53{
54 using Constraints = GetPropType<TypeTag, Properties::Constraints>; // TODO: make TPSA constraints
66
67 using MaterialState = MaterialStateTPSA<Evaluation>;
68
69 enum { numEq = getPropValue<TypeTag, Properties::NumEqTPSA>() };
70
71 using ADVectorBlock = Dune::FieldVector<Evaluation, numEq>;
72 using MatrixBlock = typename SparseMatrixAdapter::MatrixBlock;
73 using VectorBlock = Dune::FieldVector<Scalar, numEq>;
74
75public:
76 // ///
77 // Public functions
78 // ///
83 {
84 simulatorPtr_ = nullptr;
85 }
86
95 void init(Simulator& simulator)
96 {
97 simulatorPtr_ = &simulator;
99 }
100
108 {
109 jacobian_.reset();
110 }
111
115 void finalize()
116 {
117 jacobian_->finalize();
118 }
119
125 template <class SubDomainType>
126 void resetSystem_(const SubDomainType& domain)
127 {
128 if (!jacobian_) {
129 initFirstIteration_();
130 }
131 for (int globI : domain.cells) {
132 residual_[globI] = 0.0;
133 jacobian_->clearRow(globI, 0.0);
134 }
135 }
136
141 {
144 }
145
153 {
154 int succeeded;
155 try {
156 linearizeDomain(fullDomain_);
157 succeeded = 1;
158 }
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;
163 succeeded = 0;
164 }
165 catch (...) {
166 std::cout << "rank " << simulator_().gridView().comm().rank()
167 << " caught an exception while linearizing TPSA system"
168 << "\n" << std::flush;
169 succeeded = 0;
170 }
171 succeeded = simulator_().gridView().comm().min(succeeded);
172
173 if (!succeeded) {
174 throw NumericalProblem("A process did not succeed in linearizing the TPSA system");
175 }
176 }
177
186 template <class SubDomainType>
187 void linearizeDomain(const SubDomainType& domain)
188 {
189 OPM_TIMEBLOCK(linearizeDomain);
190
191 // We defer the initialization of the Jacobian matrix until here because the auxiliary modules usually assume
192 // the problem, model and grid to be fully initialized'
193 if (!jacobian_) {
194 initFirstIteration_();
195 }
196
197 if (domain.cells.size() == flowModel_().numTotalDof()) {
198 // We are on the full domain.
199 resetSystem_();
200 }
201 else {
202 resetSystem_(domain);
203 }
204
205 linearize_(domain);
206 }
207
212 { }
213
221 void setResAndJacobi(VectorBlock& res, MatrixBlock& bMat, const ADVectorBlock& resid) const
222 {
223 // Scalar local residual
224 for (unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
225 res[eqIdx] = resid[eqIdx].value();
226 }
227
228 // A[dofIdx][focusDofIdx][eqIdx][pvIdx] is the partial derivative of the residual function 'eqIdx' for the
229 // degree of freedom 'dofIdx' with regard to the focus variable 'pvIdx' of the degree of freedom 'focusDofIdx'
230 for (unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
231 for (unsigned pvIdx = 0; pvIdx < numEq; ++pvIdx) {
232 bMat[eqIdx][pvIdx] = resid[eqIdx].derivative(pvIdx);
233 }
234 }
235 }
236
241 {
242 for (auto& bdyInfo : boundaryInfo_) {
243 // Get boundary information from problem
244 const auto [type, displacementAD] = problem_().mechBoundaryCondition(bdyInfo.cell, bdyInfo.dir);
245
246 // Strip the unnecessary (and zero anyway) derivatives off displacement
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();
250 }
251
252 // Update boundary information
253 bdyInfo.bcdata.type = type;
254 bdyInfo.bcdata.displacement = displacement;
255 }
256 }
257
258 // ///
259 // Public get and set functions
260 // ///
266 const SparseMatrixAdapter& jacobian() const
267 { return *jacobian_; }
268
274 SparseMatrixAdapter& jacobian()
275 { return *jacobian_; }
276
282 const GlobalEqVector& residual() const
283 { return residual_; }
284
290 GlobalEqVector& residual()
291 { return residual_; }
292
302 { return linearizationType_; }
303
311 std::map<unsigned, Constraints> constraintsMap() const
312 { return {}; }
313
320 { linearizationType_ = linearizationType; }
321
322private:
323 // ///
324 // Private functions
325 // ///
332 void createMatrix_()
333 {
334 OPM_TIMEBLOCK(createMatrixTPSA);
335
336 // If the Jacobian has been initialize before, we jump out
337 if (!neighborInfo_.empty()) {
338 return;
339 }
340
341 // Init. the stencil
342 const auto& flowModel = flowModel_();
343 Stencil stencil(gridView_(), flowModel.dofMapper());
344
345 // Build up sparsity patterns and neighboring information for Jacobian and linearization
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_())) {
351 // Loop over primary dofs in the element
352 stencil.update(elem);
353
354 for (unsigned primaryDofIdx = 0; primaryDofIdx < stencil.numPrimaryDof(); ++primaryDofIdx) {
355 // Build up neighboring information for curret primary dof
356 const unsigned myIdx = stencil.globalSpaceIndex(primaryDofIdx);
357 loc_nbinfo.resize(stencil.numDof() - 1);
358
359 for (unsigned dofIdx = 0; dofIdx < stencil.numDof(); ++dofIdx) {
360 // NOTE: NeighborInfo could/should be expanded with cell face parameters located in problem_()
361 // needed when computing face terms in LocalResidual
362 const unsigned neighborIdx = stencil.globalSpaceIndex(dofIdx);
363 sparsityPattern[myIdx].insert(neighborIdx);
364 if (dofIdx > 0) {
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 };
369 }
370 }
371
372 // Insert local neighbor info
373 neighborInfo_.appendRow(loc_nbinfo.begin(), loc_nbinfo.end());
374
375 // Boundary condition information
376 unsigned bfIndex = 0;
377 for (const auto& intersection : intersections(gridView_(), elem)) {
378 if (intersection.boundary()) {
379 // Get boundary face direction
380 const auto& bf = stencil.boundaryFace(bfIndex);
381 const int dir_id = bf.dirId();
382
383 // Skip NNCs
384 if (dir_id < 0) {
385 continue;
386 }
387
388 // Get boundary information from problem()
389 const auto [type, displacementAD] = problem_().mechBoundaryCondition(myIdx, dir_id);
390
391 // Strip the unnecessary (and zero anyway) derivatives off displacement
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();
395 }
396
397 // Insert boundary condition data in container
398 BoundaryConditionData bcdata { type, displacement, bfIndex, bf.area() };
399 boundaryInfo_.push_back( { myIdx, dir_id, bfIndex, bcdata } );
400 ++bfIndex;
401 continue;
402 }
403 if (!intersection.neighbor()) {
404 ++bfIndex;
405 continue;
406 }
407 }
408 }
409 }
410
411 // Allocate Jacobian matrix and pointers to its sub-blocks
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);
420 }
421 }
422
423 // Create full domain
424 fullDomain_.cells.resize(numCells);
425 std::iota(fullDomain_.cells.begin(), fullDomain_.cells.end(), 0);
426 }
427
433 template <class SubDomainType>
434 void linearize_(const SubDomainType& domain)
435 {
436 OPM_TIMEBLOCK(linearizeTPSA);
437
438 // Extract misc. variables used in linearization
439 const auto& flowModel = flowModel_();
440 const auto& geoMechModel = geoMechModel_();
441 auto& problem = problem_();
442 const unsigned int numCells = domain.cells.size();
443
444#ifdef _OPENMP
445#pragma omp parallel for
446#endif
447 // Loop over cells in the domain and compute local residual and jacobian
448 for (unsigned ii = 0; ii < numCells; ++ii) {
449 OPM_TIMEBLOCK_LOCAL(linearizationForEachCellTPSA, Subsystem::Assembly);
450
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, /*timeIdx=*/0);
457
458 // ///
459 // Face term
460 // ///
461 {
462 OPM_TIMEBLOCK_LOCAL(faceCalculationForEachCellTPSA, Subsystem::Assembly);
463
464 // Loop over neighboring cells
465 for (const auto& nbInfo : nbInfos) {
466 OPM_TIMEBLOCK_LOCAL(calculationForEachFaceTPSA, Subsystem::Assembly);
467
468 // Reset local residual and Jacobian
469 res = 0.0;
470 bMat = 0.0;
471 adres = 0.0;
472
473 // Neighbor information
474 const unsigned globJ = nbInfo.neighbor;
475 assert(globJ != globI);
476 const MaterialState& materialStateEx = geoMechModel.materialState(globJ, /*timeIdx=*/0);
477
478 // Compute local face term
479 LocalResidual::computeFaceTerm(adres,
480 materialStateIn,
481 materialStateEx,
482 problem,
483 globI,
484 globJ);
485 adres *= nbInfo.faceArea;
486
487 // Extract residual and sub-block Jacobian entries from computed AD residual
488 setResAndJacobi(res, bMat, adres);
489
490 // Insert into global residual
491 residual_[globI] += res;
492
493 // Insert contribution to (globI, globI) sub-block
494 // SparseAdapter syntax: jacobian_->addToBlock(globI, globI, bMat);
495 *diagMatAddress_[globI] += bMat;
496
497 // Insert contribution to (globJ, globI) sub-block
498 // Note: since LocalResidual::computeFaceTerm have Evaluation on globI primary variables, it is
499 // natural to insert Jacobian entries for (globJ, globI) here, since we only need to flip the signs
500 // in the calculated face terms
501 // SparseAdapter syntax: jacobian_->addToBlock(globJ, globI, bMat);
502 bMat *= -1.0;
503 *nbInfo.matBlockAddress += bMat;
504 }
505 }
506
507 // ///
508 // Volume term
509 // //
510 adres = 0.0;
511 {
512 OPM_TIMEBLOCK_LOCAL(computeVolumeTerm, Subsystem::Assembly);
513
514 // Compute local volume term
515 LocalResidual::computeVolumeTerm(adres,
516 materialStateIn,
517 problem,
518 globI);
519 }
520 const double volume = flowModel.dofTotalVolume(globI);
521 adres *= volume;
522
523 // Extract residual and sub-block Jacobian entries from computed AD residual
524 setResAndJacobi(res, bMat, adres);
525
526 // Insert in global residual
527 residual_[globI] += res;
528
529 // Insert contribution to (globI, globI) sub-block
530 // SparseAdapter syntax: jacobian_->addToBlock(globI, globI, bMat);
531 *diagMatAddress_[globI] += bMat;
532
533 // ///
534 // Source term
535 // ///
536 res = 0.0;
537 bMat = 0.0;
538 adres = 0.0;
539
540 // Compute local source term
541 LocalResidual::computeSourceTerm(adres,
542 problem,
543 globI,
544 0);
545 adres *= -volume;
546
547 // Extract residual and sub-block Jacobian entries from computed AD residual
548 setResAndJacobi(res, bMat, adres);
549
550 // Insert into global residual
551 residual_[globI] += res;
552
553 // Insert contribution to (globI, globI) sub-block
554 // SparseAdapter syntax: jacobian_->addToBlock(globI, globI, bMat);
555 *diagMatAddress_[globI] += bMat;
556 } // globI loop
557
558 // ///
559 // Boundary term
560 // ///
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, /*timeIdx=*/0);
567
568 // Compute local boundary condition
569 LocalResidual::computeBoundaryTerm(adres,
570 materialStateIn,
571 bdyInfo.bcdata,
572 problem,
573 globI);
574 adres *= bdyInfo.bcdata.faceArea;
575
576 // Extract residual and sub-block Jacobian entries from computed AD residual
577 setResAndJacobi(res, bMat, adres);
578
579 // Insert in global residual
580 residual_[globI] += res;
581
582 // Insert contribution to (globI, globI) sub-block
583 // SparseAdapter syntax: jacobian_->addToBlock(globI, globI, bMat);
584 *diagMatAddress_[globI] += bMat;
585 }
586 }
587
591 void initFirstIteration_()
592 {
593 // initialize the BCRS matrix for the Jacobian of the residual function
594 createMatrix_();
595
596 // initialize the Jacobian matrix and the vector for the residual function
597 residual_.resize(flowModel_().numTotalDof());
598 resetSystem_();
599 }
600
604 void resetSystem_()
605 {
606 // Set residual vector entries to zero
607 residual_ = 0.0;
608
609 // Set all Jacobian matrix entries to zero
610 jacobian_->clear();
611 }
612
613 // ///
614 // Private get functions
615 // ///
621 Simulator& simulator_()
622 { return *simulatorPtr_; }
623
629 const Simulator& simulator_() const
630 { return *simulatorPtr_; }
631
637 Problem& problem_()
638 { return simulator_().problem(); }
639
645 const Problem& problem_() const
646 { return simulator_().problem(); }
647
653 FlowModel& flowModel_()
654 { return simulator_().model(); }
655
661 const FlowModel& flowModel_() const
662 { return simulator_().model(); }
663
669 GeomechModel& geoMechModel_()
670 { return problem_().geoMechModel(); }
671
677 const GeomechModel& geoMechModel_() const
678 { return problem_().geoMechModel(); }
679
685 const GridView& gridView_() const
686 { return problem_().gridView(); }
687
688 Simulator* simulatorPtr_{};
689 LinearizationType linearizationType_{};
690
691 std::vector<MatrixBlock*> diagMatAddress_{};
692 std::unique_ptr<SparseMatrixAdapter> jacobian_{};
693 GlobalEqVector residual_;
694
695 //
696 // Helper structs
697 //
701 struct NeighborInfo
702 {
703 unsigned int neighbor;
704 double faceArea;
705 MatrixBlock* matBlockAddress;
706 };
707 SparseTable<NeighborInfo> neighborInfo_{};
708
712 struct BoundaryConditionData
713 {
714 BCMECHType type;
715 std::vector<double> displacement;
716 unsigned boundaryFaceIndex;
717 double faceArea;
718 };
719
723 struct BoundaryInfo
724 {
725 unsigned int cell;
726 int dir;
727 unsigned int bfIndex;
728 BoundaryConditionData bcdata;
729 };
730 std::vector<BoundaryInfo> boundaryInfo_;
731
735 struct FullDomain
736 {
737 std::vector<int> cells;
738 std::vector<bool> interior;
739 };
740 FullDomain fullDomain_;
741}; // class TpsaLinearizer
742
743} // namespace Opm
744
745#endif
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