28#ifndef TPFA_LINEARIZER_HH
29#define TPFA_LINEARIZER_HH
31#include <dune/common/version.hh>
32#include <dune/common/fvector.hh>
33#include <dune/common/fmatrix.hh>
35#include <opm/common/Exceptions.hpp>
36#include <opm/common/TimingMacros.hpp>
38#include <opm/grid/utility/SparseTable.hpp>
40#include <opm/material/common/ConditionalStorage.hpp>
42#include <opm/input/eclipse/EclipseState/Grid/FaceDir.hpp>
43#include <opm/input/eclipse/Schedule/BCProp.hpp>
61#include <unordered_map>
64#include <fmt/format.h>
66#include <opm/common/utility/gpuDecorators.hpp>
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>
90template<
class TypeTag>
94template<
class Storage = std::vector<
int>>
100#if HAVE_CUDA && OPM_IS_COMPILING_WITH_GPU_COMPILER
103 if (CPUDomain.
cells.size() == 0) {
104 OPM_THROW(std::runtime_error,
"Cannot copy empty full domain to GPU.");
106 return FullDomain<gpuistl::GpuBuffer<int>>{
107 gpuistl::GpuBuffer<int>(CPUDomain.
cells)
111 FullDomain<gpuistl::GpuView<int>>
make_view(FullDomain<gpuistl::GpuBuffer<int>>& buffer)
113 if (buffer.cells.size() == 0) {
114 OPM_THROW(std::runtime_error,
"Cannot make view of empty full domain buffer.");
116 return FullDomain<gpuistl::GpuView<int>>{
122template <
class Res
idualNBInfoType,
class BlockType>
129 template <
class OtherBlockType>
140 template <
class PtrType>
157#if HAVE_CUDA && OPM_IS_COMPILING_WITH_GPU_COMPILER
159 template<
class MiniMatrixType,
class MatrixBlockType,
class Res
idualNBInfoType>
160 auto copy_to_gpu(
const SparseTable<NeighborInfoStruct<ResidualNBInfoType, MatrixBlockType>>& cpuNeighborInfoTable)
163 using StructWithMinimatrix = NeighborInfoStruct<ResidualNBInfoType, MiniMatrixType>;
164 std::vector<StructWithMinimatrix> minimatrices(cpuNeighborInfoTable.dataSize());
166 for (
auto e : cpuNeighborInfoTable.dataStorage()) {
167 minimatrices[idx++] = StructWithMinimatrix(e);
170 return SparseTable<StructWithMinimatrix, gpuistl::GpuBuffer>(
171 gpuistl::GpuBuffer<StructWithMinimatrix>(minimatrices),
172 gpuistl::GpuBuffer<int>(cpuNeighborInfoTable.rowStarts())
187template<
class TypeTag>
208 using Element =
typename GridView::template Codim<0>::Entity;
209 using ElementIterator =
typename GridView::template Codim<0>::Iterator;
211 using Vector = GlobalEqVector;
213 enum { numEq = getPropValue<TypeTag, Properties::NumEq>() };
214 enum { historySize = getPropValue<TypeTag, Properties::TimeDiscHistorySize>() };
215 enum { dimWorld = GridView::dimensionworld };
217 using MatrixBlock =
typename SparseMatrixAdapter::MatrixBlock;
218 using VectorBlock = Dune::FieldVector<Scalar, numEq>;
221#if HAVE_CUDA && OPM_IS_COMPILING_WITH_GPU_COMPILER
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>();
240 simulatorPtr_ =
nullptr;
241 separateSparseSourceTerms_ = Parameters::Get<Parameters::SeparateSparseSourceTerms>();
251 Parameters::Register<Parameters::SeparateSparseSourceTerms>
252 (
"Treat well source terms all in one go, instead of on a cell by cell basis.");
266 simulatorPtr_ = &simulator;
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;
321 std::cout <<
"rank " << simulator_().
gridView().comm().rank()
322 <<
" caught an exception while linearizing"
323 <<
"\n" << std::flush;
326 succeeded = simulator_().
gridView().comm().min(succeeded);
329 throw NumericalProblem(
"A process did not succeed in linearizing the system");
345 template <
class SubDomainType>
353 initFirstIteration_();
357 if (problem_().iterationContext().inLocalSolve()) {
368 { jacobian_->finalize(); }
376 OPM_TIMEBLOCK(linearizeAuxilaryEquations);
380 auto& model = model_();
381 const auto& comm = simulator_().
gridView().comm();
382 for (
unsigned auxModIdx = 0; auxModIdx < model.numAuxiliaryModules(); ++auxModIdx) {
383 bool succeeded =
true;
385 model.auxiliaryModule(auxModIdx)->linearize(*jacobian_, residual_);
387 catch (
const std::exception& e) {
390 std::cout <<
"rank " << simulator_().
gridView().comm().rank()
391 <<
" caught an exception while linearizing:" << e.what()
392 <<
"\n" << std::flush;
395 succeeded = comm.min(succeeded);
398 throw NumericalProblem(
"linearization of an auxiliary equation failed");
407 {
return *jacobian_; }
410 {
return *jacobian_; }
416 {
return residual_; }
419 {
return residual_; }
424 void exportSystem(
const int idx, std::string& tag,
const char *path=
"export")
426 const bool export_sparsity = exportIndex_ == -1;
429 exportCount_ = exportIndex_ == idx ? ++exportCount_ : 0;
431 tag = fmt::format(fmt::runtime(
"_{:03d}_{:02d}"), exportIndex_, exportCount_);
433 fmt::print(fmt::runtime(
"index = {:d}\n"), exportIndex_);
434 fmt::print(fmt::runtime(
"count = {:d}\n"), exportCount_);
436 Opm::exportSystem(jacobian_->istlMatrix(), residual_, export_sparsity, tag.c_str(), path);
440 { linearizationType_ = linearizationType; }
443 {
return linearizationType_; }
451 {
return flowsInfo_; }
459 {
return floresInfo_; }
468 {
return velocityInfo_; }
471 return neighborInfo_;
477 updateStoredTransmissibilities();
482 for (
auto& bdyInfo : boundaryInfo_) {
483 const auto [type, massrateAD] = problem_().boundaryCondition(bdyInfo.cell, bdyInfo.dir);
486 VectorBlock massrate(0.0);
487 for (std::size_t ii = 0; ii < massrate.size(); ++ii) {
488 massrate[ii] = massrateAD[ii].value();
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;
507 template <
class SubDomainType>
511 initFirstIteration_();
513 for (
int globI : domain.cells) {
514 residual_[globI] = 0.0;
515 jacobian_->clearRow(globI, 0.0);
521 {
return *simulatorPtr_; }
523 const Simulator& simulator_()
const
524 {
return *simulatorPtr_; }
527 {
return simulator_().
problem(); }
529 const Problem& problem_()
const
530 {
return simulator_().
problem(); }
533 {
return simulator_().
model(); }
535 const Model& model_()
const
536 {
return simulator_().
model(); }
538 const GridView& gridView_()
const
539 {
return problem_().gridView(); }
541 void initFirstIteration_()
547 residual_.resize(model_().numTotalDof());
558 if (!neighborInfo_.empty()) {
563 const auto& model = model_();
564 Stencil stencil(gridView_(), model_().dofMapper());
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);
577 for (
unsigned primaryDofIdx = 0; primaryDofIdx < stencil.numPrimaryDof(); ++primaryDofIdx) {
578 const unsigned myIdx = stencil.globalSpaceIndex(primaryDofIdx);
579 loc_nbinfo.resize(stencil.numDof() - 1);
581 for (
unsigned dofIdx = 0; dofIdx < stencil.numDof(); ++dofIdx) {
582 const unsigned neighborIdx = stencil.globalSpaceIndex(dofIdx);
583 sparsityPattern[myIdx].insert(neighborIdx);
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);
603 if constexpr (enableDiffusion) {
604 nbinfo.diffusivity = problem_().diffusivity(myIdx, neighborIdx);
606 if constexpr (enableDispersion) {
607 nbinfo.dispersivity = problem_().dispersivity(myIdx, neighborIdx);
609 loc_nbinfo[dofIdx - 1] = NeighborInfoCPU{neighborIdx, nbinfo,
nullptr};
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();
621 const auto [type, massrateAD] = problem_().boundaryCondition(myIdx, dir_id);
623 VectorBlock massrate(0.0);
624 for (std::size_t ii = 0; ii < massrate.size(); ++ii) {
625 massrate[ii] = massrateAD[ii].value();
627 const auto& exFluidState = problem_().boundaryFluidState(myIdx, dir_id);
628 BoundaryConditionData bcdata{type,
630 exFluidState.pvtRegionIndex(),
633 bf.integrationPos()[dimWorld - 1],
635 boundaryInfo_.push_back({myIdx, dir_id, bfIndex, bcdata});
643 const std::size_t numAuxMod = model.numAuxiliaryModules();
644 for (
unsigned auxModIdx = 0; auxModIdx < numAuxMod; ++auxModIdx) {
645 model.auxiliaryModule(auxModIdx)->addNeighbors(sparsityPattern);
649 jacobian_ = std::make_unique<SparseMatrixAdapter>(simulator_());
650 diagMatAddress_.resize(numCells);
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);
662 fullDomain_.
cells.resize(numCells);
663 std::iota(fullDomain_.
cells.begin(), fullDomain_.
cells.end(), 0);
677 OPM_TIMEBLOCK(createFlows);
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())) {
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);
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));
709 flowsInfo_.reserve(numCells, 6 * numCells);
711 else if (!blockFlows.empty()) {
712 flowsInfo_.reserve(numCells, 6 * blockFlows.size());
715 floresInfo_.reserve(numCells, 6 * numCells);
717 if (dispersionActive || enableBioeffects) {
718 velocityInfo_.reserve(numCells, 6 * numCells);
720 else if (!blockVelocity.empty()) {
721 velocityInfo_.reserve(numCells, 6 * blockVelocity.size());
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;
736 flowsInfo_.appendRow(loc_flinfo.begin(), loc_flinfo.begin());
737 if (!dispersionActive && !enableBioeffects && !anyFlores && blockVelocity.empty()) {
742 if (!blockVelocity.empty() && !(dispersionActive || enableBioeffects)) {
743 if (std::ranges::binary_search(blockVelocity,
744 simulator_().vanguard().cartesianIndex(myIdx))) {
745 blockVelocityFound =
true;
748 velocityInfo_.appendRow(loc_vlinfo.begin(), loc_vlinfo.begin());
749 if (!anyFlows && blockFlows.empty() && !anyFlores) {
754 const int numFaces = stencil.numBoundaryFaces() + stencil.numInteriorFaces();
755 loc_flinfo.resize(numFaces);
756 loc_vlinfo.resize(stencil.numDof() - 1);
758 for (
unsigned dofIdx = 0; dofIdx < stencil.numDof(); ++dofIdx) {
759 const unsigned neighborIdx = stencil.globalSpaceIndex(dofIdx);
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){
772 nncId = it->second.second;
775 loc_flinfo[dofIdx - 1] = FlowInfo{faceId, flow, nncId};
776 loc_vlinfo[dofIdx - 1] = VelocityInfo{faceId, flow};
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};
786 if (anyFlows || blockFlowFound) {
787 flowsInfo_.appendRow(loc_flinfo.begin(), loc_flinfo.end());
790 floresInfo_.appendRow(loc_flinfo.begin(), loc_flinfo.end());
792 if (dispersionActive || enableBioeffects || blockVelocityFound) {
793 velocityInfo_.appendRow(loc_vlinfo.begin(), loc_vlinfo.end());
802 for (
unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
803 res[eqIdx] = resid[eqIdx].value();
806 for (
unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
807 for (
unsigned pvIdx = 0; pvIdx < numEq; ++pvIdx) {
812 bMat[eqIdx][pvIdx] = resid[eqIdx].derivative(pvIdx);
819 OPM_TIMEBLOCK(updateFlows);
820 const bool enableFlows = simulator_().
problem().eclWriter().outputModule().getFlows().hasFlows();
821 const auto& blockFlows = simulator_().
problem().eclWriter().outputModule().getFlows().blockFlows();
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()) {
828 const unsigned int numCells = model_().numTotalDof();
830#pragma omp parallel for
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, 0);
840 OPM_TIMEBLOCK_LOCAL(fluxCalculationForEachCell, Subsystem::Assembly);
842 for (
const auto& nbInfo : nbInfos) {
843 OPM_TIMEBLOCK_LOCAL(fluxCalculationForEachFace, Subsystem::Assembly);
844 const unsigned globJ = nbInfo.neighbor;
845 assert(globJ != globI);
848 const IntensiveQuantities& intQuantsEx = model_().intensiveQuantities(globJ, 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();
860 else if (enableFlows) {
861 for (
unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
862 flowsInfo_[globI][loc].flow[eqIdx] = adres[eqIdx].value();
866 for (
unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
867 floresInfo_[globI][loc].flow[eqIdx] = darcyFlux[eqIdx].value();
876 for (
const auto& bdyInfo : boundaryInfo_) {
881 ADVectorBlock adres(0.0);
882 const unsigned globI = bdyInfo.cell;
883 const auto& nbInfos = neighborInfo_[globI];
884 const IntensiveQuantities& insideIntQuants = model_().intensiveQuantities(globI, 0);
885 LocalResidual::computeBoundaryFlux(adres, problem_(), bdyInfo.bcdata, insideIntQuants, globI);
886 adres *= bdyInfo.bcdata.faceArea;
887 const unsigned bfIndex = bdyInfo.bfIndex;
889 for (
unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
890 flowsInfo_[globI][nbInfos.size() + bfIndex].flow[eqIdx] = adres[eqIdx].value();
898 template <
class SubDomainType>
899 void linearize_(
const SubDomainType& domain)
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.");
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();
923#pragma omp parallel for
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, 0);
937 OPM_TIMEBLOCK_LOCAL(fluxCalculationForEachCell, Subsystem::Assembly);
939 for (
const auto& nbInfo : nbInfos) {
940 OPM_TIMEBLOCK_LOCAL(fluxCalculationForEachFace, Subsystem::Assembly);
941 const unsigned globJ = nbInfo.neighbor;
942 assert(globJ != globI);
947 const IntensiveQuantities& intQuantsEx = model_().intensiveQuantities(globJ, 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;
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;
967 residual_[globI] += res;
969 *diagMatAddress_[globI] += bMat;
972 *nbInfo.matBlockAddress += bMat;
978 const double volume = model_().dofTotalVolume(globI);
979 const Scalar storefac = volume / dt;
982 OPM_TIMEBLOCK_LOCAL(computeStorage, Subsystem::Assembly);
983 LocalResidual::template computeStorage<Evaluation>(adres, intQuantsIn);
987 if (model_().enableStorageCache()) {
991 model_().updateCachedStorage(globI, 0, res);
998 if (problem_().iterationContext().isFirstGlobalIteration()) {
1000 if (problem_().recycleFirstIterationStorage()) {
1003 model_().updateCachedStorage(globI, 1, res);
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, 1, tmp);
1012 res -= model_().cachedStorage(globI, 1);
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);
1024 residual_[globI] += res;
1026 *diagMatAddress_[globI] += bMat;
1033 if (separateSparseSourceTerms_) {
1034 LocalResidual::computeSourceDense(adres, problem_(), intQuantsIn, globI, 0);
1037 LocalResidual::computeSource(adres, problem_(), intQuantsIn, globI, 0);
1041 residual_[globI] += res;
1043 *diagMatAddress_[globI] += bMat;
1047 if (separateSparseSourceTerms_) {
1048 problem_().wellModel().addReservoirSourceTerms(residual_, diagMatAddress_);
1052 for (
const auto& bdyInfo : boundaryInfo_) {
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, 0);
1062 LocalResidual::computeBoundaryFlux(adres, problem_(), bdyInfo.bcdata, insideIntQuants, globI);
1063 adres *= bdyInfo.bcdata.faceArea;
1065 residual_[globI] += res;
1067 *diagMatAddress_[globI] += bMat;
1071 void updateStoredTransmissibilities()
1073 if (neighborInfo_.empty()) {
1077 initFirstIteration_();
1080 const unsigned numCells = model_().numTotalDof();
1082#pragma omp parallel for
1084 for (
unsigned globI = 0; globI < numCells; globI++) {
1085 auto nbInfos = neighborInfo_[globI];
1086 for (
auto& nbInfo : nbInfos) {
1087 const unsigned globJ = nbInfo.neighbor;
1088 nbInfo.res_nbinfo.trans = problem_().transmissibility(globI, globJ);
1093 Simulator* simulatorPtr_{};
1096 std::unique_ptr<SparseMatrixAdapter> jacobian_{};
1099 GlobalEqVector residual_;
1101 LinearizationType linearizationType_{};
1103 using ResidualNBInfo =
typename LocalResidual::ResidualNBInfo;
1104 using NeighborInfoCPU = NeighborInfoStruct<ResidualNBInfo, MatrixBlock>;
1106 SparseTable<NeighborInfoCPU> neighborInfo_{};
1107 std::vector<MatrixBlock*> diagMatAddress_{};
1115 SparseTable<FlowInfo> flowsInfo_;
1116 SparseTable<FlowInfo> floresInfo_;
1121 VectorBlock velocity;
1123 SparseTable<VelocityInfo> velocityInfo_;
1125 using ScalarFluidState =
typename IntensiveQuantities::ScalarFluidState;
1126 struct BoundaryConditionData
1129 VectorBlock massRate;
1130 unsigned pvtRegionIdx;
1131 unsigned boundaryFaceIndex;
1134 ScalarFluidState exFluidState;
1141 unsigned int bfIndex;
1142 BoundaryConditionData bcdata;
1144 std::vector<BoundaryInfo> boundaryInfo_;
1146 bool separateSparseSourceTerms_ =
false;
1148 FullDomain<> fullDomain_;
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