tpfalinearizer.hh
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 This file is part of the Open Porous Media project (OPM).
5
6 OPM is free software: you can redistribute it and/or modify
7 it under the terms of the GNU General Public License as published by
8 the Free Software Foundation, either version 2 of the License, or
9 (at your option) any later version.
10
11 OPM is distributed in the hope that it will be useful,
12 but WITHOUT ANY WARRANTY; without even the implied warranty of
13 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14 GNU General Public License for more details.
15
16 You should have received a copy of the GNU General Public License
17 along with OPM. If not, see <http://www.gnu.org/licenses/>.
18
19 Consult the COPYING file in the top-level source directory of this
20 module for the precise wording of the license and the list of
21 copyright holders.
22*/
28#ifndef TPFA_LINEARIZER_HH
29#define TPFA_LINEARIZER_HH
30
31#include <dune/common/version.hh>
32#include <dune/common/fvector.hh>
33#include <dune/common/fmatrix.hh>
34
35#include <opm/common/Exceptions.hpp>
36#include <opm/common/TimingMacros.hpp>
37
38#include <opm/grid/utility/SparseTable.hpp>
39
40#include <opm/material/common/ConditionalStorage.hpp>
41
42#include <opm/input/eclipse/EclipseState/Grid/FaceDir.hpp>
43#include <opm/input/eclipse/Schedule/BCProp.hpp>
44
51
52#include <cassert>
53#include <cstddef>
54#include <exception> // current_exception, rethrow_exception
55#include <iostream>
56#include <map>
57#include <memory>
58#include <numeric>
59#include <set>
60#include <stdexcept>
61#include <unordered_map>
62#include <vector>
63
64#include <fmt/format.h>
65
66#include <opm/common/utility/gpuDecorators.hpp>
67#if HAVE_CUDA
68#if USE_HIP
69#include <opm/simulators/linalg/gpuistl_hip/GpuBuffer.hpp>
70#include <opm/simulators/linalg/gpuistl_hip/GpuView.hpp>
71#include <opm/simulators/linalg/gpuistl_hip/MiniMatrix.hpp>
72#include <opm/simulators/linalg/gpuistl_hip/MiniVector.hpp>
73#else
78#endif
79#endif
80
81namespace Opm::Parameters {
82
83struct SeparateSparseSourceTerms { static constexpr bool value = false; };
84
85} // namespace Opm::Parameters
86
87namespace Opm {
88
89// forward declarations
90template<class TypeTag>
92
93// Moved these structs out of the class to make them visible in the GPU code.
94template<class Storage = std::vector<int>>
96{
97 Storage cells;
98};
99
100#if HAVE_CUDA && OPM_IS_COMPILING_WITH_GPU_COMPILER
102 {
103 if (CPUDomain.cells.size() == 0) {
104 OPM_THROW(std::runtime_error, "Cannot copy empty full domain to GPU.");
105 }
106 return FullDomain<gpuistl::GpuBuffer<int>>{
107 gpuistl::GpuBuffer<int>(CPUDomain.cells)
108 };
109 };
110
111 FullDomain<gpuistl::GpuView<int>> make_view(FullDomain<gpuistl::GpuBuffer<int>>& buffer)
112 {
113 if (buffer.cells.size() == 0) {
114 OPM_THROW(std::runtime_error, "Cannot make view of empty full domain buffer.");
115 }
116 return FullDomain<gpuistl::GpuView<int>>{
117 gpuistl::make_view(buffer.cells)
118 };
119 };
120#endif
121
122template <class ResidualNBInfoType,class BlockType>
124{
125 unsigned int neighbor;
126 ResidualNBInfoType res_nbinfo;
127 BlockType* matBlockAddress;
128
129 template <class OtherBlockType>
131 : neighbor(other.neighbor)
132 , res_nbinfo(other.res_nbinfo)
133 , matBlockAddress(nullptr)
134 {
135 if (other.matBlockAddress) {
136 matBlockAddress = reinterpret_cast<BlockType*>(other.matBlockAddress);
137 }
138 }
139
140 template <class PtrType>
141 NeighborInfoStruct(unsigned int n, const ResidualNBInfoType& r, PtrType ptr)
142 : neighbor(n)
143 , res_nbinfo(r)
144 , matBlockAddress(static_cast<BlockType*>(ptr))
145 {
146 }
147
148 // Add a default constructor
150 : neighbor(0)
151 , res_nbinfo()
152 , matBlockAddress(nullptr)
153 {
154 }
155};
156
157#if HAVE_CUDA && OPM_IS_COMPILING_WITH_GPU_COMPILER
158namespace gpuistl {
159 template<class MiniMatrixType, class MatrixBlockType, class ResidualNBInfoType>
160 auto copy_to_gpu(const SparseTable<NeighborInfoStruct<ResidualNBInfoType, MatrixBlockType>>& cpuNeighborInfoTable)
161 {
162 // Convert the DUNE FieldMatrix/MatrixBlock to MiniMatrix types
163 using StructWithMinimatrix = NeighborInfoStruct<ResidualNBInfoType, MiniMatrixType>;
164 std::vector<StructWithMinimatrix> minimatrices(cpuNeighborInfoTable.dataSize());
165 size_t idx = 0;
166 for (auto e : cpuNeighborInfoTable.dataStorage()) {
167 minimatrices[idx++] = StructWithMinimatrix(e);
168 }
169
170 return SparseTable<StructWithMinimatrix, gpuistl::GpuBuffer>(
171 gpuistl::GpuBuffer<StructWithMinimatrix>(minimatrices),
172 gpuistl::GpuBuffer<int>(cpuNeighborInfoTable.rowStarts())
173 );
174 }
175}
176#endif
177
187template<class TypeTag>
189{
197
207
208 using Element = typename GridView::template Codim<0>::Entity;
209 using ElementIterator = typename GridView::template Codim<0>::Iterator;
210
211 using Vector = GlobalEqVector;
212
213 enum { numEq = getPropValue<TypeTag, Properties::NumEq>() };
214 enum { historySize = getPropValue<TypeTag, Properties::TimeDiscHistorySize>() };
215 enum { dimWorld = GridView::dimensionworld };
216
217 using MatrixBlock = typename SparseMatrixAdapter::MatrixBlock;
218 using VectorBlock = Dune::FieldVector<Scalar, numEq>;
220
221#if HAVE_CUDA && OPM_IS_COMPILING_WITH_GPU_COMPILER
222 using MatrixBlockGPU = gpuistl::MiniMatrix<Scalar, numEq * numEq>;
223 using VectorBlockGPU = gpuistl::MiniVector<Scalar, numEq>;
224#endif
225
226 static constexpr bool linearizeNonLocalElements =
227 getPropValue<TypeTag, Properties::LinearizeNonLocalElements>();
228 static constexpr bool enableFullyImplicitThermal = (getPropValue<TypeTag, Properties::EnergyModuleType>() == EnergyModules::FullyImplicitThermal);
229 static constexpr bool enableDiffusion = getPropValue<TypeTag, Properties::EnableDiffusion>();
230 static constexpr bool enableDispersion = getPropValue<TypeTag, Properties::EnableDispersion>();
231 static const bool enableBioeffects = getPropValue<TypeTag, Properties::EnableBioeffects>();
232
233 // copying the linearizer is not a good idea
234 TpfaLinearizer(const TpfaLinearizer&) = delete;
236
237public:
239 {
240 simulatorPtr_ = nullptr;
241 separateSparseSourceTerms_ = Parameters::Get<Parameters::SeparateSparseSourceTerms>();
242 exportIndex_=-1;
243 exportCount_=-1;
244 }
245
249 static void registerParameters()
250 {
251 Parameters::Register<Parameters::SeparateSparseSourceTerms>
252 ("Treat well source terms all in one go, instead of on a cell by cell basis.");
253 }
254
264 void init(Simulator& simulator)
265 {
266 simulatorPtr_ = &simulator;
267 eraseMatrix();
268 }
269
278 {
279 jacobian_.reset();
280 }
281
292 {
295 }
296
308 {
309 int succeeded;
310 try {
311 linearizeDomain(fullDomain_);
312 succeeded = 1;
313 }
314 catch (const std::exception& e) {
315 std::cout << "rank " << simulator_().gridView().comm().rank()
316 << " caught an exception while linearizing:" << e.what()
317 << "\n" << std::flush;
318 succeeded = 0;
319 }
320 catch (...) {
321 std::cout << "rank " << simulator_().gridView().comm().rank()
322 << " caught an exception while linearizing"
323 << "\n" << std::flush;
324 succeeded = 0;
325 }
326 succeeded = simulator_().gridView().comm().min(succeeded);
327
328 if (!succeeded) {
329 throw NumericalProblem("A process did not succeed in linearizing the system");
330 }
331 }
332
346 template <class SubDomainType>
347 void linearizeDomain(const SubDomainType& domain, const bool isNlddLocalSolve = false)
348 {
349 OPM_TIMEBLOCK(linearizeDomain);
350 // we defer the initialization of the Jacobian matrix until here because the
351 // auxiliary modules usually assume the problem, model and grid to be fully
352 // initialized...
353 if (!jacobian_) {
354 initFirstIteration_();
355 }
356
357 // Called here because it is no longer called from linearize_().
358 if (isNlddLocalSolve) {
359 resetSystem_(domain);
360 }
361 else {
362 resetSystem_();
363 }
364
365 linearize_(domain, isNlddLocalSolve);
366 }
367
368 void finalize()
369 { jacobian_->finalize(); }
370
376 {
377 OPM_TIMEBLOCK(linearizeAuxilaryEquations);
378 // flush possible local caches into matrix structure
379 jacobian_->commit();
380
381 auto& model = model_();
382 const auto& comm = simulator_().gridView().comm();
383 for (unsigned auxModIdx = 0; auxModIdx < model.numAuxiliaryModules(); ++auxModIdx) {
384 bool succeeded = true;
385 try {
386 model.auxiliaryModule(auxModIdx)->linearize(*jacobian_, residual_);
387 }
388 catch (const std::exception& e) {
389 succeeded = false;
390
391 std::cout << "rank " << simulator_().gridView().comm().rank()
392 << " caught an exception while linearizing:" << e.what()
393 << "\n" << std::flush;
394 }
395
396 succeeded = comm.min(succeeded);
397
398 if (!succeeded) {
399 throw NumericalProblem("linearization of an auxiliary equation failed");
400 }
401 }
402 }
403
407 const SparseMatrixAdapter& jacobian() const
408 { return *jacobian_; }
409
410 SparseMatrixAdapter& jacobian()
411 { return *jacobian_; }
412
416 const GlobalEqVector& residual() const
417 { return residual_; }
418
419 GlobalEqVector& residual()
420 { return residual_; }
421
425 void exportSystem(const int idx, std::string& tag, const char *path="export")
426 {
427 const bool export_sparsity = exportIndex_ == -1;
428
429 // increment indices and generate tag
430 exportCount_ = exportIndex_ == idx ? ++exportCount_ : 0;
431 exportIndex_ = idx;
432 tag = fmt::format(fmt::runtime("_{:03d}_{:02d}"), exportIndex_, exportCount_);
433
434 fmt::print(fmt::runtime("index = {:d}\n"), exportIndex_);
435 fmt::print(fmt::runtime("count = {:d}\n"), exportCount_);
436
437 Opm::exportSystem(jacobian_->istlMatrix(), residual_, export_sparsity, tag.c_str(), path);
438 }
439
441 { linearizationType_ = linearizationType; }
442
444 { return linearizationType_; }
445
451 const auto& getFlowsInfo() const
452 { return flowsInfo_; }
453
459 const auto& getFloresInfo() const
460 { return floresInfo_; }
461
468 const auto& getVelocityInfo() const
469 { return velocityInfo_; }
470
471 const auto& getNeighborInfo() const {
472 return neighborInfo_;
473 }
474
475
477 {
478 updateStoredTransmissibilities();
479 }
480
482 {
483 for (auto& bdyInfo : boundaryInfo_) {
484 const auto [type, massrateAD] = problem_().boundaryCondition(bdyInfo.cell, bdyInfo.dir);
485
486 // Strip the unnecessary (and zero anyway) derivatives off massrate.
487 VectorBlock massrate(0.0);
488 for (std::size_t ii = 0; ii < massrate.size(); ++ii) {
489 massrate[ii] = massrateAD[ii].value();
490 }
491 if (type != BCType::NONE) {
492 const auto& exFluidState = problem_().boundaryFluidState(bdyInfo.cell, bdyInfo.dir);
493 bdyInfo.bcdata.type = type;
494 bdyInfo.bcdata.massRate = massrate;
495 bdyInfo.bcdata.exFluidState = exFluidState;
496 }
497 }
498 }
499
505 std::map<unsigned, Constraints> constraintsMap() const
506 { return {}; }
507
508 template <class SubDomainType>
509 void resetSystem_(const SubDomainType& domain)
510 {
511 if (!jacobian_) {
512 initFirstIteration_();
513 }
514 for (int globI : domain.cells) {
515 residual_[globI] = 0.0;
516 jacobian_->clearRow(globI, 0.0);
517 }
518 }
519
520private:
521 Simulator& simulator_()
522 { return *simulatorPtr_; }
523
524 const Simulator& simulator_() const
525 { return *simulatorPtr_; }
526
527 Problem& problem_()
528 { return simulator_().problem(); }
529
530 const Problem& problem_() const
531 { return simulator_().problem(); }
532
533 Model& model_()
534 { return simulator_().model(); }
535
536 const Model& model_() const
537 { return simulator_().model(); }
538
539 const GridView& gridView_() const
540 { return problem_().gridView(); }
541
542 void initFirstIteration_()
543 {
544 // initialize the BCRS matrix for the Jacobian of the residual function
545 createMatrix_();
546
547 // initialize the Jacobian matrix and the vector for the residual function
548 residual_.resize(model_().numTotalDof());
549 resetSystem_();
550
551 // initialize the sparse tables for Flows and Flores
552 createFlows_();
553 }
554
555 // Construct the BCRS matrix for the Jacobian of the residual function
556 void createMatrix_()
557 {
558 OPM_TIMEBLOCK(createMatrix);
559 if (!neighborInfo_.empty()) {
560 // It is ok to call this function multiple times, but it
561 // should not do anything if already called.
562 return;
563 }
564 const auto& model = model_();
565 Stencil stencil(gridView_(), model_().dofMapper());
566
567 // for the main model, find out the global indices of the neighboring degrees of
568 // freedom of each primary degree of freedom
569 using NeighborSet = std::set<unsigned>;
570 std::vector<NeighborSet> sparsityPattern(model.numTotalDof());
571 const Scalar gravity = problem_().gravity()[dimWorld - 1];
572 unsigned numCells = model.numTotalDof();
573 neighborInfo_.reserve(numCells, 6 * numCells);
574 std::vector<NeighborInfoCPU> loc_nbinfo;
575 for (const auto& elem : elements(gridView_())) {
576 stencil.update(elem);
577
578 for (unsigned primaryDofIdx = 0; primaryDofIdx < stencil.numPrimaryDof(); ++primaryDofIdx) {
579 const unsigned myIdx = stencil.globalSpaceIndex(primaryDofIdx);
580 loc_nbinfo.resize(stencil.numDof() - 1); // Do not include the primary dof in neighborInfo_
581
582 for (unsigned dofIdx = 0; dofIdx < stencil.numDof(); ++dofIdx) {
583 const unsigned neighborIdx = stencil.globalSpaceIndex(dofIdx);
584 sparsityPattern[myIdx].insert(neighborIdx);
585 if (dofIdx > 0) {
586 const Scalar trans = problem_().transmissibility(myIdx, neighborIdx);
587 const auto scvfIdx = dofIdx - 1;
588 const auto& scvf = stencil.interiorFace(scvfIdx);
589 const Scalar area = scvf.area();
590 const Scalar Vin = problem_().model().dofTotalVolume(myIdx);
591 const Scalar Vex = problem_().model().dofTotalVolume(neighborIdx);
592 const Scalar zIn = problem_().dofCenterDepth(myIdx);
593 const Scalar zEx = problem_().dofCenterDepth(neighborIdx);
594 const Scalar dZg = (zIn - zEx)*gravity;
595 const Scalar thpres = problem_().thresholdPressure(myIdx, neighborIdx);
596 const auto dirId = scvf.dirId();
597 auto faceDir = dirId < 0 ? FaceDir::DirEnum::Unknown
598 : FaceDir::FromIntersectionIndex(dirId);
599 ResidualNBInfo nbinfo{trans, area, thpres, dZg, faceDir, Vin, Vex, {}, {}, {}, {}};
600 if constexpr (enableFullyImplicitThermal) {
601 nbinfo.inAlpha = problem_().thermalHalfTransmissibility(myIdx, neighborIdx);
602 nbinfo.outAlpha = problem_().thermalHalfTransmissibility(neighborIdx, myIdx);
603 }
604 if constexpr (enableDiffusion) {
605 nbinfo.diffusivity = problem_().diffusivity(myIdx, neighborIdx);
606 }
607 if constexpr (enableDispersion) {
608 nbinfo.dispersivity = problem_().dispersivity(myIdx, neighborIdx);
609 }
610 loc_nbinfo[dofIdx - 1] = NeighborInfoCPU{neighborIdx, nbinfo, nullptr};
611 }
612 }
613 neighborInfo_.appendRow(loc_nbinfo.begin(), loc_nbinfo.end());
614 if (problem_().nonTrivialBoundaryConditions()) {
615 for (unsigned bfIndex = 0; bfIndex < stencil.numBoundaryFaces(); ++bfIndex) {
616 const auto& bf = stencil.boundaryFace(bfIndex);
617 const int dir_id = bf.dirId();
618 // not for NNCs
619 if (dir_id < 0) {
620 continue;
621 }
622 const auto [type, massrateAD] = problem_().boundaryCondition(myIdx, dir_id);
623 // Strip the unnecessary (and zero anyway) derivatives off massrate.
624 VectorBlock massrate(0.0);
625 for (std::size_t ii = 0; ii < massrate.size(); ++ii) {
626 massrate[ii] = massrateAD[ii].value();
627 }
628 const auto& exFluidState = problem_().boundaryFluidState(myIdx, dir_id);
629 BoundaryConditionData bcdata{type,
630 massrate,
631 exFluidState.pvtRegionIndex(),
632 bfIndex,
633 bf.area(),
634 bf.integrationPos()[dimWorld - 1],
635 exFluidState};
636 boundaryInfo_.push_back({myIdx, dir_id, bfIndex, bcdata});
637 }
638 }
639 }
640 }
641
642 // add the additional neighbors and degrees of freedom caused by the auxiliary
643 // equations
644 const std::size_t numAuxMod = model.numAuxiliaryModules();
645 for (unsigned auxModIdx = 0; auxModIdx < numAuxMod; ++auxModIdx) {
646 model.auxiliaryModule(auxModIdx)->addNeighbors(sparsityPattern);
647 }
648
649 // allocate raw matrix
650 jacobian_ = std::make_unique<SparseMatrixAdapter>(simulator_());
651 diagMatAddress_.resize(numCells);
652 // create matrix structure based on sparsity pattern
653 jacobian_->reserve(sparsityPattern);
654 for (unsigned globI = 0; globI < numCells; globI++) {
655 const auto& nbInfos = neighborInfo_[globI];
656 diagMatAddress_[globI] = jacobian_->blockAddress(globI, globI);
657 for (auto& nbInfo : nbInfos) {
658 nbInfo.matBlockAddress = jacobian_->blockAddress(nbInfo.neighbor, globI);
659 }
660 }
661
662 // Create dummy full domain.
663 fullDomain_.cells.resize(numCells);
664 std::iota(fullDomain_.cells.begin(), fullDomain_.cells.end(), 0);
665 }
666
667 // reset the global linear system of equations.
668 void resetSystem_()
669 {
670 residual_ = 0.0;
671 // zero all matrix entries
672 jacobian_->clear();
673 }
674
675 // Initialize the flows, flores, and velocity sparse tables
676 void createFlows_()
677 {
678 OPM_TIMEBLOCK(createFlows);
679 // If FLOWS/FLORES is set in any RPTRST in the schedule, then we initializate the sparse tables.
680 // If DISPERC is in the deck, we initialize the sparse table here as well.
681 const bool anyFlows = simulator_().problem().eclWriter().outputModule().getFlows().anyFlows();
682 const auto& blockFlows = simulator_().problem().eclWriter().outputModule().getFlows().blockFlows();
683 const auto& blockVelocity = simulator_().problem().eclWriter().outputModule().getFlows().blockVelocity();
684 const bool isTemp = simulator_().vanguard().eclState().getSimulationConfig().isTemp();
685 const bool anyFlores = simulator_().problem().eclWriter().outputModule().getFlows().anyFlores() || isTemp;
686 const bool dispersionActive = simulator_().vanguard().eclState().getSimulationConfig().rock_config().dispersion();
687 if (!dispersionActive && !enableBioeffects && blockVelocity.empty()
688 && !((anyFlows || !blockFlows.empty()) && flowsInfo_.empty())
689 && !(anyFlores && floresInfo_.empty())) {
690 return;
691 }
692 const auto& model = model_();
693 const auto& nncOutput = simulator_().problem().eclWriter().getOutputNnc();
694 Stencil stencil(gridView_(), model_().dofMapper());
695 const unsigned numCells = model.numTotalDof();
696 std::unordered_multimap<int, std::pair<int, int>> nncIndices;
697 std::vector<FlowInfo> loc_flinfo;
698 std::vector<VelocityInfo> loc_vlinfo;
699 unsigned int nncId = 0;
700 VectorBlock flow(0.0);
701
702 // Create a nnc structure to use fast lookup
703 for (unsigned nncIdx = 0; nncIdx < nncOutput.size(); ++nncIdx) {
704 const int ci1 = nncOutput[nncIdx].cell1;
705 const int ci2 = nncOutput[nncIdx].cell2;
706 nncIndices.emplace(ci1, std::make_pair(ci2, nncIdx));
707 }
708
709 if (anyFlows) {
710 flowsInfo_.reserve(numCells, 6 * numCells);
711 }
712 else if (!blockFlows.empty()) {
713 flowsInfo_.reserve(numCells, 6 * blockFlows.size());
714 }
715 if (anyFlores) {
716 floresInfo_.reserve(numCells, 6 * numCells);
717 }
718 if (dispersionActive || enableBioeffects) {
719 velocityInfo_.reserve(numCells, 6 * numCells);
720 }
721 else if (!blockVelocity.empty()) {
722 velocityInfo_.reserve(numCells, 6 * blockVelocity.size());
723 }
724
725 for (const auto& elem : elements(gridView_())) {
726 stencil.update(elem);
727 for (unsigned primaryDofIdx = 0; primaryDofIdx < stencil.numPrimaryDof(); ++primaryDofIdx) {
728 const unsigned myIdx = stencil.globalSpaceIndex(primaryDofIdx);
729 bool blockFlowFound = false;
730 bool blockVelocityFound = false;
731 if (!blockFlows.empty()) {
732 if (std::ranges::binary_search(blockFlows,
733 simulator_().vanguard().cartesianIndex(myIdx))) {
734 blockFlowFound = true;
735 }
736 else {
737 flowsInfo_.appendRow(loc_flinfo.begin(), loc_flinfo.begin());
738 if (!dispersionActive && !enableBioeffects && !anyFlores && blockVelocity.empty()) {
739 continue;
740 }
741 }
742 }
743 if (!blockVelocity.empty() && !(dispersionActive || enableBioeffects)) {
744 if (std::ranges::binary_search(blockVelocity,
745 simulator_().vanguard().cartesianIndex(myIdx))) {
746 blockVelocityFound = true;
747 }
748 else {
749 velocityInfo_.appendRow(loc_vlinfo.begin(), loc_vlinfo.begin());
750 if (!anyFlows && blockFlows.empty() && !anyFlores) {
751 continue;
752 }
753 }
754 }
755 const int numFaces = stencil.numBoundaryFaces() + stencil.numInteriorFaces();
756 loc_flinfo.resize(numFaces);
757 loc_vlinfo.resize(stencil.numDof() - 1);
758
759 for (unsigned dofIdx = 0; dofIdx < stencil.numDof(); ++dofIdx) {
760 const unsigned neighborIdx = stencil.globalSpaceIndex(dofIdx);
761 if (dofIdx > 0) {
762 const auto scvfIdx = dofIdx - 1;
763 const auto& scvf = stencil.interiorFace(scvfIdx);
764 int faceId = scvf.dirId();
765 const int cartMyIdx = simulator_().vanguard().cartesianIndex(myIdx);
766 const int cartNeighborIdx = simulator_().vanguard().cartesianIndex(neighborIdx);
767 const auto& range = nncIndices.equal_range(cartMyIdx);
768 for (auto it = range.first; it != range.second; ++it) {
769 if (it->second.first == cartNeighborIdx){
770 // -1 gives problem since is used for the nncInput from the deck
771 faceId = -2;
772 // the index is stored to be used for writting the outputs
773 nncId = it->second.second;
774 }
775 }
776 loc_flinfo[dofIdx - 1] = FlowInfo{faceId, flow, nncId};
777 loc_vlinfo[dofIdx - 1] = VelocityInfo{faceId, flow};
778 }
779 }
780
781 for (unsigned bdfIdx = 0; bdfIdx < stencil.numBoundaryFaces(); ++bdfIdx) {
782 const auto& scvf = stencil.boundaryFace(bdfIdx);
783 const int faceId = scvf.dirId();
784 loc_flinfo[stencil.numInteriorFaces() + bdfIdx] = FlowInfo{faceId, flow, nncId};
785 }
786
787 if (anyFlows || blockFlowFound) {
788 flowsInfo_.appendRow(loc_flinfo.begin(), loc_flinfo.end());
789 }
790 if (anyFlores) {
791 floresInfo_.appendRow(loc_flinfo.begin(), loc_flinfo.end());
792 }
793 if (dispersionActive || enableBioeffects || blockVelocityFound) {
794 velocityInfo_.appendRow(loc_vlinfo.begin(), loc_vlinfo.end());
795 }
796 }
797 }
798 }
799
800public:
801 void setResAndJacobi(VectorBlock& res, MatrixBlock& bMat, const ADVectorBlock& resid) const
802 {
803 for (unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
804 res[eqIdx] = resid[eqIdx].value();
805 }
806
807 for (unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
808 for (unsigned pvIdx = 0; pvIdx < numEq; ++pvIdx) {
809 // A[dofIdx][focusDofIdx][eqIdx][pvIdx] is the partial derivative of
810 // the residual function 'eqIdx' for the degree of freedom 'dofIdx'
811 // with regard to the focus variable 'pvIdx' of the degree of freedom
812 // 'focusDofIdx'
813 bMat[eqIdx][pvIdx] = resid[eqIdx].derivative(pvIdx);
814 }
815 }
816 }
817
819 {
820 OPM_TIMEBLOCK(updateFlows);
821 const bool enableFlows = simulator_().problem().eclWriter().outputModule().getFlows().hasFlows();
822 const auto& blockFlows = simulator_().problem().eclWriter().outputModule().getFlows().blockFlows();
823 // We reuse the fluxes in the TEMP option
824 const bool isTemp = simulator_().vanguard().eclState().getSimulationConfig().isTemp();
825 const bool enableFlores = simulator_().problem().eclWriter().outputModule().getFlows().hasFlores() || isTemp;
826 if (!enableFlows && !enableFlores && blockFlows.empty()) {
827 return;
828 }
829 const unsigned int numCells = model_().numTotalDof();
830#ifdef _OPENMP
831#pragma omp parallel for
832#endif
833 for (unsigned globI = 0; globI < numCells; ++globI) {
834 OPM_TIMEBLOCK_LOCAL(linearizationForEachCell, Subsystem::Assembly);
835 const auto& nbInfos = neighborInfo_[globI];
836 ADVectorBlock adres(0.0);
837 ADVectorBlock darcyFlux(0.0);
838 const IntensiveQuantities& intQuantsIn = model_().intensiveQuantities(globI, /*timeIdx*/ 0);
839 // Flux term.
840 {
841 OPM_TIMEBLOCK_LOCAL(fluxCalculationForEachCell, Subsystem::Assembly);
842 short loc = 0;
843 for (const auto& nbInfo : nbInfos) {
844 OPM_TIMEBLOCK_LOCAL(fluxCalculationForEachFace, Subsystem::Assembly);
845 const unsigned globJ = nbInfo.neighbor;
846 assert(globJ != globI);
847 adres = 0.0;
848 darcyFlux = 0.0;
849 const IntensiveQuantities& intQuantsEx = model_().intensiveQuantities(globJ, /*timeIdx*/ 0);
850 LocalResidual::computeFlux(adres, darcyFlux, globI, globJ, intQuantsIn,
851 intQuantsEx, nbInfo.res_nbinfo, problem_().moduleParams());
852 adres *= nbInfo.res_nbinfo.faceArea;
853 if (!blockFlows.empty()) {
854 if (std::ranges::binary_search(blockFlows,
855 simulator_().vanguard().cartesianIndex(globI))) {
856 for (unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
857 flowsInfo_[globI][loc].flow[eqIdx] = adres[eqIdx].value();
858 }
859 }
860 }
861 else if (enableFlows) {
862 for (unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
863 flowsInfo_[globI][loc].flow[eqIdx] = adres[eqIdx].value();
864 }
865 }
866 if (enableFlores) {
867 for (unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
868 floresInfo_[globI][loc].flow[eqIdx] = darcyFlux[eqIdx].value();
869 }
870 }
871 ++loc;
872 }
873 }
874 }
875
876 // Boundary terms. Only looping over cells with nontrivial bcs.
877 for (const auto& bdyInfo : boundaryInfo_) {
878 if (bdyInfo.bcdata.type == BCType::NONE) {
879 continue;
880 }
881
882 ADVectorBlock adres(0.0);
883 const unsigned globI = bdyInfo.cell;
884 const auto& nbInfos = neighborInfo_[globI];
885 const IntensiveQuantities& insideIntQuants = model_().intensiveQuantities(globI, /*timeIdx*/ 0);
886 LocalResidual::computeBoundaryFlux(adres, problem_(), bdyInfo.bcdata, insideIntQuants, globI);
887 adres *= bdyInfo.bcdata.faceArea;
888 const unsigned bfIndex = bdyInfo.bfIndex;
889 if (enableFlows) {
890 for (unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
891 flowsInfo_[globI][nbInfos.size() + bfIndex].flow[eqIdx] = adres[eqIdx].value();
892 }
893 }
894 // TODO also store Flores?
895 }
896 }
897
898private:
899 template <class SubDomainType>
900 void linearize_(const SubDomainType& domain, bool isNlddLocalSolve)
901 {
902 // This check should be removed once this is addressed by
903 // for example storing the previous timesteps' values for
904 // rsmax (for DRSDT) and similar.
905 if (!problem_().recycleFirstIterationStorage()) {
906 if (!model_().storeIntensiveQuantities() && !model_().enableStorageCache()) {
907 OPM_THROW(std::runtime_error, "Must have cached either IQs or storage when we cannot recycle.");
908 }
909 }
910
911 OPM_TIMEBLOCK(linearize);
912
913 // We do not call resetSystem_() here, since that will set
914 // the full system to zero, not just our part.
915 // Instead, that must be called before starting the linearization.
916 const bool dispersionActive = simulator_().vanguard().eclState().getSimulationConfig().rock_config().dispersion();
917 const auto& blockVelocity = simulator_().problem().eclWriter().outputModule().getFlows().blockVelocity();
918 const unsigned int numCells = domain.cells.size();
919
920 // Fetch timestepsize used later in accumulation term.
921 const double dt = simulator_().timeStepSize();
922
923#ifdef _OPENMP
924#pragma omp parallel for
925#endif
926 for (unsigned ii = 0; ii < numCells; ++ii) {
927 OPM_TIMEBLOCK_LOCAL(linearizationForEachCell, Subsystem::Assembly);
928 const unsigned globI = domain.cells[ii];
929 const auto& nbInfos = neighborInfo_[globI];
930 VectorBlock res(0.0);
931 MatrixBlock bMat(0.0);
932 ADVectorBlock adres(0.0);
933 ADVectorBlock darcyFlux(0.0);
934 const IntensiveQuantities& intQuantsIn = model_().intensiveQuantities(globI, /*timeIdx*/ 0);
935
936 // Flux term.
937 {
938 OPM_TIMEBLOCK_LOCAL(fluxCalculationForEachCell, Subsystem::Assembly);
939 short loc = 0;
940 for (const auto& nbInfo : nbInfos) {
941 OPM_TIMEBLOCK_LOCAL(fluxCalculationForEachFace, Subsystem::Assembly);
942 const unsigned globJ = nbInfo.neighbor;
943 assert(globJ != globI);
944 res = 0.0;
945 bMat = 0.0;
946 adres = 0.0;
947 darcyFlux = 0.0;
948 const IntensiveQuantities& intQuantsEx = model_().intensiveQuantities(globJ, /*timeIdx*/ 0);
949 LocalResidual::computeFlux(adres, darcyFlux, globI, globJ, intQuantsIn, intQuantsEx,
950 nbInfo.res_nbinfo, problem_().moduleParams());
951 adres *= nbInfo.res_nbinfo.faceArea;
952 if (dispersionActive || enableBioeffects) {
953 for (unsigned phaseIdx = 0; phaseIdx < numEq; ++phaseIdx) {
954 velocityInfo_[globI][loc].velocity[phaseIdx] =
955 darcyFlux[phaseIdx].value() / nbInfo.res_nbinfo.faceArea;
956 }
957 }
958 else if (!blockVelocity.empty()) {
959 if (std::ranges::binary_search(blockVelocity,
960 simulator_().vanguard().cartesianIndex(globI))) {
961 for (unsigned phaseIdx = 0; phaseIdx < numEq; ++phaseIdx) {
962 velocityInfo_[globI][loc].velocity[phaseIdx] =
963 darcyFlux[phaseIdx].value() / nbInfo.res_nbinfo.faceArea;
964 }
965 }
966 }
967 setResAndJacobi(res, bMat, adres);
968 residual_[globI] += res;
969 //SparseAdapter syntax: jacobian_->addToBlock(globI, globI, bMat);
970 *diagMatAddress_[globI] += bMat;
971 bMat *= -1.0;
972 //SparseAdapter syntax: jacobian_->addToBlock(globJ, globI, bMat);
973 *nbInfo.matBlockAddress += bMat;
974 ++loc;
975 }
976 }
977
978 // Accumulation term.
979 const double volume = model_().dofTotalVolume(globI);
980 const Scalar storefac = volume / dt;
981 adres = 0.0;
982 {
983 OPM_TIMEBLOCK_LOCAL(computeStorage, Subsystem::Assembly);
984 LocalResidual::template computeStorage<Evaluation>(adres, intQuantsIn);
985 }
986 setResAndJacobi(res, bMat, adres);
987 // Either use cached storage term, or compute it on the fly.
988 if (model_().enableStorageCache()) {
989 // The cached storage for timeIdx 0 (current time) is not
990 // used, but after storage cache is shifted at the end of the
991 // timestep, it will become cached storage for timeIdx 1.
992 model_().updateCachedStorage(globI, /*timeIdx=*/0, res);
993 // We should not update the storage cache here for NLDD local solves.
994 // This will reset the start-of-step storage to incorrect numbers when
995 // we do local solves, where the iteration number will start from 0,
996 // but the starting state may not be identical to the start-of-step state.
997 // Note that a full assembly must be done before local solves
998 // otherwise this will be left un-updated.
999 if (model_().newtonMethod().numIterations() == 0 && !isNlddLocalSolve) {
1000 // Need to update the storage cache.
1001 if (problem_().recycleFirstIterationStorage()) {
1002 // Assumes nothing have changed in the system which
1003 // affects masses calculated from primary variables.
1004 model_().updateCachedStorage(globI, /*timeIdx=*/1, res);
1005 }
1006 else {
1007 Dune::FieldVector<Scalar, numEq> tmp;
1008 const IntensiveQuantities intQuantOld = model_().intensiveQuantities(globI, 1);
1009 LocalResidual::template computeStorage<Scalar>(tmp, intQuantOld);
1010 model_().updateCachedStorage(globI, /*timeIdx=*/1, tmp);
1011 }
1012 }
1013 res -= model_().cachedStorage(globI, 1);
1014 }
1015 else {
1016 OPM_TIMEBLOCK_LOCAL(computeStorage0, Subsystem::Assembly);
1017 Dune::FieldVector<Scalar, numEq> tmp;
1018 const IntensiveQuantities intQuantOld = model_().intensiveQuantities(globI, 1);
1019 LocalResidual::template computeStorage<Scalar>(tmp, intQuantOld);
1020 // assume volume do not change
1021 res -= tmp;
1022 }
1023 res *= storefac;
1024 bMat *= storefac;
1025 residual_[globI] += res;
1026 //SparseAdapter syntax: jacobian_->addToBlock(globI, globI, bMat);
1027 *diagMatAddress_[globI] += bMat;
1028
1029 // Cell-wise source terms.
1030 // This will include well sources if SeparateSparseSourceTerms is false.
1031 res = 0.0;
1032 bMat = 0.0;
1033 adres = 0.0;
1034 if (separateSparseSourceTerms_) {
1035 LocalResidual::computeSourceDense(adres, problem_(), intQuantsIn, globI, 0);
1036 }
1037 else {
1038 LocalResidual::computeSource(adres, problem_(), intQuantsIn, globI, 0);
1039 }
1040 adres *= -volume;
1041 setResAndJacobi(res, bMat, adres);
1042 residual_[globI] += res;
1043 //SparseAdapter syntax: jacobian_->addToBlock(globI, globI, bMat);
1044 *diagMatAddress_[globI] += bMat;
1045 } // end of loop for cell globI.
1046
1047 // Add sparse source terms. For now only wells.
1048 if (separateSparseSourceTerms_) {
1049 problem_().wellModel().addReservoirSourceTerms(residual_, diagMatAddress_);
1050 }
1051
1052 // Boundary terms. Only looping over cells with nontrivial bcs.
1053 for (const auto& bdyInfo : boundaryInfo_) {
1054 if (bdyInfo.bcdata.type == BCType::NONE) {
1055 continue;
1056 }
1057
1058 VectorBlock res(0.0);
1059 MatrixBlock bMat(0.0);
1060 ADVectorBlock adres(0.0);
1061 const unsigned globI = bdyInfo.cell;
1062 const IntensiveQuantities& insideIntQuants = model_().intensiveQuantities(globI, /*timeIdx*/ 0);
1063 LocalResidual::computeBoundaryFlux(adres, problem_(), bdyInfo.bcdata, insideIntQuants, globI);
1064 adres *= bdyInfo.bcdata.faceArea;
1065 setResAndJacobi(res, bMat, adres);
1066 residual_[globI] += res;
1068 *diagMatAddress_[globI] += bMat;
1069 }
1070 }
1071
1072 void updateStoredTransmissibilities()
1073 {
1074 if (neighborInfo_.empty()) {
1075 // This function was called before createMatrix_() was called.
1076 // We call initFirstIteration_(), not createMatrix_(), because
1077 // that will also initialize the residual consistently.
1078 initFirstIteration_();
1079 }
1080
1081 const unsigned numCells = model_().numTotalDof();
1082#ifdef _OPENMP
1083#pragma omp parallel for
1084#endif
1085 for (unsigned globI = 0; globI < numCells; globI++) {
1086 auto nbInfos = neighborInfo_[globI]; // nbInfos will be a SparseTable<...>::mutable_iterator_range.
1087 for (auto& nbInfo : nbInfos) {
1088 const unsigned globJ = nbInfo.neighbor;
1089 nbInfo.res_nbinfo.trans = problem_().transmissibility(globI, globJ);
1090 }
1091 }
1092 }
1093
1094 Simulator* simulatorPtr_{};
1095
1096 // the jacobian matrix
1097 std::unique_ptr<SparseMatrixAdapter> jacobian_{};
1098
1099 // the right-hand side
1100 GlobalEqVector residual_;
1101
1102 LinearizationType linearizationType_{};
1103
1104 using ResidualNBInfo = typename LocalResidual::ResidualNBInfo;
1105 using NeighborInfoCPU = NeighborInfoStruct<ResidualNBInfo, MatrixBlock>;
1106
1107 SparseTable<NeighborInfoCPU> neighborInfo_{};
1108 std::vector<MatrixBlock*> diagMatAddress_{};
1109
1110 struct FlowInfo
1111 {
1112 int faceId;
1113 VectorBlock flow;
1114 unsigned int nncId;
1115 };
1116 SparseTable<FlowInfo> flowsInfo_;
1117 SparseTable<FlowInfo> floresInfo_;
1118
1119 struct VelocityInfo
1120 {
1121 int faceId;
1122 VectorBlock velocity;
1123 };
1124 SparseTable<VelocityInfo> velocityInfo_;
1125
1126 using ScalarFluidState = typename IntensiveQuantities::ScalarFluidState;
1127 struct BoundaryConditionData
1128 {
1129 BCType type;
1130 VectorBlock massRate;
1131 unsigned pvtRegionIdx;
1132 unsigned boundaryFaceIndex;
1133 double faceArea;
1134 double faceZCoord;
1135 ScalarFluidState exFluidState;
1136 };
1137
1138 struct BoundaryInfo
1139 {
1140 unsigned int cell;
1141 int dir;
1142 unsigned int bfIndex;
1143 BoundaryConditionData bcdata;
1144 };
1145 std::vector<BoundaryInfo> boundaryInfo_;
1146
1147 bool separateSparseSourceTerms_ = false;
1148
1149 FullDomain<> fullDomain_;
1150
1151 int exportIndex_;
1152 int exportCount_;
1153};
1154} // namespace Opm
1155
1156#endif // TPFA_LINEARIZER
Declares the properties required by the black oil model.
The base class for the element-centered finite-volume discretization scheme.
Definition: ecfvdiscretization.hh:160
Definition: matrixblock.hh:229
Manages the initializing and running of time dependent problems.
Definition: simulator.hh:84
Scalar timeStepSize() const
Returns the time step length so that we don't miss the beginning of the next episode or cross the en...
Definition: simulator.hh:413
Vanguard & vanguard()
Return a reference to the grid manager of simulation.
Definition: simulator.hh:234
Problem & problem()
Return the object which specifies the pysical setup of the simulation.
Definition: simulator.hh:265
const GridView & gridView() const
Return the grid view for which the simulation is done.
Definition: simulator.hh:246
Model & model()
Return the physical model used in the simulation.
Definition: simulator.hh:252
The common code for the linearizers of non-linear systems of equations.
Definition: tpfalinearizer.hh:189
const auto & getFloresInfo() const
Return constant reference to the floresInfo.
Definition: tpfalinearizer.hh:459
const LinearizationType & getLinearizationType() const
Definition: tpfalinearizer.hh:443
const auto & getNeighborInfo() const
Definition: tpfalinearizer.hh:471
void updateBoundaryConditionData()
Definition: tpfalinearizer.hh:481
void linearize()
Linearize the full system of non-linear equations.
Definition: tpfalinearizer.hh:291
SparseMatrixAdapter & jacobian()
Definition: tpfalinearizer.hh:410
const auto & getFlowsInfo() const
Return constant reference to the flowsInfo.
Definition: tpfalinearizer.hh:451
std::map< unsigned, Constraints > constraintsMap() const
Returns the map of constraint degrees of freedom.
Definition: tpfalinearizer.hh:505
TpfaLinearizer()
Definition: tpfalinearizer.hh:238
void finalize()
Definition: tpfalinearizer.hh:368
void init(Simulator &simulator)
Initialize the linearizer.
Definition: tpfalinearizer.hh:264
void updateFlowsInfo()
Definition: tpfalinearizer.hh:818
void setResAndJacobi(VectorBlock &res, MatrixBlock &bMat, const ADVectorBlock &resid) const
Definition: tpfalinearizer.hh:801
static void registerParameters()
Register all run-time parameters for the Jacobian linearizer.
Definition: tpfalinearizer.hh:249
void exportSystem(const int idx, std::string &tag, const char *path="export")
Export block sparse linear system.
Definition: tpfalinearizer.hh:425
void linearizeDomain()
Linearize the part of the non-linear system of equations that is associated with the spatial domain.
Definition: tpfalinearizer.hh:307
const SparseMatrixAdapter & jacobian() const
Return constant reference to global Jacobian matrix backend.
Definition: tpfalinearizer.hh:407
void setLinearizationType(LinearizationType linearizationType)
Definition: tpfalinearizer.hh:440
void linearizeDomain(const SubDomainType &domain, const bool isNlddLocalSolve=false)
Linearize the part of the non-linear system of equations that is associated with a part of the spatia...
Definition: tpfalinearizer.hh:347
GlobalEqVector & residual()
Definition: tpfalinearizer.hh:419
void eraseMatrix()
Causes the Jacobian matrix to be recreated from scratch before the next iteration.
Definition: tpfalinearizer.hh:277
const auto & getVelocityInfo() const
Return constant reference to the velocityInfo.
Definition: tpfalinearizer.hh:468
void updateDiscretizationParameters()
Definition: tpfalinearizer.hh:476
void resetSystem_(const SubDomainType &domain)
Definition: tpfalinearizer.hh:509
void linearizeAuxiliaryEquations()
Linearize the part of the non-linear system of equations that is associated with the spatial domain.
Definition: tpfalinearizer.hh:375
const GlobalEqVector & residual() const
Return constant reference to global residual vector.
Definition: tpfalinearizer.hh:416
A small fixed-size square matrix class for use in CUDA kernels.
Definition: MiniMatrix.hpp:37
Definition: MiniVector.hpp:52
Declare the properties used by the infrastructure code of the finite volume discretizations.
Defines the common properties required by the porous medium multi-phase models.
@ NONE
Definition: DeferredLogger.hpp:46
Definition: blackoilnewtonmethodparams.hpp:31
HYPRE_IJMatrix createMatrix(HYPRE_Int N, HYPRE_Int dof_offset, const CommType &comm)
Create Hypre matrix.
Definition: HypreSetup.hpp:170
PointerView< T > make_view(const std::shared_ptr< T > &ptr)
Definition: gpu_smart_pointer.hpp:318
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
void exportSystem(const IstlMatrix &jacobian, const GlobalEqVector &residual, const bool export_sparsity, const char *tag, const char *path="export")
Export blocks-sparse linear system.
Definition: exportSystem.hpp:42
Definition: tpfalinearizer.hh:96
Storage cells
Definition: tpfalinearizer.hh:97
Definition: linearizationtype.hh:34
Definition: tpfalinearizer.hh:124
NeighborInfoStruct(unsigned int n, const ResidualNBInfoType &r, PtrType ptr)
Definition: tpfalinearizer.hh:141
NeighborInfoStruct()
Definition: tpfalinearizer.hh:149
NeighborInfoStruct(const NeighborInfoStruct< ResidualNBInfoType, OtherBlockType > &other)
Definition: tpfalinearizer.hh:130
BlockType * matBlockAddress
Definition: tpfalinearizer.hh:127
ResidualNBInfoType res_nbinfo
Definition: tpfalinearizer.hh:126
unsigned int neighbor
Definition: tpfalinearizer.hh:125
Definition: tpfalinearizer.hh:83
static constexpr bool value
Definition: tpfalinearizer.hh:83