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>
75 static constexpr bool enableBrine = getPropValue<TypeTag, Properties::EnableBrine>();
76 enum { enableVapwat = getPropValue<TypeTag, Properties::EnableVapwat>() };
77 enum { enableDisgasInWater = getPropValue<TypeTag, Properties::EnableDisgasInWater>() };
78 enum { enableSaltPrecipitation = getPropValue<TypeTag, Properties::EnableSaltPrecipitation>() };
79 static constexpr bool enableSolvent = getPropValue<TypeTag, Properties::EnableSolvent>();
80 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
81 static constexpr bool compositionSwitchEnabled = Indices::compositionSwitchIdx >= 0;
89 compositionSwitchEnabled,
92 enableSaltPrecipitation,
100 unsigned globalDofIdx,
103 const EvaluationTemp T = EvaluationTemp::createVariable(problem.temperature(globalDofIdx, timeIdx), 0);
108 const unsigned globalSpaceIdx,
109 const unsigned timeIdx)
113 for (
int phaseIdx = 0; phaseIdx < numPhases; ++ phaseIdx) {
114 if (!FluidSystem::phaseIsActive(phaseIdx)) {
118 const auto& h = FluidSystem::enthalpy(
fluidState_, phaseIdx, problem.pvtRegionIndex(globalSpaceIdx));
122 const auto& solidEnergyLawParams = problem.solidEnergyLawParams(globalSpaceIdx, timeIdx);
125 const auto& thermalConductionLawParams = problem.thermalConductionLawParams(globalSpaceIdx, timeIdx);
133 rockFraction_ = problem.rockFraction(globalSpaceIdx, timeIdx);
148 template <
class Flu
idState>
152 fluidState_.setPvtRegionIndex(fs.pvtRegionIndex());
155 for (
int phaseIdx = 0; phaseIdx < numPhases; ++ phaseIdx) {
156 fluidState_.setPressure(phaseIdx, getValue(fs.pressure(phaseIdx)));
157 fluidState_.setDensity(phaseIdx, getValue(fs.density(phaseIdx)));
158 fluidState_.setSaturation(phaseIdx, getValue(fs.saturation(phaseIdx)));
159 fluidState_.setInvB(phaseIdx, getValue(fs.invB(phaseIdx)));
175template <class TypeTag, bool enableTempV = getPropValue<TypeTag, Properties::EnergyModuleType>() == EnergyModules::SequentialImplicitThermal >
177 GetPropType<TypeTag, Properties::GridView>,
178 GetPropType<TypeTag, Properties::DofMapper>,
179 GetPropType<TypeTag, Properties::Stencil>,
180 GetPropType<TypeTag, Properties::FluidSystem>,
181 GetPropType<TypeTag, Properties::Scalar>>
200 using IndexTraits =
typename FluidSystem::IndexTraitsType;
213 struct ResidualNBInfo
222 enum { numEq = getPropValue<TypeTag, Properties::NumEq>() };
223 enum { numPhases = FluidSystem::numPhases };
224 enum { waterPhaseIdx = FluidSystem::waterPhaseIdx };
225 enum { oilPhaseIdx = FluidSystem::oilPhaseIdx };
226 enum { gasPhaseIdx = FluidSystem::gasPhaseIdx };
227 static constexpr int temperatureIdx = 0;
231 :
BaseType(simulator.vanguard().gridView(),
232 simulator.vanguard().eclState(),
233 simulator.vanguard().cartesianIndexMapper(),
234 simulator.model().dofMapper())
240 const unsigned int numCells =
simulator_.model().numTotalDof();
250 for (
unsigned globI = 0; globI < numCells; ++globI) {
258 const auto& elemMapper =
simulator_.model().elementMapper();
262 scalingFactor_ = getPropValue<TypeTag, Properties::BlackOilEnergyScalingFactor>();
265 using NeighborSet = std::set<unsigned>;
266 std::vector<NeighborSet> neighbors(numCells);
269 std::vector<NeighborInfoCPU> loc_nbinfo;
270 for (
const auto& elem : elements(this->
gridView_)) {
271 stencil.update(elem);
272 for (
unsigned primaryDofIdx = 0; primaryDofIdx < stencil.numPrimaryDof(); ++primaryDofIdx) {
273 const unsigned myIdx = stencil.globalSpaceIndex(primaryDofIdx);
274 loc_nbinfo.resize(stencil.numDof() - 1);
275 for (
unsigned dofIdx = 0; dofIdx < stencil.numDof(); ++dofIdx) {
276 const unsigned neighborIdx = stencil.globalSpaceIndex(dofIdx);
277 neighbors[myIdx].insert(neighborIdx);
279 const auto scvfIdx = dofIdx - 1;
280 const auto& scvf = stencil.interiorFace(scvfIdx);
281 const Scalar area = scvf.area();
282 Scalar inAlpha =
simulator_.problem().thermalHalfTransmissibility(myIdx, neighborIdx);
283 Scalar outAlpha =
simulator_.problem().thermalHalfTransmissibility(neighborIdx, myIdx);
284 ResidualNBInfo nbinfo{area, inAlpha, outAlpha};
285 loc_nbinfo[dofIdx - 1] =
NeighborInfoCPU{neighborIdx, nbinfo,
nullptr};
288 neighborInfo_.appendRow(loc_nbinfo.begin(), loc_nbinfo.end());
295 for (
unsigned globI = 0; globI < numCells; globI++) {
298 for (
auto& nbInfo : nbInfos) {
299 nbInfo.matBlockAddress =
energyMatrix_->blockAddress(nbInfo.neighbor, globI);
311 const unsigned int numCells =
simulator_.model().numTotalDof();
313 #pragma omp parallel for
315 for (
unsigned globI = 0; globI < numCells; ++globI) {
322 const int nw =
simulator_.problem().wellModel().wellState().numWells();
334 OPM_TIMEBLOCK(TemperatureModel_endTimeStep);
337 const unsigned int numCells =
simulator_.model().numTotalDof();
339 #pragma omp parallel for
341 for (
unsigned globI = 0; globI < numCells; ++globI) {
349 const int nw = wellState.
numWells();
350 for (
auto wellID = 0*nw; wellID < nw; ++wellID) {
351 auto& ws = wellState.
well(wellID);
356 const auto& wellPtrs =
simulator_.problem().wellModel().localNonshutWells();
357 for (
const auto& wellPtr : wellPtrs) {
358 auto& ws = wellState.
well(wellPtr->name());
367 template <
class Restarter>
377 template <
class Restarter>
385 const unsigned int numCells =
simulator_.model().numTotalDof();
387 #pragma omp parallel for
389 for (
unsigned globI = 0; globI < numCells; ++globI) {
390 Scalar storage = 0.0;
398 OPM_TIMEBLOCK(TemperatureModel_advanceTemperatureFields);
399 const int max_iter = 20;
400 const int min_iter = 1;
401 bool is_converged =
false;
403 for (
int iter = 0; iter < max_iter; ++iter) {
405 if (iter >= min_iter &&
converged(iter)) {
413 fmt::format(fmt::runtime(
"Temperature model (TEMP): Newton did not converge after {} iterations. \n"
414 "The Simulator will continue to the next step with an unconverged solution."),
422 OPM_TIMEBLOCK(TemperatureModel_solveAndUpdate);
423 const unsigned int numCells =
simulator_.model().numTotalDof();
424 EnergyVector dx(numCells);
427 if (
simulator_.gridView().comm().rank() == 0) {
428 OpmLog::warning(
"Temp model: Linear solver did not converge. Temperature values not updated.");
432 OPM_TIMEBLOCK(TemperatureModel_solveAndUpdate_update);
434 #pragma omp parallel for
436 for (
unsigned globI = 0; globI < numCells; ++globI) {
446 OPM_TIMEBLOCK(TemperatureModel_converged);
448 Scalar maxNorm = 0.0;
449 Scalar sumNorm = 0.0;
450 const auto tolerance_cnv_energy_strict = Parameters::Get<Parameters::ToleranceCnvEnergy<Scalar>>();
451 const auto& elemMapper =
simulator_.model().elementMapper();
454 Scalar sum_pv_not_converged = 0.0;
455 for (
const auto& elem : elements(
simulator_.gridView(), Dune::Partitions::interior)) {
456 unsigned globI = elemMapper.index(elem);
457 const auto pvValue =
simulator_.problem().referencePorosity(globI, 0)
460 const Scalar scaled_norm = dt * std::abs(this->
energyVector_[globI])/ pvValue;
461 maxNorm = max(maxNorm, scaled_norm);
462 sumNorm += scaled_norm;
463 if (!isNumericalAquiferCell(elem)) {
464 if (scaled_norm > tolerance_cnv_energy_strict) {
465 sum_pv_not_converged += pvValue;
471 OPM_TIMEBLOCK(TemperatureModel_converged_communicate);
473 maxNorm =
simulator_.gridView().comm().max(maxNorm);
474 sumNorm =
simulator_.gridView().comm().sum(sumNorm);
475 sum_pv =
simulator_.gridView().comm().sum(sum_pv);
479 sum_pv_not_converged =
simulator_.gridView().comm().sum(sum_pv_not_converged);
481 Scalar relaxed_max_pv_fraction = Parameters::Get<Parameters::RelaxedMaxPvFraction<Scalar>>();
482 const bool relax = (sum_pv_not_converged / sum_pv) < relaxed_max_pv_fraction;
483 const auto tolerance_energy_balance = relax? Parameters::Get<Parameters::ToleranceEnergyBalanceRelaxed<Scalar>>():
485 const bool tolerance_cnv_energy = relax? Parameters::Get<Parameters::ToleranceCnvEnergyRelaxed<Scalar>>():
486 tolerance_cnv_energy_strict;
488 const auto msg = fmt::format(fmt::runtime(
"Temperature model (TEMP): Newton iter {}: "
489 "CNV(E): {:.1e}, EB: {:.1e}"),
490 iter, maxNorm, sumNorm);
492 if (maxNorm < tolerance_cnv_energy && sumNorm < tolerance_energy_balance) {
494 fmt::format(fmt::runtime(
"Temperature model (TEMP): Newton converged after {} iterations"),
502 template<
class LhsEval>
506 const auto& poro = getValue(
simulator_.model().intensiveQuantities(globI, 0).porosity());
508 const auto& fs = intQuants.fluidStateTemp();
509 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++ phaseIdx) {
510 if (!FluidSystem::phaseIsActive(phaseIdx)) {
514 const auto& u = decay<LhsEval>(fs.internalEnergy(phaseIdx));
515 const auto& S = decay<LhsEval>(fs.saturation(phaseIdx));
516 const auto& rho = decay<LhsEval>(fs.density(phaseIdx));
518 storage += poro*S*u*rho;
522 const Scalar rockFraction = intQuants.rockFraction();
523 const auto& uRock = decay<LhsEval>(intQuants.rockInternalEnergy());
524 storage += rockFraction*uRock;
528 template <
class RateVector>
530 const FluidStateTemp& fsEx,
531 const RateVector& darcyFlux,
534 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++ phaseIdx) {
535 if (!FluidSystem::phaseIsActive(phaseIdx)) {
539 const unsigned activeCompIdx =
540 FluidSystem::canonicalToActiveCompIdx(FluidSystem::solventComponentIndex(phaseIdx));
541 bool inIsUp = darcyFlux[activeCompIdx] > 0;
542 const auto& fs = inIsUp ? fsIn : fsEx;
544 flux += fs.enthalpy(phaseIdx)
545 * fs.density(phaseIdx)
546 * getValue(darcyFlux[activeCompIdx]);
549 flux += getValue(fs.enthalpy(phaseIdx))
550 * getValue(fs.density(phaseIdx))
551 * getValue(darcyFlux[activeCompIdx]);
557 template <
class Res
idualNBInfo>
560 const ResidualNBInfo& res_nbinfo,
561 Evaluation& heatFlux)
563 short interiorDofIdx = 0;
564 short exteriorDofIdx = 1;
565 EnergyModule::ExtensiveQuantities::updateEnergy(heatFlux,
575 res_nbinfo.faceArea);
581 OPM_TIMEBLOCK(TemperatureModel_assembleEquations);
583 const unsigned int numCells =
simulator_.model().numTotalDof();
584 for (
unsigned globI = 0; globI < numCells; ++globI) {
592 OPM_TIMEBLOCK(TemperatureModel_assembleEquations_storage);
594 #pragma omp parallel for
596 for (
unsigned globI = 0; globI < numCells; ++globI) {
597 MatrixBlockTemp bMat;
598 Scalar volume =
simulator_.model().dofTotalVolume(globI);
599 Scalar storefac = volume / dt;
600 Evaluation storage = 0.0;
603 bMat[0][0] = storefac * storage.derivative(temperatureIdx);
610 OPM_TIMEBLOCK(TemperatureModel_assembleEquations_flux);
611 const auto& floresInfo = this->
simulator_.problem().model().linearizer().getFloresInfo();
612 const bool enableDriftCompensation = Parameters::Get<Parameters::EnableDriftCompensationTemp>();
616 #pragma omp parallel for
618 for (
unsigned globI = 0; globI < numCells; ++globI) {
620 const auto& floresInfos = floresInfo[globI];
623 MatrixBlockTemp bMat;
624 for (
const auto& nbInfo : nbInfos) {
625 unsigned globJ = nbInfo.neighbor;
627 assert(globJ != globI);
628 const auto& darcyflux = floresInfos[loc].flow;
630 Evaluation flux = 0.0;
631 computeFluxTerm(intQuantsIn.fluidStateTemp(), intQuantsEx.fluidStateTemp(), darcyflux, flux);
633 Evaluation heatFlux = 0.0;
637 bMat[0][0] = heatFlux.derivative(temperatureIdx);
641 *nbInfo.matBlockAddress += bMat;
645 if (enableDriftCompensation) {
646 auto dofDriftRate = problem.drift()[globI]/dt;
647 const auto& fs = intQuantsIn.fluidStateTemp();
648 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++ phaseIdx) {
649 const unsigned activeCompIdx =
650 FluidSystem::canonicalToActiveCompIdx(FluidSystem::solventComponentIndex(phaseIdx));
651 auto drift_hrate = dofDriftRate[activeCompIdx]*getValue(fs.enthalpy(phaseIdx)) * getValue(fs.density(phaseIdx)) / getValue(fs.invB(phaseIdx));
660 OPM_TIMEBLOCK(TemperatureModel_assembleEquations_wells);
661 const auto& wellPtrs =
simulator_.problem().wellModel().localNonshutWells();
662 for (
const auto& wellPtr : wellPtrs) {
685 const auto& eclWell = well.wellEcl();
686 std::size_t well_index =
simulator_.problem().wellModel().wellState().index(well.name()).value();
687 const auto& ws =
simulator_.problem().wellModel().wellState().well(well_index);
689 MatrixBlockTemp bMat;
690 for (std::size_t i = 0; i < ws.perf_data.size(); ++i) {
691 const auto globI = ws.perf_data.cell_index[i];
693 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++ phaseIdx) {
694 if (!FluidSystem::phaseIsActive(phaseIdx)) {
698 Evaluation rate = well.volumetricSurfaceRateForConnection(globI, phaseIdx);
699 if (rate > 0 && eclWell.isInjector()) {
700 fs.setTemperature(eclWell.inj_temperature());
701 const auto& rho = FluidSystem::density(fs, phaseIdx, fs.pvtRegionIndex());
702 fs.setDensity(phaseIdx, rho);
703 const auto& h = FluidSystem::enthalpy(fs, phaseIdx, fs.pvtRegionIndex());
704 fs.setEnthalpy(phaseIdx, h);
705 rate *= getValue(fs.enthalpy(phaseIdx)) * getValue(fs.density(phaseIdx)) / getValue(fs.invB(phaseIdx));
707 const Evaluation d = 1.0 - fs.Rv() * fs.Rs();
708 if (phaseIdx == gasPhaseIdx && d > 0) {
709 const auto& oilrate = well.volumetricSurfaceRateForConnection(globI, oilPhaseIdx);
710 rate -= oilrate * getValue(fs.Rs());
713 if (phaseIdx == oilPhaseIdx && d > 0) {
714 const auto& gasrate = well.volumetricSurfaceRateForConnection(globI, gasPhaseIdx);
715 rate -= gasrate * getValue(fs.Rv());
718 rate *= fs.enthalpy(phaseIdx) * getValue(fs.density(phaseIdx)) / getValue(fs.invB(phaseIdx));
723 bMat[0][0] = -rate.derivative(temperatureIdx);
729 template<
class Well,
class SingleWellState>
732 if (well.isInjector()) {
733 if (ws.
status != WellStatus::STOP) {
738 const int np =
simulator_.problem().wellModel().wellState().numPhases();
739 std::array<Scalar,2> weighted{0.0,0.0};
740 auto& [weighted_temperature, total_weight] = weighted;
741 for (std::size_t i = 0; i < ws.
perf_data.size(); ++i) {
742 const auto globI = ws.
perf_data.cell_index[i];
743 const auto& fs =
intQuants_[globI].fluidStateTemp();
744 Scalar weight_factor =
simulator_.problem().wellModel().computeTemperatureWeightFactor(i, np, fs, ws);
745 total_weight += weight_factor;
746 weighted_temperature += weight_factor * fs.temperature(0).value();
749 ws.
temperature = weighted_temperature / total_weight;
764template <
class TypeTag>
779 template <
class Restarter>
789 template <
class Restarter>
Definition: TemperatureModel.hpp:67
FluidStateTemp fluidState_
Definition: TemperatureModel.hpp:166
OPM_HOST_DEVICE const Scalar & rockFraction() const
Definition: TemperatureModel.hpp:142
OPM_HOST_DEVICE void setFluidState(const FluidState &fs)
Definition: TemperatureModel.hpp:149
Scalar rockFraction_
Definition: TemperatureModel.hpp:167
EvaluationTemp rockInternalEnergy_
Definition: TemperatureModel.hpp:164
EvaluationTemp totalThermalConductivity_
Definition: TemperatureModel.hpp:165
OPM_HOST_DEVICE void updateTemperature_(const Problem &problem, unsigned globalDofIdx, unsigned timeIdx)
Definition: TemperatureModel.hpp:99
OPM_HOST_DEVICE const EvaluationTemp & totalThermalConductivity() const
Definition: TemperatureModel.hpp:139
OPM_HOST_DEVICE void updateEnergyQuantities_(const Problem &problem, const unsigned globalSpaceIdx, const unsigned timeIdx)
Definition: TemperatureModel.hpp:107
OPM_HOST_DEVICE const EvaluationTemp & rockInternalEnergy() const
Definition: TemperatureModel.hpp:136
DenseAd::Evaluation< Scalar, 1 > EvaluationTemp
Definition: TemperatureModel.hpp:84
OPM_HOST_DEVICE const FluidStateTemp & fluidStateTemp() const
Definition: TemperatureModel.hpp:145
BlackOilFluidState< EvaluationTemp, FluidSystem, true, true, compositionSwitchEnabled, enableVapwat, enableBrine, enableSaltPrecipitation, enableDisgasInWater, enableSolvent, Indices::numPhases > FluidStateTemp
Definition: TemperatureModel.hpp:95
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:88
void beginTimeStep()
Definition: TemperatureModel.hpp:794
void init()
Definition: TemperatureModel.hpp:793
const Scalar temperature(size_t) const
Definition: TemperatureModel.hpp:795
void deserialize(Restarter &)
This method restores the complete state of the temperature from disk.
Definition: TemperatureModel.hpp:790
void serialize(Restarter &)
This method writes the complete state of all temperature to the hard disk.
Definition: TemperatureModel.hpp:780
TemperatureModel(Simulator &)
Definition: TemperatureModel.hpp:772
A class which handles sequential implicit solution of the energy equation as specified in by TEMP.
Definition: TemperatureModel.hpp:182
void serialize(Restarter &)
This method writes the complete state of all temperature to the hard disk.
Definition: TemperatureModel.hpp:368
void assembleEquationWell(const Well &well)
Definition: TemperatureModel.hpp:683
void computeFluxTerm(const FluidStateTemp &fsIn, const FluidStateTemp &fsEx, const RateVector &darcyFlux, Evaluation &flux)
Definition: TemperatureModel.hpp:529
std::vector< int > interiorRows_
Definition: TemperatureModel.hpp:759
EnergyVector storage1_
Definition: TemperatureModel.hpp:753
bool converged(const int iter)
Definition: TemperatureModel.hpp:444
void computeStorageTerm(unsigned globI, LhsEval &storage)
Definition: TemperatureModel.hpp:503
void deserialize(Restarter &)
This method restores the complete state of the temperature from disk.
Definition: TemperatureModel.hpp:378
void endTimeStep(WellStateType &wellState)
Informs the temperature model that a time step has just been finished.
Definition: TemperatureModel.hpp:329
const Simulator & simulator_
Definition: TemperatureModel.hpp:752
void computeHeatFluxTerm(const IntensiveQuantitiesTemp &intQuantsIn, const IntensiveQuantitiesTemp &intQuantsEx, const ResidualNBInfo &res_nbinfo, Evaluation &heatFlux)
Definition: TemperatureModel.hpp:558
void beginTimeStep()
Definition: TemperatureModel.hpp:304
void computeWellTemperature(const Well &well, SingleWellState &ws)
Definition: TemperatureModel.hpp:730
void advanceTemperatureFields()
Definition: TemperatureModel.hpp:396
void updateStorageCache()
Definition: TemperatureModel.hpp:382
Scalar scalingFactor_
Definition: TemperatureModel.hpp:760
std::unique_ptr< SpareMatrixEnergyAdapter > energyMatrix_
Definition: TemperatureModel.hpp:757
void assembleEquations()
Definition: TemperatureModel.hpp:579
void solveAndUpdate()
Definition: TemperatureModel.hpp:420
std::vector< MatrixBlockTemp * > diagMatAddress_
Definition: TemperatureModel.hpp:756
std::vector< IntensiveQuantitiesTemp > intQuants_
Definition: TemperatureModel.hpp:754
std::vector< int > overlapRows_
Definition: TemperatureModel.hpp:758
SparseTable< NeighborInfoCPU > neighborInfo_
Definition: TemperatureModel.hpp:755
TemperatureModel(Simulator &simulator)
Definition: TemperatureModel.hpp:230
void init()
Definition: TemperatureModel.hpp:238
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:77
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