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
345 template <class SubDomainType>
346 void linearizeDomain(const SubDomainType& domain)
347 {
348 OPM_TIMEBLOCK(linearizeDomain);
349 // we defer the initialization of the Jacobian matrix until here because the
350 // auxiliary modules usually assume the problem, model and grid to be fully
351 // initialized...
352 if (!jacobian_) {
353 initFirstIteration_();
354 }
355
356 // Called here because it is no longer called from linearize_().
357 if (problem_().iterationContext().inLocalSolve()) {
358 resetSystem_(domain);
359 }
360 else {
361 resetSystem_();
362 }
363
364 linearize_(domain);
365 }
366
367 void finalize()
368 { jacobian_->finalize(); }
369
375 {
376 OPM_TIMEBLOCK(linearizeAuxilaryEquations);
377 // flush possible local caches into matrix structure
378 jacobian_->commit();
379
380 auto& model = model_();
381 const auto& comm = simulator_().gridView().comm();
382 for (unsigned auxModIdx = 0; auxModIdx < model.numAuxiliaryModules(); ++auxModIdx) {
383 bool succeeded = true;
384 try {
385 model.auxiliaryModule(auxModIdx)->linearize(*jacobian_, residual_);
386 }
387 catch (const std::exception& e) {
388 succeeded = false;
389
390 std::cout << "rank " << simulator_().gridView().comm().rank()
391 << " caught an exception while linearizing:" << e.what()
392 << "\n" << std::flush;
393 }
394
395 succeeded = comm.min(succeeded);
396
397 if (!succeeded) {
398 throw NumericalProblem("linearization of an auxiliary equation failed");
399 }
400 }
401 }
402
406 const SparseMatrixAdapter& jacobian() const
407 { return *jacobian_; }
408
409 SparseMatrixAdapter& jacobian()
410 { return *jacobian_; }
411
415 const GlobalEqVector& residual() const
416 { return residual_; }
417
418 GlobalEqVector& residual()
419 { return residual_; }
420
424 void exportSystem(const int idx, std::string& tag, const char *path="export")
425 {
426 const bool export_sparsity = exportIndex_ == -1;
427
428 // increment indices and generate tag
429 exportCount_ = exportIndex_ == idx ? ++exportCount_ : 0;
430 exportIndex_ = idx;
431 tag = fmt::format(fmt::runtime("_{:03d}_{:02d}"), exportIndex_, exportCount_);
432
433 fmt::print(fmt::runtime("index = {:d}\n"), exportIndex_);
434 fmt::print(fmt::runtime("count = {:d}\n"), exportCount_);
435
436 Opm::exportSystem(jacobian_->istlMatrix(), residual_, export_sparsity, tag.c_str(), path);
437 }
438
440 { linearizationType_ = linearizationType; }
441
443 { return linearizationType_; }
444
450 const auto& getFlowsInfo() const
451 { return flowsInfo_; }
452
458 const auto& getFloresInfo() const
459 { return floresInfo_; }
460
467 const auto& getVelocityInfo() const
468 { return velocityInfo_; }
469
470 const auto& getNeighborInfo() const {
471 return neighborInfo_;
472 }
473
474
476 {
477 updateStoredTransmissibilities();
478 }
479
481 {
482 for (auto& bdyInfo : boundaryInfo_) {
483 const auto [type, massrateAD] = problem_().boundaryCondition(bdyInfo.cell, bdyInfo.dir);
484
485 // Strip the unnecessary (and zero anyway) derivatives off massrate.
486 VectorBlock massrate(0.0);
487 for (std::size_t ii = 0; ii < massrate.size(); ++ii) {
488 massrate[ii] = massrateAD[ii].value();
489 }
490 if (type != BCType::NONE) {
491 const auto& exFluidState = problem_().boundaryFluidState(bdyInfo.cell, bdyInfo.dir);
492 bdyInfo.bcdata.type = type;
493 bdyInfo.bcdata.massRate = massrate;
494 bdyInfo.bcdata.exFluidState = exFluidState;
495 }
496 }
497 }
498
504 std::map<unsigned, Constraints> constraintsMap() const
505 { return {}; }
506
507 template <class SubDomainType>
508 void resetSystem_(const SubDomainType& domain)
509 {
510 if (!jacobian_) {
511 initFirstIteration_();
512 }
513 for (int globI : domain.cells) {
514 residual_[globI] = 0.0;
515 jacobian_->clearRow(globI, 0.0);
516 }
517 }
518
519private:
520 Simulator& simulator_()
521 { return *simulatorPtr_; }
522
523 const Simulator& simulator_() const
524 { return *simulatorPtr_; }
525
526 Problem& problem_()
527 { return simulator_().problem(); }
528
529 const Problem& problem_() const
530 { return simulator_().problem(); }
531
532 Model& model_()
533 { return simulator_().model(); }
534
535 const Model& model_() const
536 { return simulator_().model(); }
537
538 const GridView& gridView_() const
539 { return problem_().gridView(); }
540
541 void initFirstIteration_()
542 {
543 // initialize the BCRS matrix for the Jacobian of the residual function
544 createMatrix_();
545
546 // initialize the Jacobian matrix and the vector for the residual function
547 residual_.resize(model_().numTotalDof());
548 resetSystem_();
549
550 // initialize the sparse tables for Flows and Flores
551 createFlows_();
552 }
553
554 // Construct the BCRS matrix for the Jacobian of the residual function
555 void createMatrix_()
556 {
557 OPM_TIMEBLOCK(createMatrix);
558 if (!neighborInfo_.empty()) {
559 // It is ok to call this function multiple times, but it
560 // should not do anything if already called.
561 return;
562 }
563 const auto& model = model_();
564 Stencil stencil(gridView_(), model_().dofMapper());
565
566 // for the main model, find out the global indices of the neighboring degrees of
567 // freedom of each primary degree of freedom
568 using NeighborSet = std::set<unsigned>;
569 std::vector<NeighborSet> sparsityPattern(model.numTotalDof());
570 const Scalar gravity = problem_().gravity()[dimWorld - 1];
571 unsigned numCells = model.numTotalDof();
572 neighborInfo_.reserve(numCells, 6 * numCells);
573 std::vector<NeighborInfoCPU> loc_nbinfo;
574 for (const auto& elem : elements(gridView_())) {
575 stencil.update(elem);
576
577 for (unsigned primaryDofIdx = 0; primaryDofIdx < stencil.numPrimaryDof(); ++primaryDofIdx) {
578 const unsigned myIdx = stencil.globalSpaceIndex(primaryDofIdx);
579 loc_nbinfo.resize(stencil.numDof() - 1); // Do not include the primary dof in neighborInfo_
580
581 for (unsigned dofIdx = 0; dofIdx < stencil.numDof(); ++dofIdx) {
582 const unsigned neighborIdx = stencil.globalSpaceIndex(dofIdx);
583 sparsityPattern[myIdx].insert(neighborIdx);
584 if (dofIdx > 0) {
585 const Scalar trans = problem_().transmissibility(myIdx, neighborIdx);
586 const auto scvfIdx = dofIdx - 1;
587 const auto& scvf = stencil.interiorFace(scvfIdx);
588 const Scalar area = scvf.area();
589 const Scalar Vin = problem_().model().dofTotalVolume(myIdx);
590 const Scalar Vex = problem_().model().dofTotalVolume(neighborIdx);
591 const Scalar zIn = problem_().dofCenterDepth(myIdx);
592 const Scalar zEx = problem_().dofCenterDepth(neighborIdx);
593 const Scalar dZg = (zIn - zEx)*gravity;
594 const Scalar thpres = problem_().thresholdPressure(myIdx, neighborIdx);
595 const auto dirId = scvf.dirId();
596 auto faceDir = dirId < 0 ? FaceDir::DirEnum::Unknown
597 : FaceDir::FromIntersectionIndex(dirId);
598 ResidualNBInfo nbinfo{trans, area, thpres, dZg, faceDir, Vin, Vex, {}, {}, {}, {}};
599 if constexpr (enableFullyImplicitThermal) {
600 nbinfo.inAlpha = problem_().thermalHalfTransmissibility(myIdx, neighborIdx);
601 nbinfo.outAlpha = problem_().thermalHalfTransmissibility(neighborIdx, myIdx);
602 }
603 if constexpr (enableDiffusion) {
604 nbinfo.diffusivity = problem_().diffusivity(myIdx, neighborIdx);
605 }
606 if constexpr (enableDispersion) {
607 nbinfo.dispersivity = problem_().dispersivity(myIdx, neighborIdx);
608 }
609 loc_nbinfo[dofIdx - 1] = NeighborInfoCPU{neighborIdx, nbinfo, nullptr};
610 }
611 }
612 neighborInfo_.appendRow(loc_nbinfo.begin(), loc_nbinfo.end());
613 if (problem_().nonTrivialBoundaryConditions()) {
614 for (unsigned bfIndex = 0; bfIndex < stencil.numBoundaryFaces(); ++bfIndex) {
615 const auto& bf = stencil.boundaryFace(bfIndex);
616 const int dir_id = bf.dirId();
617 // not for NNCs
618 if (dir_id < 0) {
619 continue;
620 }
621 const auto [type, massrateAD] = problem_().boundaryCondition(myIdx, dir_id);
622 // Strip the unnecessary (and zero anyway) derivatives off massrate.
623 VectorBlock massrate(0.0);
624 for (std::size_t ii = 0; ii < massrate.size(); ++ii) {
625 massrate[ii] = massrateAD[ii].value();
626 }
627 const auto& exFluidState = problem_().boundaryFluidState(myIdx, dir_id);
628 BoundaryConditionData bcdata{type,
629 massrate,
630 exFluidState.pvtRegionIndex(),
631 bfIndex,
632 bf.area(),
633 bf.integrationPos()[dimWorld - 1],
634 exFluidState};
635 boundaryInfo_.push_back({myIdx, dir_id, bfIndex, bcdata});
636 }
637 }
638 }
639 }
640
641 // add the additional neighbors and degrees of freedom caused by the auxiliary
642 // equations
643 const std::size_t numAuxMod = model.numAuxiliaryModules();
644 for (unsigned auxModIdx = 0; auxModIdx < numAuxMod; ++auxModIdx) {
645 model.auxiliaryModule(auxModIdx)->addNeighbors(sparsityPattern);
646 }
647
648 // allocate raw matrix
649 jacobian_ = std::make_unique<SparseMatrixAdapter>(simulator_());
650 diagMatAddress_.resize(numCells);
651 // create matrix structure based on sparsity pattern
652 jacobian_->reserve(sparsityPattern);
653 for (unsigned globI = 0; globI < numCells; globI++) {
654 const auto& nbInfos = neighborInfo_[globI];
655 diagMatAddress_[globI] = jacobian_->blockAddress(globI, globI);
656 for (auto& nbInfo : nbInfos) {
657 nbInfo.matBlockAddress = jacobian_->blockAddress(nbInfo.neighbor, globI);
658 }
659 }
660
661 // Create dummy full domain.
662 fullDomain_.cells.resize(numCells);
663 std::iota(fullDomain_.cells.begin(), fullDomain_.cells.end(), 0);
664 }
665
666 // reset the global linear system of equations.
667 void resetSystem_()
668 {
669 residual_ = 0.0;
670 // zero all matrix entries
671 jacobian_->clear();
672 }
673
674 // Initialize the flows, flores, and velocity sparse tables
675 void createFlows_()
676 {
677 OPM_TIMEBLOCK(createFlows);
678 // If FLOWS/FLORES is set in any RPTRST in the schedule, then we initializate the sparse tables.
679 // If DISPERC is in the deck, we initialize the sparse table here as well.
680 const bool anyFlows = simulator_().problem().eclWriter().outputModule().getFlows().anyFlows();
681 const auto& blockFlows = simulator_().problem().eclWriter().outputModule().getFlows().blockFlows();
682 const auto& blockVelocity = simulator_().problem().eclWriter().outputModule().getFlows().blockVelocity();
683 const bool isTemp = simulator_().vanguard().eclState().getSimulationConfig().isTemp();
684 const bool anyFlores = simulator_().problem().eclWriter().outputModule().getFlows().anyFlores() || isTemp;
685 const bool dispersionActive = simulator_().vanguard().eclState().getSimulationConfig().rock_config().dispersion();
686 if (!dispersionActive && !enableBioeffects && blockVelocity.empty()
687 && !((anyFlows || !blockFlows.empty()) && flowsInfo_.empty())
688 && !(anyFlores && floresInfo_.empty())) {
689 return;
690 }
691 const auto& model = model_();
692 const auto& nncOutput = simulator_().problem().eclWriter().getOutputNnc().front();
693 Stencil stencil(gridView_(), model_().dofMapper());
694 const unsigned numCells = model.numTotalDof();
695 std::unordered_multimap<int, std::pair<int, int>> nncIndices;
696 std::vector<FlowInfo> loc_flinfo;
697 std::vector<VelocityInfo> loc_vlinfo;
698 unsigned int nncId = 0;
699 VectorBlock flow(0.0);
700
701 // Create a nnc structure to use fast lookup
702 for (unsigned nncIdx = 0; nncIdx < nncOutput.size(); ++nncIdx) {
703 const int ci1 = nncOutput[nncIdx].cell1;
704 const int ci2 = nncOutput[nncIdx].cell2;
705 nncIndices.emplace(ci1, std::make_pair(ci2, nncIdx));
706 }
707
708 if (anyFlows) {
709 flowsInfo_.reserve(numCells, 6 * numCells);
710 }
711 else if (!blockFlows.empty()) {
712 flowsInfo_.reserve(numCells, 6 * blockFlows.size());
713 }
714 if (anyFlores) {
715 floresInfo_.reserve(numCells, 6 * numCells);
716 }
717 if (dispersionActive || enableBioeffects) {
718 velocityInfo_.reserve(numCells, 6 * numCells);
719 }
720 else if (!blockVelocity.empty()) {
721 velocityInfo_.reserve(numCells, 6 * blockVelocity.size());
722 }
723
724 for (const auto& elem : elements(gridView_())) {
725 stencil.update(elem);
726 for (unsigned primaryDofIdx = 0; primaryDofIdx < stencil.numPrimaryDof(); ++primaryDofIdx) {
727 const unsigned myIdx = stencil.globalSpaceIndex(primaryDofIdx);
728 bool blockFlowFound = false;
729 bool blockVelocityFound = false;
730 if (!blockFlows.empty()) {
731 if (std::ranges::binary_search(blockFlows,
732 simulator_().vanguard().cartesianIndex(myIdx))) {
733 blockFlowFound = true;
734 }
735 else {
736 flowsInfo_.appendRow(loc_flinfo.begin(), loc_flinfo.begin());
737 if (!dispersionActive && !enableBioeffects && !anyFlores && blockVelocity.empty()) {
738 continue;
739 }
740 }
741 }
742 if (!blockVelocity.empty() && !(dispersionActive || enableBioeffects)) {
743 if (std::ranges::binary_search(blockVelocity,
744 simulator_().vanguard().cartesianIndex(myIdx))) {
745 blockVelocityFound = true;
746 }
747 else {
748 velocityInfo_.appendRow(loc_vlinfo.begin(), loc_vlinfo.begin());
749 if (!anyFlows && blockFlows.empty() && !anyFlores) {
750 continue;
751 }
752 }
753 }
754 const int numFaces = stencil.numBoundaryFaces() + stencil.numInteriorFaces();
755 loc_flinfo.resize(numFaces);
756 loc_vlinfo.resize(stencil.numDof() - 1);
757
758 for (unsigned dofIdx = 0; dofIdx < stencil.numDof(); ++dofIdx) {
759 const unsigned neighborIdx = stencil.globalSpaceIndex(dofIdx);
760 if (dofIdx > 0) {
761 const auto scvfIdx = dofIdx - 1;
762 const auto& scvf = stencil.interiorFace(scvfIdx);
763 int faceId = scvf.dirId();
764 const int cartMyIdx = simulator_().vanguard().cartesianIndex(myIdx);
765 const int cartNeighborIdx = simulator_().vanguard().cartesianIndex(neighborIdx);
766 const auto& range = nncIndices.equal_range(cartMyIdx);
767 for (auto it = range.first; it != range.second; ++it) {
768 if (it->second.first == cartNeighborIdx){
769 // -1 gives problem since is used for the nncInput from the deck
770 faceId = -2;
771 // the index is stored to be used for writting the outputs
772 nncId = it->second.second;
773 }
774 }
775 loc_flinfo[dofIdx - 1] = FlowInfo{faceId, flow, nncId};
776 loc_vlinfo[dofIdx - 1] = VelocityInfo{faceId, flow};
777 }
778 }
779
780 for (unsigned bdfIdx = 0; bdfIdx < stencil.numBoundaryFaces(); ++bdfIdx) {
781 const auto& scvf = stencil.boundaryFace(bdfIdx);
782 const int faceId = scvf.dirId();
783 loc_flinfo[stencil.numInteriorFaces() + bdfIdx] = FlowInfo{faceId, flow, nncId};
784 }
785
786 if (anyFlows || blockFlowFound) {
787 flowsInfo_.appendRow(loc_flinfo.begin(), loc_flinfo.end());
788 }
789 if (anyFlores) {
790 floresInfo_.appendRow(loc_flinfo.begin(), loc_flinfo.end());
791 }
792 if (dispersionActive || enableBioeffects || blockVelocityFound) {
793 velocityInfo_.appendRow(loc_vlinfo.begin(), loc_vlinfo.end());
794 }
795 }
796 }
797 }
798
799public:
800 void setResAndJacobi(VectorBlock& res, MatrixBlock& bMat, const ADVectorBlock& resid) const
801 {
802 for (unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
803 res[eqIdx] = resid[eqIdx].value();
804 }
805
806 for (unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
807 for (unsigned pvIdx = 0; pvIdx < numEq; ++pvIdx) {
808 // A[dofIdx][focusDofIdx][eqIdx][pvIdx] is the partial derivative of
809 // the residual function 'eqIdx' for the degree of freedom 'dofIdx'
810 // with regard to the focus variable 'pvIdx' of the degree of freedom
811 // 'focusDofIdx'
812 bMat[eqIdx][pvIdx] = resid[eqIdx].derivative(pvIdx);
813 }
814 }
815 }
816
818 {
819 OPM_TIMEBLOCK(updateFlows);
820 const bool enableFlows = simulator_().problem().eclWriter().outputModule().getFlows().hasFlows();
821 const auto& blockFlows = simulator_().problem().eclWriter().outputModule().getFlows().blockFlows();
822 // We reuse the fluxes in the TEMP option
823 const bool isTemp = simulator_().vanguard().eclState().getSimulationConfig().isTemp();
824 const bool enableFlores = simulator_().problem().eclWriter().outputModule().getFlows().hasFlores() || isTemp;
825 if (!enableFlows && !enableFlores && blockFlows.empty()) {
826 return;
827 }
828 const unsigned int numCells = model_().numTotalDof();
829#ifdef _OPENMP
830#pragma omp parallel for
831#endif
832 for (unsigned globI = 0; globI < numCells; ++globI) {
833 OPM_TIMEBLOCK_LOCAL(linearizationForEachCell, Subsystem::Assembly);
834 const auto& nbInfos = neighborInfo_[globI];
835 ADVectorBlock adres(0.0);
836 ADVectorBlock darcyFlux(0.0);
837 const IntensiveQuantities& intQuantsIn = model_().intensiveQuantities(globI, /*timeIdx*/ 0);
838 // Flux term.
839 {
840 OPM_TIMEBLOCK_LOCAL(fluxCalculationForEachCell, Subsystem::Assembly);
841 short loc = 0;
842 for (const auto& nbInfo : nbInfos) {
843 OPM_TIMEBLOCK_LOCAL(fluxCalculationForEachFace, Subsystem::Assembly);
844 const unsigned globJ = nbInfo.neighbor;
845 assert(globJ != globI);
846 adres = 0.0;
847 darcyFlux = 0.0;
848 const IntensiveQuantities& intQuantsEx = model_().intensiveQuantities(globJ, /*timeIdx*/ 0);
849 LocalResidual::computeFlux(adres, darcyFlux, globI, globJ, intQuantsIn,
850 intQuantsEx, nbInfo.res_nbinfo, problem_().moduleParams());
851 adres *= nbInfo.res_nbinfo.faceArea;
852 if (!blockFlows.empty()) {
853 if (std::ranges::binary_search(blockFlows,
854 simulator_().vanguard().cartesianIndex(globI))) {
855 for (unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
856 flowsInfo_[globI][loc].flow[eqIdx] = adres[eqIdx].value();
857 }
858 }
859 }
860 else if (enableFlows) {
861 for (unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
862 flowsInfo_[globI][loc].flow[eqIdx] = adres[eqIdx].value();
863 }
864 }
865 if (enableFlores) {
866 for (unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
867 floresInfo_[globI][loc].flow[eqIdx] = darcyFlux[eqIdx].value();
868 }
869 }
870 ++loc;
871 }
872 }
873 }
874
875 // Boundary terms. Only looping over cells with nontrivial bcs.
876 for (const auto& bdyInfo : boundaryInfo_) {
877 if (bdyInfo.bcdata.type == BCType::NONE) {
878 continue;
879 }
880
881 ADVectorBlock adres(0.0);
882 const unsigned globI = bdyInfo.cell;
883 const auto& nbInfos = neighborInfo_[globI];
884 const IntensiveQuantities& insideIntQuants = model_().intensiveQuantities(globI, /*timeIdx*/ 0);
885 LocalResidual::computeBoundaryFlux(adres, problem_(), bdyInfo.bcdata, insideIntQuants, globI);
886 adres *= bdyInfo.bcdata.faceArea;
887 const unsigned bfIndex = bdyInfo.bfIndex;
888 if (enableFlows) {
889 for (unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
890 flowsInfo_[globI][nbInfos.size() + bfIndex].flow[eqIdx] = adres[eqIdx].value();
891 }
892 }
893 // TODO also store Flores?
894 }
895 }
896
897private:
898 template <class SubDomainType>
899 void linearize_(const SubDomainType& domain)
900 {
901 // This check should be removed once this is addressed by
902 // for example storing the previous timesteps' values for
903 // rsmax (for DRSDT) and similar.
904 if (!problem_().recycleFirstIterationStorage()) {
905 if (!model_().storeIntensiveQuantities() && !model_().enableStorageCache()) {
906 OPM_THROW(std::runtime_error, "Must have cached either IQs or storage when we cannot recycle.");
907 }
908 }
909
910 OPM_TIMEBLOCK(linearize);
911
912 // We do not call resetSystem_() here, since that will set
913 // the full system to zero, not just our part.
914 // Instead, that must be called before starting the linearization.
915 const bool dispersionActive = simulator_().vanguard().eclState().getSimulationConfig().rock_config().dispersion();
916 const auto& blockVelocity = simulator_().problem().eclWriter().outputModule().getFlows().blockVelocity();
917 const unsigned int numCells = domain.cells.size();
918
919 // Fetch timestepsize used later in accumulation term.
920 const double dt = simulator_().timeStepSize();
921
922#ifdef _OPENMP
923#pragma omp parallel for
924#endif
925 for (unsigned ii = 0; ii < numCells; ++ii) {
926 OPM_TIMEBLOCK_LOCAL(linearizationForEachCell, Subsystem::Assembly);
927 const unsigned globI = domain.cells[ii];
928 const auto& nbInfos = neighborInfo_[globI];
929 VectorBlock res(0.0);
930 MatrixBlock bMat(0.0);
931 ADVectorBlock adres(0.0);
932 ADVectorBlock darcyFlux(0.0);
933 const IntensiveQuantities& intQuantsIn = model_().intensiveQuantities(globI, /*timeIdx*/ 0);
934
935 // Flux term.
936 {
937 OPM_TIMEBLOCK_LOCAL(fluxCalculationForEachCell, Subsystem::Assembly);
938 short loc = 0;
939 for (const auto& nbInfo : nbInfos) {
940 OPM_TIMEBLOCK_LOCAL(fluxCalculationForEachFace, Subsystem::Assembly);
941 const unsigned globJ = nbInfo.neighbor;
942 assert(globJ != globI);
943 res = 0.0;
944 bMat = 0.0;
945 adres = 0.0;
946 darcyFlux = 0.0;
947 const IntensiveQuantities& intQuantsEx = model_().intensiveQuantities(globJ, /*timeIdx*/ 0);
948 LocalResidual::computeFlux(adres, darcyFlux, globI, globJ, intQuantsIn, intQuantsEx,
949 nbInfo.res_nbinfo, problem_().moduleParams());
950 adres *= nbInfo.res_nbinfo.faceArea;
951 if (dispersionActive || enableBioeffects) {
952 for (unsigned phaseIdx = 0; phaseIdx < numEq; ++phaseIdx) {
953 velocityInfo_[globI][loc].velocity[phaseIdx] =
954 darcyFlux[phaseIdx].value() / nbInfo.res_nbinfo.faceArea;
955 }
956 }
957 else if (!blockVelocity.empty()) {
958 if (std::ranges::binary_search(blockVelocity,
959 simulator_().vanguard().cartesianIndex(globI))) {
960 for (unsigned phaseIdx = 0; phaseIdx < numEq; ++phaseIdx) {
961 velocityInfo_[globI][loc].velocity[phaseIdx] =
962 darcyFlux[phaseIdx].value() / nbInfo.res_nbinfo.faceArea;
963 }
964 }
965 }
966 setResAndJacobi(res, bMat, adres);
967 residual_[globI] += res;
968 //SparseAdapter syntax: jacobian_->addToBlock(globI, globI, bMat);
969 *diagMatAddress_[globI] += bMat;
970 bMat *= -1.0;
971 //SparseAdapter syntax: jacobian_->addToBlock(globJ, globI, bMat);
972 *nbInfo.matBlockAddress += bMat;
973 ++loc;
974 }
975 }
976
977 // Accumulation term.
978 const double volume = model_().dofTotalVolume(globI);
979 const Scalar storefac = volume / dt;
980 adres = 0.0;
981 {
982 OPM_TIMEBLOCK_LOCAL(computeStorage, Subsystem::Assembly);
983 LocalResidual::template computeStorage<Evaluation>(adres, intQuantsIn);
984 }
985 setResAndJacobi(res, bMat, adres);
986 // Either use cached storage term, or compute it on the fly.
987 if (model_().enableStorageCache()) {
988 // The cached storage for timeIdx 0 (current time) is not
989 // used, but after storage cache is shifted at the end of the
990 // timestep, it will become cached storage for timeIdx 1.
991 model_().updateCachedStorage(globI, /*timeIdx=*/0, res);
992 // We should not update the storage cache here for NLDD local solves.
993 // This will reset the start-of-step storage to incorrect numbers when
994 // we do local solves, where the iteration number will start from 0,
995 // but the starting state may not be identical to the start-of-step state.
996 // Note that a full assembly must be done before local solves
997 // otherwise this will be left un-updated.
998 if (problem_().iterationContext().isFirstGlobalIteration()) {
999 // Need to update the storage cache.
1000 if (problem_().recycleFirstIterationStorage()) {
1001 // Assumes nothing have changed in the system which
1002 // affects masses calculated from primary variables.
1003 model_().updateCachedStorage(globI, /*timeIdx=*/1, res);
1004 }
1005 else {
1006 Dune::FieldVector<Scalar, numEq> tmp;
1007 const IntensiveQuantities intQuantOld = model_().intensiveQuantities(globI, 1);
1008 LocalResidual::template computeStorage<Scalar>(tmp, intQuantOld);
1009 model_().updateCachedStorage(globI, /*timeIdx=*/1, tmp);
1010 }
1011 }
1012 res -= model_().cachedStorage(globI, 1);
1013 }
1014 else {
1015 OPM_TIMEBLOCK_LOCAL(computeStorage0, Subsystem::Assembly);
1016 Dune::FieldVector<Scalar, numEq> tmp;
1017 const IntensiveQuantities intQuantOld = model_().intensiveQuantities(globI, 1);
1018 LocalResidual::template computeStorage<Scalar>(tmp, intQuantOld);
1019 // assume volume do not change
1020 res -= tmp;
1021 }
1022 res *= storefac;
1023 bMat *= storefac;
1024 residual_[globI] += res;
1025 //SparseAdapter syntax: jacobian_->addToBlock(globI, globI, bMat);
1026 *diagMatAddress_[globI] += bMat;
1027
1028 // Cell-wise source terms.
1029 // This will include well sources if SeparateSparseSourceTerms is false.
1030 res = 0.0;
1031 bMat = 0.0;
1032 adres = 0.0;
1033 if (separateSparseSourceTerms_) {
1034 LocalResidual::computeSourceDense(adres, problem_(), intQuantsIn, globI, 0);
1035 }
1036 else {
1037 LocalResidual::computeSource(adres, problem_(), intQuantsIn, globI, 0);
1038 }
1039 adres *= -volume;
1040 setResAndJacobi(res, bMat, adres);
1041 residual_[globI] += res;
1042 //SparseAdapter syntax: jacobian_->addToBlock(globI, globI, bMat);
1043 *diagMatAddress_[globI] += bMat;
1044 } // end of loop for cell globI.
1045
1046 // Add sparse source terms. For now only wells.
1047 if (separateSparseSourceTerms_) {
1048 problem_().wellModel().addReservoirSourceTerms(residual_, diagMatAddress_);
1049 }
1050
1051 // Boundary terms. Only looping over cells with nontrivial bcs.
1052 for (const auto& bdyInfo : boundaryInfo_) {
1053 if (bdyInfo.bcdata.type == BCType::NONE) {
1054 continue;
1055 }
1056
1057 VectorBlock res(0.0);
1058 MatrixBlock bMat(0.0);
1059 ADVectorBlock adres(0.0);
1060 const unsigned globI = bdyInfo.cell;
1061 const IntensiveQuantities& insideIntQuants = model_().intensiveQuantities(globI, /*timeIdx*/ 0);
1062 LocalResidual::computeBoundaryFlux(adres, problem_(), bdyInfo.bcdata, insideIntQuants, globI);
1063 adres *= bdyInfo.bcdata.faceArea;
1064 setResAndJacobi(res, bMat, adres);
1065 residual_[globI] += res;
1067 *diagMatAddress_[globI] += bMat;
1068 }
1069 }
1070
1071 void updateStoredTransmissibilities()
1072 {
1073 if (neighborInfo_.empty()) {
1074 // This function was called before createMatrix_() was called.
1075 // We call initFirstIteration_(), not createMatrix_(), because
1076 // that will also initialize the residual consistently.
1077 initFirstIteration_();
1078 }
1079
1080 const unsigned numCells = model_().numTotalDof();
1081#ifdef _OPENMP
1082#pragma omp parallel for
1083#endif
1084 for (unsigned globI = 0; globI < numCells; globI++) {
1085 auto nbInfos = neighborInfo_[globI]; // nbInfos will be a SparseTable<...>::mutable_iterator_range.
1086 for (auto& nbInfo : nbInfos) {
1087 const unsigned globJ = nbInfo.neighbor;
1088 nbInfo.res_nbinfo.trans = problem_().transmissibility(globI, globJ);
1089 }
1090 }
1091 }
1092
1093 Simulator* simulatorPtr_{};
1094
1095 // the jacobian matrix
1096 std::unique_ptr<SparseMatrixAdapter> jacobian_{};
1097
1098 // the right-hand side
1099 GlobalEqVector residual_;
1100
1101 LinearizationType linearizationType_{};
1102
1103 using ResidualNBInfo = typename LocalResidual::ResidualNBInfo;
1104 using NeighborInfoCPU = NeighborInfoStruct<ResidualNBInfo, MatrixBlock>;
1105
1106 SparseTable<NeighborInfoCPU> neighborInfo_{};
1107 std::vector<MatrixBlock*> diagMatAddress_{};
1108
1109 struct FlowInfo
1110 {
1111 int faceId;
1112 VectorBlock flow;
1113 unsigned int nncId;
1114 };
1115 SparseTable<FlowInfo> flowsInfo_;
1116 SparseTable<FlowInfo> floresInfo_;
1117
1118 struct VelocityInfo
1119 {
1120 int faceId;
1121 VectorBlock velocity;
1122 };
1123 SparseTable<VelocityInfo> velocityInfo_;
1124
1125 using ScalarFluidState = typename IntensiveQuantities::ScalarFluidState;
1126 struct BoundaryConditionData
1127 {
1128 BCType type;
1129 VectorBlock massRate;
1130 unsigned pvtRegionIdx;
1131 unsigned boundaryFaceIndex;
1132 double faceArea;
1133 double faceZCoord;
1134 ScalarFluidState exFluidState;
1135 };
1136
1137 struct BoundaryInfo
1138 {
1139 unsigned int cell;
1140 int dir;
1141 unsigned int bfIndex;
1142 BoundaryConditionData bcdata;
1143 };
1144 std::vector<BoundaryInfo> boundaryInfo_;
1145
1146 bool separateSparseSourceTerms_ = false;
1147
1148 FullDomain<> fullDomain_;
1149
1150 int exportIndex_;
1151 int exportCount_;
1152};
1153} // namespace Opm
1154
1155#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:458
const LinearizationType & getLinearizationType() const
Definition: tpfalinearizer.hh:442
const auto & getNeighborInfo() const
Definition: tpfalinearizer.hh:470
void updateBoundaryConditionData()
Definition: tpfalinearizer.hh:480
void linearize()
Linearize the full system of non-linear equations.
Definition: tpfalinearizer.hh:291
SparseMatrixAdapter & jacobian()
Definition: tpfalinearizer.hh:409
const auto & getFlowsInfo() const
Return constant reference to the flowsInfo.
Definition: tpfalinearizer.hh:450
std::map< unsigned, Constraints > constraintsMap() const
Returns the map of constraint degrees of freedom.
Definition: tpfalinearizer.hh:504
TpfaLinearizer()
Definition: tpfalinearizer.hh:238
void linearizeDomain(const SubDomainType &domain)
Linearize the part of the non-linear system of equations that is associated with a part of the spatia...
Definition: tpfalinearizer.hh:346
void finalize()
Definition: tpfalinearizer.hh:367
void init(Simulator &simulator)
Initialize the linearizer.
Definition: tpfalinearizer.hh:264
void updateFlowsInfo()
Definition: tpfalinearizer.hh:817
void setResAndJacobi(VectorBlock &res, MatrixBlock &bMat, const ADVectorBlock &resid) const
Definition: tpfalinearizer.hh:800
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:424
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:406
void setLinearizationType(LinearizationType linearizationType)
Definition: tpfalinearizer.hh:439
GlobalEqVector & residual()
Definition: tpfalinearizer.hh:418
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:467
void updateDiscretizationParameters()
Definition: tpfalinearizer.hh:475
void resetSystem_(const SubDomainType &domain)
Definition: tpfalinearizer.hh:508
void linearizeAuxiliaryEquations()
Linearize the part of the non-linear system of equations that is associated with the spatial domain.
Definition: tpfalinearizer.hh:374
const GlobalEqVector & residual() const
Return constant reference to global residual vector.
Definition: tpfalinearizer.hh:415
A small fixed-size square matrix class for use in CUDA kernels.
Definition: MiniMatrix.hpp:37
Definition: MiniVector.hpp:50
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