28#ifndef OPM_TEMPERATURE_MODEL_HPP
29#define OPM_TEMPERATURE_MODEL_HPP
31#include <opm/common/OpmLog/OpmLog.hpp>
32#include <opm/common/utility/gpuDecorators.hpp>
53template<
class TypeTag,
class MyTypeTag>
62template<
typename Scalar,
typename IndexTraits>
class WellState;
65template <
class TypeTag>
77 enum { enableBrine = getPropValue<TypeTag, Properties::EnableBrine>() };
78 enum { enableVapwat = getPropValue<TypeTag, Properties::EnableVapwat>() };
79 enum { enableDisgasInWater = getPropValue<TypeTag, Properties::EnableDisgasInWater>() };
80 enum { enableSaltPrecipitation = getPropValue<TypeTag, Properties::EnableSaltPrecipitation>() };
81 enum { enableSolvent = getPropValue<TypeTag, Properties::EnableSolvent>() };
82 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
83 static constexpr bool compositionSwitchEnabled = Indices::compositionSwitchIdx >= 0;
91 compositionSwitchEnabled,
94 enableSaltPrecipitation,
102 unsigned globalDofIdx,
105 const EvaluationTemp T = EvaluationTemp::createVariable(problem.temperature(globalDofIdx, timeIdx), 0);
110 const unsigned globalSpaceIdx,
111 const unsigned timeIdx)
115 for (
int phaseIdx = 0; phaseIdx < numPhases; ++ phaseIdx) {
116 if (!FluidSystem::phaseIsActive(phaseIdx)) {
120 const auto& h = FluidSystem::enthalpy(
fluidState_, phaseIdx, problem.pvtRegionIndex(globalSpaceIdx));
124 const auto& solidEnergyLawParams = problem.solidEnergyLawParams(globalSpaceIdx, timeIdx);
127 const auto& thermalConductionLawParams = problem.thermalConductionLawParams(globalSpaceIdx, timeIdx);
135 rockFraction_ = problem.rockFraction(globalSpaceIdx, timeIdx);
150 template <
class Flu
idState>
154 fluidState_.setPvtRegionIndex(fs.pvtRegionIndex());
157 for (
int phaseIdx = 0; phaseIdx < numPhases; ++ phaseIdx) {
158 fluidState_.setPressure(phaseIdx, getValue(fs.pressure(phaseIdx)));
159 fluidState_.setDensity(phaseIdx, getValue(fs.density(phaseIdx)));
160 fluidState_.setSaturation(phaseIdx, getValue(fs.saturation(phaseIdx)));
161 fluidState_.setInvB(phaseIdx, getValue(fs.invB(phaseIdx)));
177template <class TypeTag, bool enableTempV = getPropValue<TypeTag, Properties::EnergyModuleType>() == EnergyModules::SequentialImplicitThermal >
179 GetPropType<TypeTag, Properties::GridView>,
180 GetPropType<TypeTag, Properties::DofMapper>,
181 GetPropType<TypeTag, Properties::Stencil>,
182 GetPropType<TypeTag, Properties::FluidSystem>,
183 GetPropType<TypeTag, Properties::Scalar>>
202 using IndexTraits =
typename FluidSystem::IndexTraitsType;
215 struct ResidualNBInfo
224 enum { numEq = getPropValue<TypeTag, Properties::NumEq>() };
225 enum { numPhases = FluidSystem::numPhases };
226 enum { waterPhaseIdx = FluidSystem::waterPhaseIdx };
227 enum { oilPhaseIdx = FluidSystem::oilPhaseIdx };
228 enum { gasPhaseIdx = FluidSystem::gasPhaseIdx };
229 static constexpr int temperatureIdx = 0;
233 :
BaseType(simulator.vanguard().gridView(),
234 simulator.vanguard().eclState(),
235 simulator.vanguard().cartesianIndexMapper(),
236 simulator.model().dofMapper())
242 const unsigned int numCells =
simulator_.model().numTotalDof();
252 for (
unsigned globI = 0; globI < numCells; ++globI) {
260 const auto& elemMapper =
simulator_.model().elementMapper();
264 scalingFactor_ = getPropValue<TypeTag, Properties::BlackOilEnergyScalingFactor>();
267 using NeighborSet = std::set<unsigned>;
268 std::vector<NeighborSet> neighbors(numCells);
271 std::vector<NeighborInfoCPU> loc_nbinfo;
272 for (
const auto& elem : elements(this->
gridView_)) {
273 stencil.update(elem);
274 for (
unsigned primaryDofIdx = 0; primaryDofIdx < stencil.numPrimaryDof(); ++primaryDofIdx) {
275 const unsigned myIdx = stencil.globalSpaceIndex(primaryDofIdx);
276 loc_nbinfo.resize(stencil.numDof() - 1);
277 for (
unsigned dofIdx = 0; dofIdx < stencil.numDof(); ++dofIdx) {
278 const unsigned neighborIdx = stencil.globalSpaceIndex(dofIdx);
279 neighbors[myIdx].insert(neighborIdx);
281 const auto scvfIdx = dofIdx - 1;
282 const auto& scvf = stencil.interiorFace(scvfIdx);
283 const Scalar area = scvf.area();
284 Scalar inAlpha =
simulator_.problem().thermalHalfTransmissibility(myIdx, neighborIdx);
285 Scalar outAlpha =
simulator_.problem().thermalHalfTransmissibility(neighborIdx, myIdx);
286 ResidualNBInfo nbinfo{area, inAlpha, outAlpha};
287 loc_nbinfo[dofIdx - 1] =
NeighborInfoCPU{neighborIdx, nbinfo,
nullptr};
290 neighborInfo_.appendRow(loc_nbinfo.begin(), loc_nbinfo.end());
297 for (
unsigned globI = 0; globI < numCells; globI++) {
300 for (
auto& nbInfo : nbInfos) {
301 nbInfo.matBlockAddress =
energyMatrix_->blockAddress(nbInfo.neighbor, globI);
313 const unsigned int numCells =
simulator_.model().numTotalDof();
315 #pragma omp parallel for
317 for (
unsigned globI = 0; globI < numCells; ++globI) {
324 const int nw =
simulator_.problem().wellModel().wellState().numWells();
336 OPM_TIMEBLOCK(TemperatureModel_endTimeStep);
339 const unsigned int numCells =
simulator_.model().numTotalDof();
341 #pragma omp parallel for
343 for (
unsigned globI = 0; globI < numCells; ++globI) {
351 const int nw = wellState.
numWells();
352 for (
auto wellID = 0*nw; wellID < nw; ++wellID) {
353 auto& ws = wellState.
well(wellID);
358 const auto& wellPtrs =
simulator_.problem().wellModel().localNonshutWells();
359 for (
const auto& wellPtr : wellPtrs) {
360 auto& ws = wellState.
well(wellPtr->name());
369 template <
class Restarter>
379 template <
class Restarter>
387 const unsigned int numCells =
simulator_.model().numTotalDof();
389 #pragma omp parallel for
391 for (
unsigned globI = 0; globI < numCells; ++globI) {
392 Scalar storage = 0.0;
400 OPM_TIMEBLOCK(TemperatureModel_advanceTemperatureFields);
401 const int max_iter = 20;
402 const int min_iter = 1;
403 bool is_converged =
false;
405 for (
int iter = 0; iter < max_iter; ++iter) {
407 if (iter >= min_iter &&
converged(iter)) {
415 fmt::format(fmt::runtime(
"Temperature model (TEMP): Newton did not converge after {} iterations. \n"
416 "The Simulator will continue to the next step with an unconverged solution."),
424 OPM_TIMEBLOCK(TemperatureModel_solveAndUpdate);
425 const unsigned int numCells =
simulator_.model().numTotalDof();
426 EnergyVector dx(numCells);
429 if (
simulator_.gridView().comm().rank() == 0) {
430 OpmLog::warning(
"Temp model: Linear solver did not converge. Temperature values not updated.");
434 OPM_TIMEBLOCK(TemperatureModel_solveAndUpdate_update);
436 #pragma omp parallel for
438 for (
unsigned globI = 0; globI < numCells; ++globI) {
448 OPM_TIMEBLOCK(TemperatureModel_converged);
450 Scalar maxNorm = 0.0;
451 Scalar sumNorm = 0.0;
452 const auto tolerance_cnv_energy_strict = Parameters::Get<Parameters::ToleranceCnvEnergy<Scalar>>();
453 const auto& elemMapper =
simulator_.model().elementMapper();
456 Scalar sum_pv_not_converged = 0.0;
457 for (
const auto& elem : elements(
simulator_.gridView(), Dune::Partitions::interior)) {
458 unsigned globI = elemMapper.index(elem);
459 const auto pvValue =
simulator_.problem().referencePorosity(globI, 0)
462 const Scalar scaled_norm = dt * std::abs(this->
energyVector_[globI])/ pvValue;
463 maxNorm = max(maxNorm, scaled_norm);
464 sumNorm += scaled_norm;
465 if (!isNumericalAquiferCell(elem)) {
466 if (scaled_norm > tolerance_cnv_energy_strict) {
467 sum_pv_not_converged += pvValue;
473 OPM_TIMEBLOCK(TemperatureModel_converged_communicate);
475 maxNorm =
simulator_.gridView().comm().max(maxNorm);
476 sumNorm =
simulator_.gridView().comm().sum(sumNorm);
477 sum_pv =
simulator_.gridView().comm().sum(sum_pv);
481 sum_pv_not_converged =
simulator_.gridView().comm().sum(sum_pv_not_converged);
483 Scalar relaxed_max_pv_fraction = Parameters::Get<Parameters::RelaxedMaxPvFraction<Scalar>>();
484 const bool relax = (sum_pv_not_converged / sum_pv) < relaxed_max_pv_fraction;
485 const auto tolerance_energy_balance = relax? Parameters::Get<Parameters::ToleranceEnergyBalanceRelaxed<Scalar>>():
487 const bool tolerance_cnv_energy = relax? Parameters::Get<Parameters::ToleranceCnvEnergyRelaxed<Scalar>>():
488 tolerance_cnv_energy_strict;
490 const auto msg = fmt::format(fmt::runtime(
"Temperature model (TEMP): Newton iter {}: "
491 "CNV(E): {:.1e}, EB: {:.1e}"),
492 iter, maxNorm, sumNorm);
494 if (maxNorm < tolerance_cnv_energy && sumNorm < tolerance_energy_balance) {
496 fmt::format(fmt::runtime(
"Temperature model (TEMP): Newton converged after {} iterations"),
504 template<
class LhsEval>
508 const auto& poro = getValue(
simulator_.model().intensiveQuantities(globI, 0).porosity());
510 const auto& fs = intQuants.fluidStateTemp();
511 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++ phaseIdx) {
512 if (!FluidSystem::phaseIsActive(phaseIdx)) {
516 const auto& u = decay<LhsEval>(fs.internalEnergy(phaseIdx));
517 const auto& S = decay<LhsEval>(fs.saturation(phaseIdx));
518 const auto& rho = decay<LhsEval>(fs.density(phaseIdx));
520 storage += poro*S*u*rho;
524 const Scalar rockFraction = intQuants.rockFraction();
525 const auto& uRock = decay<LhsEval>(intQuants.rockInternalEnergy());
526 storage += rockFraction*uRock;
530 template <
class RateVector>
532 const FluidStateTemp& fsEx,
533 const RateVector& darcyFlux,
536 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++ phaseIdx) {
537 if (!FluidSystem::phaseIsActive(phaseIdx)) {
541 const unsigned activeCompIdx =
542 FluidSystem::canonicalToActiveCompIdx(FluidSystem::solventComponentIndex(phaseIdx));
543 bool inIsUp = darcyFlux[activeCompIdx] > 0;
544 const auto& fs = inIsUp ? fsIn : fsEx;
546 flux += fs.enthalpy(phaseIdx)
547 * fs.density(phaseIdx)
548 * getValue(darcyFlux[activeCompIdx]);
551 flux += getValue(fs.enthalpy(phaseIdx))
552 * getValue(fs.density(phaseIdx))
553 * getValue(darcyFlux[activeCompIdx]);
559 template <
class Res
idualNBInfo>
562 const ResidualNBInfo& res_nbinfo,
563 Evaluation& heatFlux)
565 short interiorDofIdx = 0;
566 short exteriorDofIdx = 1;
567 EnergyModule::ExtensiveQuantities::updateEnergy(heatFlux,
577 res_nbinfo.faceArea);
583 OPM_TIMEBLOCK(TemperatureModel_assembleEquations);
585 const unsigned int numCells =
simulator_.model().numTotalDof();
586 for (
unsigned globI = 0; globI < numCells; ++globI) {
594 OPM_TIMEBLOCK(TemperatureModel_assembleEquations_storage);
596 #pragma omp parallel for
598 for (
unsigned globI = 0; globI < numCells; ++globI) {
599 MatrixBlockTemp bMat;
600 Scalar volume =
simulator_.model().dofTotalVolume(globI);
601 Scalar storefac = volume / dt;
602 Evaluation storage = 0.0;
605 bMat[0][0] = storefac * storage.derivative(temperatureIdx);
612 OPM_TIMEBLOCK(TemperatureModel_assembleEquations_flux);
613 const auto& floresInfo = this->
simulator_.problem().model().linearizer().getFloresInfo();
614 const bool enableDriftCompensation = Parameters::Get<Parameters::EnableDriftCompensationTemp>();
618 #pragma omp parallel for
620 for (
unsigned globI = 0; globI < numCells; ++globI) {
622 const auto& floresInfos = floresInfo[globI];
625 MatrixBlockTemp bMat;
626 for (
const auto& nbInfo : nbInfos) {
627 unsigned globJ = nbInfo.neighbor;
629 assert(globJ != globI);
630 const auto& darcyflux = floresInfos[loc].flow;
632 Evaluation flux = 0.0;
633 computeFluxTerm(intQuantsIn.fluidStateTemp(), intQuantsEx.fluidStateTemp(), darcyflux, flux);
635 Evaluation heatFlux = 0.0;
639 bMat[0][0] = heatFlux.derivative(temperatureIdx);
643 *nbInfo.matBlockAddress += bMat;
647 if (enableDriftCompensation) {
648 auto dofDriftRate = problem.drift()[globI]/dt;
649 const auto& fs = intQuantsIn.fluidStateTemp();
650 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++ phaseIdx) {
651 const unsigned activeCompIdx =
652 FluidSystem::canonicalToActiveCompIdx(FluidSystem::solventComponentIndex(phaseIdx));
653 auto drift_hrate = dofDriftRate[activeCompIdx]*getValue(fs.enthalpy(phaseIdx)) * getValue(fs.density(phaseIdx)) / getValue(fs.invB(phaseIdx));
662 OPM_TIMEBLOCK(TemperatureModel_assembleEquations_wells);
663 const auto& wellPtrs =
simulator_.problem().wellModel().localNonshutWells();
664 for (
const auto& wellPtr : wellPtrs) {
687 const auto& eclWell = well.wellEcl();
688 std::size_t well_index =
simulator_.problem().wellModel().wellState().index(well.name()).value();
689 const auto& ws =
simulator_.problem().wellModel().wellState().well(well_index);
691 MatrixBlockTemp bMat;
692 for (std::size_t i = 0; i < ws.perf_data.size(); ++i) {
693 const auto globI = ws.perf_data.cell_index[i];
695 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++ phaseIdx) {
696 if (!FluidSystem::phaseIsActive(phaseIdx)) {
700 Evaluation rate = well.volumetricSurfaceRateForConnection(globI, phaseIdx);
701 if (rate > 0 && eclWell.isInjector()) {
702 fs.setTemperature(eclWell.inj_temperature());
703 const auto& rho = FluidSystem::density(fs, phaseIdx, fs.pvtRegionIndex());
704 fs.setDensity(phaseIdx, rho);
705 const auto& h = FluidSystem::enthalpy(fs, phaseIdx, fs.pvtRegionIndex());
706 fs.setEnthalpy(phaseIdx, h);
707 rate *= getValue(fs.enthalpy(phaseIdx)) * getValue(fs.density(phaseIdx)) / getValue(fs.invB(phaseIdx));
709 const Evaluation d = 1.0 - fs.Rv() * fs.Rs();
710 if (phaseIdx == gasPhaseIdx && d > 0) {
711 const auto& oilrate = well.volumetricSurfaceRateForConnection(globI, oilPhaseIdx);
712 rate -= oilrate * getValue(fs.Rs());
715 if (phaseIdx == oilPhaseIdx && d > 0) {
716 const auto& gasrate = well.volumetricSurfaceRateForConnection(globI, gasPhaseIdx);
717 rate -= gasrate * getValue(fs.Rv());
720 rate *= fs.enthalpy(phaseIdx) * getValue(fs.density(phaseIdx)) / getValue(fs.invB(phaseIdx));
725 bMat[0][0] = -rate.derivative(temperatureIdx);
731 template<
class Well,
class SingleWellState>
734 if (well.isInjector()) {
735 if (ws.
status != WellStatus::STOP) {
740 const int np =
simulator_.problem().wellModel().wellState().numPhases();
741 std::array<Scalar,2> weighted{0.0,0.0};
742 auto& [weighted_temperature, total_weight] = weighted;
743 for (std::size_t i = 0; i < ws.
perf_data.size(); ++i) {
744 const auto globI = ws.
perf_data.cell_index[i];
745 const auto& fs =
intQuants_[globI].fluidStateTemp();
746 Scalar weight_factor =
simulator_.problem().wellModel().computeTemperatureWeightFactor(i, np, fs, ws);
747 total_weight += weight_factor;
748 weighted_temperature += weight_factor * fs.temperature(0).value();
751 ws.
temperature = weighted_temperature / total_weight;
766template <
class TypeTag>
781 template <
class Restarter>
791 template <
class Restarter>
Definition: TemperatureModel.hpp:67
FluidStateTemp fluidState_
Definition: TemperatureModel.hpp:168
OPM_HOST_DEVICE const Scalar & rockFraction() const
Definition: TemperatureModel.hpp:144
OPM_HOST_DEVICE void setFluidState(const FluidState &fs)
Definition: TemperatureModel.hpp:151
Scalar rockFraction_
Definition: TemperatureModel.hpp:169
EvaluationTemp rockInternalEnergy_
Definition: TemperatureModel.hpp:166
EvaluationTemp totalThermalConductivity_
Definition: TemperatureModel.hpp:167
OPM_HOST_DEVICE void updateTemperature_(const Problem &problem, unsigned globalDofIdx, unsigned timeIdx)
Definition: TemperatureModel.hpp:101
OPM_HOST_DEVICE const EvaluationTemp & totalThermalConductivity() const
Definition: TemperatureModel.hpp:141
OPM_HOST_DEVICE void updateEnergyQuantities_(const Problem &problem, const unsigned globalSpaceIdx, const unsigned timeIdx)
Definition: TemperatureModel.hpp:109
OPM_HOST_DEVICE const EvaluationTemp & rockInternalEnergy() const
Definition: TemperatureModel.hpp:138
DenseAd::Evaluation< Scalar, 1 > EvaluationTemp
Definition: TemperatureModel.hpp:86
OPM_HOST_DEVICE const FluidStateTemp & fluidStateTemp() const
Definition: TemperatureModel.hpp:147
BlackOilFluidState< EvaluationTemp, FluidSystem, true, true, compositionSwitchEnabled, enableVapwat, enableBrine, enableSaltPrecipitation, enableDisgasInWater, enableSolvent, Indices::numPhases > FluidStateTemp
Definition: TemperatureModel.hpp:97
Contains the high level supplements required to extend the black oil model by energy.
Definition: blackoilenergymodules.hh:64
Definition: GenericTemperatureModel.hpp:53
Opm::GenericTemperatureModel< GetPropType< TypeTag, Properties::Grid >, GetPropType< TypeTag, Properties::GridView >, GetPropType< TypeTag, Properties::DofMapper >, GetPropType< TypeTag, Properties::Stencil >, GetPropType< TypeTag, Properties::FluidSystem >, GetPropType< TypeTag, Properties::Scalar > >::energyVector_ EnergyVector energyVector_
Definition: GenericTemperatureModel.hpp:91
Opm::GenericTemperatureModel< GetPropType< TypeTag, Properties::Grid >, GetPropType< TypeTag, Properties::GridView >, GetPropType< TypeTag, Properties::DofMapper >, GetPropType< TypeTag, Properties::Stencil >, GetPropType< TypeTag, Properties::FluidSystem >, GetPropType< TypeTag, Properties::Scalar > >::dofMapper_ const GetPropType< TypeTag, Properties::DofMapper > & dofMapper_
Definition: GenericTemperatureModel.hpp:89
Opm::GenericTemperatureModel< GetPropType< TypeTag, Properties::Grid >, GetPropType< TypeTag, Properties::GridView >, GetPropType< TypeTag, Properties::DofMapper >, GetPropType< TypeTag, Properties::Stencil >, GetPropType< TypeTag, Properties::FluidSystem >, GetPropType< TypeTag, Properties::Scalar > >::energy_rates_ std::vector< GetPropType< TypeTag, Properties::Scalar > > energy_rates_
Definition: GenericTemperatureModel.hpp:93
Opm::GenericTemperatureModel< GetPropType< TypeTag, Properties::Grid >, GetPropType< TypeTag, Properties::GridView >, GetPropType< TypeTag, Properties::DofMapper >, GetPropType< TypeTag, Properties::Stencil >, GetPropType< TypeTag, Properties::FluidSystem >, GetPropType< TypeTag, Properties::Scalar > >::MatrixBlockTemp MatrixBlock< GetPropType< TypeTag, Properties::Scalar >, 1, 1 > MatrixBlockTemp
Definition: GenericTemperatureModel.hpp:56
Opm::GenericTemperatureModel< GetPropType< TypeTag, Properties::Grid >, GetPropType< TypeTag, Properties::GridView >, GetPropType< TypeTag, Properties::DofMapper >, GetPropType< TypeTag, Properties::Stencil >, GetPropType< TypeTag, Properties::FluidSystem >, GetPropType< TypeTag, Properties::Scalar > >::EnergyMatrix Dune::BCRSMatrix< MatrixBlockTemp > EnergyMatrix
Definition: GenericTemperatureModel.hpp:57
Opm::GenericTemperatureModel< GetPropType< TypeTag, Properties::Grid >, GetPropType< TypeTag, Properties::GridView >, GetPropType< TypeTag, Properties::DofMapper >, GetPropType< TypeTag, Properties::Stencil >, GetPropType< TypeTag, Properties::FluidSystem >, GetPropType< TypeTag, Properties::Scalar > >::linearSolve_ bool linearSolve_(const EnergyMatrix &M, EnergyVector &x, EnergyVector &b)
Definition: GenericTemperatureModel_impl.hpp:166
Opm::GenericTemperatureModel< GetPropType< TypeTag, Properties::Grid >, GetPropType< TypeTag, Properties::GridView >, GetPropType< TypeTag, Properties::DofMapper >, GetPropType< TypeTag, Properties::Stencil >, GetPropType< TypeTag, Properties::FluidSystem >, GetPropType< TypeTag, Properties::Scalar > >::EnergyVector Dune::BlockVector< Dune::FieldVector< GetPropType< TypeTag, Properties::Scalar >, 1 > > EnergyVector
Definition: GenericTemperatureModel.hpp:58
Opm::GenericTemperatureModel< GetPropType< TypeTag, Properties::Grid >, GetPropType< TypeTag, Properties::GridView >, GetPropType< TypeTag, Properties::DofMapper >, GetPropType< TypeTag, Properties::Stencil >, GetPropType< TypeTag, Properties::FluidSystem >, GetPropType< TypeTag, Properties::Scalar > >::doTemp bool doTemp()
Definition: GenericTemperatureModel.hpp:62
Opm::GenericTemperatureModel< GetPropType< TypeTag, Properties::Grid >, GetPropType< TypeTag, Properties::GridView >, GetPropType< TypeTag, Properties::DofMapper >, GetPropType< TypeTag, Properties::Stencil >, GetPropType< TypeTag, Properties::FluidSystem >, GetPropType< TypeTag, Properties::Scalar > >::doInit void doInit(std::size_t numGridDof)
Initialize all internal data structures needed by the temperature module.
Definition: GenericTemperatureModel_impl.hpp:115
Opm::GenericTemperatureModel< GetPropType< TypeTag, Properties::Grid >, GetPropType< TypeTag, Properties::GridView >, GetPropType< TypeTag, Properties::DofMapper >, GetPropType< TypeTag, Properties::Stencil >, GetPropType< TypeTag, Properties::FluidSystem >, GetPropType< TypeTag, Properties::Scalar > >::maxTempChange_ GetPropType< TypeTag, Properties::Scalar > maxTempChange_
Definition: GenericTemperatureModel.hpp:95
Opm::GenericTemperatureModel< GetPropType< TypeTag, Properties::Grid >, GetPropType< TypeTag, Properties::GridView >, GetPropType< TypeTag, Properties::DofMapper >, GetPropType< TypeTag, Properties::Stencil >, GetPropType< TypeTag, Properties::FluidSystem >, GetPropType< TypeTag, Properties::Scalar > >::temperature_ std::vector< GetPropType< TypeTag, Properties::Scalar > > temperature_
Definition: GenericTemperatureModel.hpp:92
Opm::GenericTemperatureModel< GetPropType< TypeTag, Properties::Grid >, GetPropType< TypeTag, Properties::GridView >, GetPropType< TypeTag, Properties::DofMapper >, GetPropType< TypeTag, Properties::Stencil >, GetPropType< TypeTag, Properties::FluidSystem >, GetPropType< TypeTag, Properties::Scalar > >::gridView_ const GetPropType< TypeTag, Properties::GridView > & gridView_
Definition: GenericTemperatureModel.hpp:86
A sparse matrix interface backend for BCRSMatrix from dune-istl.
Definition: istlsparsematrixadapter.hh:43
Definition: SingleWellState.hpp:44
WellStatus status
Definition: SingleWellState.hpp:100
Scalar temperature
Definition: SingleWellState.hpp:108
PerfData< Scalar > perf_data
Definition: SingleWellState.hpp:156
Definition: BlackoilWellModel.hpp:87
void beginTimeStep()
Definition: TemperatureModel.hpp:796
void init()
Definition: TemperatureModel.hpp:795
const Scalar temperature(size_t) const
Definition: TemperatureModel.hpp:797
void deserialize(Restarter &)
This method restores the complete state of the temperature from disk.
Definition: TemperatureModel.hpp:792
void serialize(Restarter &)
This method writes the complete state of all temperature to the hard disk.
Definition: TemperatureModel.hpp:782
TemperatureModel(Simulator &)
Definition: TemperatureModel.hpp:774
A class which handles sequential implicit solution of the energy equation as specified in by TEMP.
Definition: TemperatureModel.hpp:184
void serialize(Restarter &)
This method writes the complete state of all temperature to the hard disk.
Definition: TemperatureModel.hpp:370
void assembleEquationWell(const Well &well)
Definition: TemperatureModel.hpp:685
void computeFluxTerm(const FluidStateTemp &fsIn, const FluidStateTemp &fsEx, const RateVector &darcyFlux, Evaluation &flux)
Definition: TemperatureModel.hpp:531
std::vector< int > interiorRows_
Definition: TemperatureModel.hpp:761
EnergyVector storage1_
Definition: TemperatureModel.hpp:755
bool converged(const int iter)
Definition: TemperatureModel.hpp:446
void computeStorageTerm(unsigned globI, LhsEval &storage)
Definition: TemperatureModel.hpp:505
void deserialize(Restarter &)
This method restores the complete state of the temperature from disk.
Definition: TemperatureModel.hpp:380
void endTimeStep(WellStateType &wellState)
Informs the temperature model that a time step has just been finished.
Definition: TemperatureModel.hpp:331
const Simulator & simulator_
Definition: TemperatureModel.hpp:754
void computeHeatFluxTerm(const IntensiveQuantitiesTemp &intQuantsIn, const IntensiveQuantitiesTemp &intQuantsEx, const ResidualNBInfo &res_nbinfo, Evaluation &heatFlux)
Definition: TemperatureModel.hpp:560
void beginTimeStep()
Definition: TemperatureModel.hpp:306
void computeWellTemperature(const Well &well, SingleWellState &ws)
Definition: TemperatureModel.hpp:732
void advanceTemperatureFields()
Definition: TemperatureModel.hpp:398
void updateStorageCache()
Definition: TemperatureModel.hpp:384
Scalar scalingFactor_
Definition: TemperatureModel.hpp:762
std::unique_ptr< SpareMatrixEnergyAdapter > energyMatrix_
Definition: TemperatureModel.hpp:759
void assembleEquations()
Definition: TemperatureModel.hpp:581
void solveAndUpdate()
Definition: TemperatureModel.hpp:422
std::vector< MatrixBlockTemp * > diagMatAddress_
Definition: TemperatureModel.hpp:758
std::vector< IntensiveQuantitiesTemp > intQuants_
Definition: TemperatureModel.hpp:756
std::vector< int > overlapRows_
Definition: TemperatureModel.hpp:760
SparseTable< NeighborInfoCPU > neighborInfo_
Definition: TemperatureModel.hpp:757
TemperatureModel(Simulator &simulator)
Definition: TemperatureModel.hpp:232
void init()
Definition: TemperatureModel.hpp:240
Definition: WellState.hpp:66
int numWells() const
Definition: WellState.hpp:114
const SingleWellState< Scalar, IndexTraits > & well(std::size_t well_index) const
Definition: WellState.hpp:306
Provides data handles for parallel communication which operate on DOFs.
auto Get(bool errorIfNotRegistered=true)
Retrieve a runtime parameter.
Definition: parametersystem.hpp:187
Definition: blackoilmodel.hh:80
void findOverlapAndInterior(const Grid &grid, const Mapper &mapper, std::vector< int > &overlapRows, std::vector< int > &interiorRows)
Find the rows corresponding to overlap cells.
Definition: findOverlapRowsAndColumns.hpp:92
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
The Opm property system, traits with inheritance.
Definition: AquiferGridUtils.hpp:35
Definition: tpfalinearizer.hh:124
Definition: BlackoilModelParameters.hpp:61
Definition: TemperatureModel.hpp:54
a tag to mark properties as undefined
Definition: propertysystem.hh:38