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");
346 template <
class SubDomainType>
347 void linearizeDomain(
const SubDomainType& domain,
const bool isNlddLocalSolve =
false)
354 initFirstIteration_();
358 if (isNlddLocalSolve) {
365 linearize_(domain, isNlddLocalSolve);
369 { jacobian_->finalize(); }
377 OPM_TIMEBLOCK(linearizeAuxilaryEquations);
381 auto& model = model_();
382 const auto& comm = simulator_().
gridView().comm();
383 for (
unsigned auxModIdx = 0; auxModIdx < model.numAuxiliaryModules(); ++auxModIdx) {
384 bool succeeded =
true;
386 model.auxiliaryModule(auxModIdx)->linearize(*jacobian_, residual_);
388 catch (
const std::exception& e) {
391 std::cout <<
"rank " << simulator_().
gridView().comm().rank()
392 <<
" caught an exception while linearizing:" << e.what()
393 <<
"\n" << std::flush;
396 succeeded = comm.min(succeeded);
399 throw NumericalProblem(
"linearization of an auxiliary equation failed");
408 {
return *jacobian_; }
411 {
return *jacobian_; }
417 {
return residual_; }
420 {
return residual_; }
425 void exportSystem(
const int idx, std::string& tag,
const char *path=
"export")
427 const bool export_sparsity = exportIndex_ == -1;
430 exportCount_ = exportIndex_ == idx ? ++exportCount_ : 0;
432 tag = fmt::format(fmt::runtime(
"_{:03d}_{:02d}"), exportIndex_, exportCount_);
434 fmt::print(fmt::runtime(
"index = {:d}\n"), exportIndex_);
435 fmt::print(fmt::runtime(
"count = {:d}\n"), exportCount_);
437 Opm::exportSystem(jacobian_->istlMatrix(), residual_, export_sparsity, tag.c_str(), path);
441 { linearizationType_ = linearizationType; }
444 {
return linearizationType_; }
452 {
return flowsInfo_; }
460 {
return floresInfo_; }
469 {
return velocityInfo_; }
472 return neighborInfo_;
478 updateStoredTransmissibilities();
483 for (
auto& bdyInfo : boundaryInfo_) {
484 const auto [type, massrateAD] = problem_().boundaryCondition(bdyInfo.cell, bdyInfo.dir);
487 VectorBlock massrate(0.0);
488 for (std::size_t ii = 0; ii < massrate.size(); ++ii) {
489 massrate[ii] = massrateAD[ii].value();
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;
508 template <
class SubDomainType>
512 initFirstIteration_();
514 for (
int globI : domain.cells) {
515 residual_[globI] = 0.0;
516 jacobian_->clearRow(globI, 0.0);
522 {
return *simulatorPtr_; }
524 const Simulator& simulator_()
const
525 {
return *simulatorPtr_; }
528 {
return simulator_().
problem(); }
530 const Problem& problem_()
const
531 {
return simulator_().
problem(); }
534 {
return simulator_().
model(); }
536 const Model& model_()
const
537 {
return simulator_().
model(); }
539 const GridView& gridView_()
const
540 {
return problem_().gridView(); }
542 void initFirstIteration_()
548 residual_.resize(model_().numTotalDof());
559 if (!neighborInfo_.empty()) {
564 const auto& model = model_();
565 Stencil stencil(gridView_(), model_().dofMapper());
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);
578 for (
unsigned primaryDofIdx = 0; primaryDofIdx < stencil.numPrimaryDof(); ++primaryDofIdx) {
579 const unsigned myIdx = stencil.globalSpaceIndex(primaryDofIdx);
580 loc_nbinfo.resize(stencil.numDof() - 1);
582 for (
unsigned dofIdx = 0; dofIdx < stencil.numDof(); ++dofIdx) {
583 const unsigned neighborIdx = stencil.globalSpaceIndex(dofIdx);
584 sparsityPattern[myIdx].insert(neighborIdx);
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);
604 if constexpr (enableDiffusion) {
605 nbinfo.diffusivity = problem_().diffusivity(myIdx, neighborIdx);
607 if constexpr (enableDispersion) {
608 nbinfo.dispersivity = problem_().dispersivity(myIdx, neighborIdx);
610 loc_nbinfo[dofIdx - 1] = NeighborInfoCPU{neighborIdx, nbinfo,
nullptr};
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();
622 const auto [type, massrateAD] = problem_().boundaryCondition(myIdx, dir_id);
624 VectorBlock massrate(0.0);
625 for (std::size_t ii = 0; ii < massrate.size(); ++ii) {
626 massrate[ii] = massrateAD[ii].value();
628 const auto& exFluidState = problem_().boundaryFluidState(myIdx, dir_id);
629 BoundaryConditionData bcdata{type,
631 exFluidState.pvtRegionIndex(),
634 bf.integrationPos()[dimWorld - 1],
636 boundaryInfo_.push_back({myIdx, dir_id, bfIndex, bcdata});
644 const std::size_t numAuxMod = model.numAuxiliaryModules();
645 for (
unsigned auxModIdx = 0; auxModIdx < numAuxMod; ++auxModIdx) {
646 model.auxiliaryModule(auxModIdx)->addNeighbors(sparsityPattern);
650 jacobian_ = std::make_unique<SparseMatrixAdapter>(simulator_());
651 diagMatAddress_.resize(numCells);
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);
663 fullDomain_.
cells.resize(numCells);
664 std::iota(fullDomain_.
cells.begin(), fullDomain_.
cells.end(), 0);
678 OPM_TIMEBLOCK(createFlows);
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())) {
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);
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));
710 flowsInfo_.reserve(numCells, 6 * numCells);
712 else if (!blockFlows.empty()) {
713 flowsInfo_.reserve(numCells, 6 * blockFlows.size());
716 floresInfo_.reserve(numCells, 6 * numCells);
718 if (dispersionActive || enableBioeffects) {
719 velocityInfo_.reserve(numCells, 6 * numCells);
721 else if (!blockVelocity.empty()) {
722 velocityInfo_.reserve(numCells, 6 * blockVelocity.size());
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;
737 flowsInfo_.appendRow(loc_flinfo.begin(), loc_flinfo.begin());
738 if (!dispersionActive && !enableBioeffects && !anyFlores && blockVelocity.empty()) {
743 if (!blockVelocity.empty() && !(dispersionActive || enableBioeffects)) {
744 if (std::ranges::binary_search(blockVelocity,
745 simulator_().vanguard().cartesianIndex(myIdx))) {
746 blockVelocityFound =
true;
749 velocityInfo_.appendRow(loc_vlinfo.begin(), loc_vlinfo.begin());
750 if (!anyFlows && blockFlows.empty() && !anyFlores) {
755 const int numFaces = stencil.numBoundaryFaces() + stencil.numInteriorFaces();
756 loc_flinfo.resize(numFaces);
757 loc_vlinfo.resize(stencil.numDof() - 1);
759 for (
unsigned dofIdx = 0; dofIdx < stencil.numDof(); ++dofIdx) {
760 const unsigned neighborIdx = stencil.globalSpaceIndex(dofIdx);
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){
773 nncId = it->second.second;
776 loc_flinfo[dofIdx - 1] = FlowInfo{faceId, flow, nncId};
777 loc_vlinfo[dofIdx - 1] = VelocityInfo{faceId, flow};
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};
787 if (anyFlows || blockFlowFound) {
788 flowsInfo_.appendRow(loc_flinfo.begin(), loc_flinfo.end());
791 floresInfo_.appendRow(loc_flinfo.begin(), loc_flinfo.end());
793 if (dispersionActive || enableBioeffects || blockVelocityFound) {
794 velocityInfo_.appendRow(loc_vlinfo.begin(), loc_vlinfo.end());
803 for (
unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
804 res[eqIdx] = resid[eqIdx].value();
807 for (
unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
808 for (
unsigned pvIdx = 0; pvIdx < numEq; ++pvIdx) {
813 bMat[eqIdx][pvIdx] = resid[eqIdx].derivative(pvIdx);
820 OPM_TIMEBLOCK(updateFlows);
821 const bool enableFlows = simulator_().
problem().eclWriter().outputModule().getFlows().hasFlows();
822 const auto& blockFlows = simulator_().
problem().eclWriter().outputModule().getFlows().blockFlows();
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()) {
829 const unsigned int numCells = model_().numTotalDof();
831#pragma omp parallel for
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, 0);
841 OPM_TIMEBLOCK_LOCAL(fluxCalculationForEachCell, Subsystem::Assembly);
843 for (
const auto& nbInfo : nbInfos) {
844 OPM_TIMEBLOCK_LOCAL(fluxCalculationForEachFace, Subsystem::Assembly);
845 const unsigned globJ = nbInfo.neighbor;
846 assert(globJ != globI);
849 const IntensiveQuantities& intQuantsEx = model_().intensiveQuantities(globJ, 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();
861 else if (enableFlows) {
862 for (
unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
863 flowsInfo_[globI][loc].flow[eqIdx] = adres[eqIdx].value();
867 for (
unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
868 floresInfo_[globI][loc].flow[eqIdx] = darcyFlux[eqIdx].value();
877 for (
const auto& bdyInfo : boundaryInfo_) {
882 ADVectorBlock adres(0.0);
883 const unsigned globI = bdyInfo.cell;
884 const auto& nbInfos = neighborInfo_[globI];
885 const IntensiveQuantities& insideIntQuants = model_().intensiveQuantities(globI, 0);
886 LocalResidual::computeBoundaryFlux(adres, problem_(), bdyInfo.bcdata, insideIntQuants, globI);
887 adres *= bdyInfo.bcdata.faceArea;
888 const unsigned bfIndex = bdyInfo.bfIndex;
890 for (
unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
891 flowsInfo_[globI][nbInfos.size() + bfIndex].flow[eqIdx] = adres[eqIdx].value();
899 template <
class SubDomainType>
900 void linearize_(
const SubDomainType& domain,
bool isNlddLocalSolve)
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.");
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();
924#pragma omp parallel for
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, 0);
938 OPM_TIMEBLOCK_LOCAL(fluxCalculationForEachCell, Subsystem::Assembly);
940 for (
const auto& nbInfo : nbInfos) {
941 OPM_TIMEBLOCK_LOCAL(fluxCalculationForEachFace, Subsystem::Assembly);
942 const unsigned globJ = nbInfo.neighbor;
943 assert(globJ != globI);
948 const IntensiveQuantities& intQuantsEx = model_().intensiveQuantities(globJ, 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;
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;
968 residual_[globI] += res;
970 *diagMatAddress_[globI] += bMat;
973 *nbInfo.matBlockAddress += bMat;
979 const double volume = model_().dofTotalVolume(globI);
980 const Scalar storefac = volume / dt;
983 OPM_TIMEBLOCK_LOCAL(computeStorage, Subsystem::Assembly);
984 LocalResidual::template computeStorage<Evaluation>(adres, intQuantsIn);
988 if (model_().enableStorageCache()) {
992 model_().updateCachedStorage(globI, 0, res);
999 if (model_().newtonMethod().numIterations() == 0 && !isNlddLocalSolve) {
1001 if (problem_().recycleFirstIterationStorage()) {
1004 model_().updateCachedStorage(globI, 1, res);
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, 1, tmp);
1013 res -= model_().cachedStorage(globI, 1);
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);
1025 residual_[globI] += res;
1027 *diagMatAddress_[globI] += bMat;
1034 if (separateSparseSourceTerms_) {
1035 LocalResidual::computeSourceDense(adres, problem_(), intQuantsIn, globI, 0);
1038 LocalResidual::computeSource(adres, problem_(), intQuantsIn, globI, 0);
1042 residual_[globI] += res;
1044 *diagMatAddress_[globI] += bMat;
1048 if (separateSparseSourceTerms_) {
1049 problem_().wellModel().addReservoirSourceTerms(residual_, diagMatAddress_);
1053 for (
const auto& bdyInfo : boundaryInfo_) {
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, 0);
1063 LocalResidual::computeBoundaryFlux(adres, problem_(), bdyInfo.bcdata, insideIntQuants, globI);
1064 adres *= bdyInfo.bcdata.faceArea;
1066 residual_[globI] += res;
1068 *diagMatAddress_[globI] += bMat;
1072 void updateStoredTransmissibilities()
1074 if (neighborInfo_.empty()) {
1078 initFirstIteration_();
1081 const unsigned numCells = model_().numTotalDof();
1083#pragma omp parallel for
1085 for (
unsigned globI = 0; globI < numCells; globI++) {
1086 auto nbInfos = neighborInfo_[globI];
1087 for (
auto& nbInfo : nbInfos) {
1088 const unsigned globJ = nbInfo.neighbor;
1089 nbInfo.res_nbinfo.trans = problem_().transmissibility(globI, globJ);
1094 Simulator* simulatorPtr_{};
1097 std::unique_ptr<SparseMatrixAdapter> jacobian_{};
1100 GlobalEqVector residual_;
1102 LinearizationType linearizationType_{};
1104 using ResidualNBInfo =
typename LocalResidual::ResidualNBInfo;
1105 using NeighborInfoCPU = NeighborInfoStruct<ResidualNBInfo, MatrixBlock>;
1107 SparseTable<NeighborInfoCPU> neighborInfo_{};
1108 std::vector<MatrixBlock*> diagMatAddress_{};
1116 SparseTable<FlowInfo> flowsInfo_;
1117 SparseTable<FlowInfo> floresInfo_;
1122 VectorBlock velocity;
1124 SparseTable<VelocityInfo> velocityInfo_;
1126 using ScalarFluidState =
typename IntensiveQuantities::ScalarFluidState;
1127 struct BoundaryConditionData
1130 VectorBlock massRate;
1131 unsigned pvtRegionIdx;
1132 unsigned boundaryFaceIndex;
1135 ScalarFluidState exFluidState;
1142 unsigned int bfIndex;
1143 BoundaryConditionData bcdata;
1145 std::vector<BoundaryInfo> boundaryInfo_;
1147 bool separateSparseSourceTerms_ =
false;
1149 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: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