28#ifndef OPM_TEMPERATURE_MODEL_HPP
29#define OPM_TEMPERATURE_MODEL_HPP
31#include <opm/common/OpmLog/OpmLog.hpp>
52template<
class TypeTag,
class MyTypeTag>
61template<
typename Scalar,
typename IndexTraits>
class WellState;
64template <
class TypeTag>
76 enum { enableBrine = getPropValue<TypeTag, Properties::EnableBrine>() };
77 enum { enableVapwat = getPropValue<TypeTag, Properties::EnableVapwat>() };
78 enum { enableDisgasInWater = getPropValue<TypeTag, Properties::EnableDisgasInWater>() };
79 enum { enableSaltPrecipitation = getPropValue<TypeTag, Properties::EnableSaltPrecipitation>() };
80 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
81 static constexpr bool compositionSwitchEnabled = Indices::compositionSwitchIdx >= 0;
89 compositionSwitchEnabled,
92 enableSaltPrecipitation,
100 const EvaluationTemp T = EvaluationTemp::createVariable(problem.temperature(globalDofIdx, timeIdx), 0);
105 const unsigned globalSpaceIdx,
106 const unsigned timeIdx)
110 for (
int phaseIdx = 0; phaseIdx < numPhases; ++ phaseIdx) {
111 if (!FluidSystem::phaseIsActive(phaseIdx)) {
115 const auto& h = FluidSystem::enthalpy(
fluidState_, phaseIdx, problem.pvtRegionIndex(globalSpaceIdx));
119 const auto& solidEnergyLawParams = problem.solidEnergyLawParams(globalSpaceIdx, timeIdx);
122 const auto& thermalConductionLawParams = problem.thermalConductionLawParams(globalSpaceIdx, timeIdx);
130 rockFraction_ = problem.rockFraction(globalSpaceIdx, timeIdx);
145 template <
class Flu
idState>
149 fluidState_.setPvtRegionIndex(fs.pvtRegionIndex());
152 for (
int phaseIdx = 0; phaseIdx < numPhases; ++ phaseIdx) {
153 fluidState_.setPressure(phaseIdx, getValue(fs.pressure(phaseIdx)));
154 fluidState_.setDensity(phaseIdx, getValue(fs.density(phaseIdx)));
155 fluidState_.setSaturation(phaseIdx, getValue(fs.saturation(phaseIdx)));
156 fluidState_.setInvB(phaseIdx, getValue(fs.invB(phaseIdx)));
172template <class TypeTag, bool enableTempV = getPropValue<TypeTag, Properties::EnergyModuleType>() == EnergyModules::SequentialImplicitThermal >
174 GetPropType<TypeTag, Properties::GridView>,
175 GetPropType<TypeTag, Properties::DofMapper>,
176 GetPropType<TypeTag, Properties::Stencil>,
177 GetPropType<TypeTag, Properties::FluidSystem>,
178 GetPropType<TypeTag, Properties::Scalar>>
197 using IndexTraits =
typename FluidSystem::IndexTraitsType;
210 struct ResidualNBInfo
219 enum { numEq = getPropValue<TypeTag, Properties::NumEq>() };
220 enum { numPhases = FluidSystem::numPhases };
221 enum { waterPhaseIdx = FluidSystem::waterPhaseIdx };
222 enum { oilPhaseIdx = FluidSystem::oilPhaseIdx };
223 enum { gasPhaseIdx = FluidSystem::gasPhaseIdx };
224 static constexpr int temperatureIdx = 0;
228 :
BaseType(simulator.vanguard().gridView(),
229 simulator.vanguard().eclState(),
230 simulator.vanguard().cartesianIndexMapper(),
231 simulator.model().dofMapper())
237 const unsigned int numCells =
simulator_.model().numTotalDof();
247 for (
unsigned globI = 0; globI < numCells; ++globI) {
255 const auto& elemMapper =
simulator_.model().elementMapper();
259 scalingFactor_ = getPropValue<TypeTag, Properties::BlackOilEnergyScalingFactor>();
262 using NeighborSet = std::set<unsigned>;
263 std::vector<NeighborSet> neighbors(numCells);
266 std::vector<NeighborInfoCPU> loc_nbinfo;
267 for (
const auto& elem : elements(this->
gridView_)) {
268 stencil.update(elem);
269 for (
unsigned primaryDofIdx = 0; primaryDofIdx < stencil.numPrimaryDof(); ++primaryDofIdx) {
270 const unsigned myIdx = stencil.globalSpaceIndex(primaryDofIdx);
271 loc_nbinfo.resize(stencil.numDof() - 1);
272 for (
unsigned dofIdx = 0; dofIdx < stencil.numDof(); ++dofIdx) {
273 const unsigned neighborIdx = stencil.globalSpaceIndex(dofIdx);
274 neighbors[myIdx].insert(neighborIdx);
276 const auto scvfIdx = dofIdx - 1;
277 const auto& scvf = stencil.interiorFace(scvfIdx);
278 const Scalar area = scvf.area();
279 Scalar inAlpha =
simulator_.problem().thermalHalfTransmissibility(myIdx, neighborIdx);
280 Scalar outAlpha =
simulator_.problem().thermalHalfTransmissibility(neighborIdx, myIdx);
281 ResidualNBInfo nbinfo{area, inAlpha, outAlpha};
282 loc_nbinfo[dofIdx - 1] =
NeighborInfoCPU{neighborIdx, nbinfo,
nullptr};
285 neighborInfo_.appendRow(loc_nbinfo.begin(), loc_nbinfo.end());
292 for (
unsigned globI = 0; globI < numCells; globI++) {
295 for (
auto& nbInfo : nbInfos) {
296 nbInfo.matBlockAddress =
energyMatrix_->blockAddress(nbInfo.neighbor, globI);
308 const unsigned int numCells =
simulator_.model().numTotalDof();
310 #pragma omp parallel for
312 for (
unsigned globI = 0; globI < numCells; ++globI) {
319 const int nw =
simulator_.problem().wellModel().wellState().numWells();
333 const unsigned int numCells =
simulator_.model().numTotalDof();
335 #pragma omp parallel for
337 for (
unsigned globI = 0; globI < numCells; ++globI) {
345 const int nw = wellState.
numWells();
346 for (
auto wellID = 0*nw; wellID < nw; ++wellID) {
347 auto& ws = wellState.
well(wellID);
352 const auto& wellPtrs =
simulator_.problem().wellModel().localNonshutWells();
353 for (
const auto& wellPtr : wellPtrs) {
354 auto& ws = wellState.
well(wellPtr->name());
363 template <
class Restarter>
373 template <
class Restarter>
381 const unsigned int numCells =
simulator_.model().numTotalDof();
383 #pragma omp parallel for
385 for (
unsigned globI = 0; globI < numCells; ++globI) {
386 Scalar storage = 0.0;
394 const int max_iter = 20;
395 const int min_iter = 1;
396 bool is_converged =
false;
398 for (
int iter = 0; iter < max_iter; ++iter) {
400 if (iter >= min_iter &&
converged(iter)) {
408 fmt::format(fmt::runtime(
"Temperature model (TEMP): Newton did not converge after {} iterations. \n"
409 "The Simulator will continue to the next step with an unconverged solution."),
417 const unsigned int numCells =
simulator_.model().numTotalDof();
418 EnergyVector dx(numCells);
421 if (
simulator_.gridView().comm().rank() == 0) {
422 OpmLog::warning(
"Temp model: Linear solver did not converge. Temperature values not updated.");
427 #pragma omp parallel for
429 for (
unsigned globI = 0; globI < numCells; ++globI) {
440 Scalar maxNorm = 0.0;
441 Scalar sumNorm = 0.0;
442 const auto tolerance_cnv_energy_strict = Parameters::Get<Parameters::ToleranceCnvEnergy<Scalar>>();
443 const auto& elemMapper =
simulator_.model().elementMapper();
446 Scalar sum_pv_not_converged = 0.0;
447 for (
const auto& elem : elements(
simulator_.gridView(), Dune::Partitions::interior)) {
448 unsigned globI = elemMapper.index(elem);
449 const auto pvValue =
simulator_.problem().referencePorosity(globI, 0)
452 const Scalar scaled_norm = dt * std::abs(this->
energyVector_[globI])/ pvValue;
453 maxNorm = max(maxNorm, scaled_norm);
454 sumNorm += scaled_norm;
455 if (!isNumericalAquiferCell(elem)) {
456 if (scaled_norm > tolerance_cnv_energy_strict) {
457 sum_pv_not_converged += pvValue;
462 maxNorm =
simulator_.gridView().comm().max(maxNorm);
463 sumNorm =
simulator_.gridView().comm().sum(sumNorm);
464 sum_pv =
simulator_.gridView().comm().sum(sum_pv);
468 sum_pv_not_converged =
simulator_.gridView().comm().sum(sum_pv_not_converged);
469 Scalar relaxed_max_pv_fraction = Parameters::Get<Parameters::RelaxedMaxPvFraction<Scalar>>();
470 const bool relax = (sum_pv_not_converged / sum_pv) < relaxed_max_pv_fraction;
471 const auto tolerance_energy_balance = relax? Parameters::Get<Parameters::ToleranceEnergyBalanceRelaxed<Scalar>>():
473 const bool tolerance_cnv_energy = relax? Parameters::Get<Parameters::ToleranceCnvEnergyRelaxed<Scalar>>():
474 tolerance_cnv_energy_strict;
476 const auto msg = fmt::format(fmt::runtime(
"Temperature model (TEMP): Newton iter {}: "
477 "CNV(E): {:.1e}, EB: {:.1e}"),
478 iter, maxNorm, sumNorm);
480 if (maxNorm < tolerance_cnv_energy && sumNorm < tolerance_energy_balance) {
482 fmt::format(fmt::runtime(
"Temperature model (TEMP): Newton converged after {} iterations"),
490 template<
class LhsEval>
494 const auto& poro = getValue(
simulator_.model().intensiveQuantities(globI, 0).porosity());
496 const auto& fs = intQuants.fluidStateTemp();
497 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++ phaseIdx) {
498 if (!FluidSystem::phaseIsActive(phaseIdx)) {
502 const auto& u = decay<LhsEval>(fs.internalEnergy(phaseIdx));
503 const auto& S = decay<LhsEval>(fs.saturation(phaseIdx));
504 const auto& rho = decay<LhsEval>(fs.density(phaseIdx));
506 storage += poro*S*u*rho;
510 const Scalar rockFraction = intQuants.rockFraction();
511 const auto& uRock = decay<LhsEval>(intQuants.rockInternalEnergy());
512 storage += rockFraction*uRock;
516 template <
class RateVector>
518 const FluidStateTemp& fsEx,
519 const RateVector& darcyFlux,
522 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++ phaseIdx) {
523 if (!FluidSystem::phaseIsActive(phaseIdx)) {
527 const unsigned activeCompIdx =
528 FluidSystem::canonicalToActiveCompIdx(FluidSystem::solventComponentIndex(phaseIdx));
529 bool inIsUp = darcyFlux[activeCompIdx] > 0;
530 const auto& fs = inIsUp ? fsIn : fsEx;
532 flux += fs.enthalpy(phaseIdx)
533 * fs.density(phaseIdx)
534 * getValue(darcyFlux[activeCompIdx]);
537 flux += getValue(fs.enthalpy(phaseIdx))
538 * getValue(fs.density(phaseIdx))
539 * getValue(darcyFlux[activeCompIdx]);
545 template <
class Res
idualNBInfo>
548 const ResidualNBInfo& res_nbinfo,
549 Evaluation& heatFlux)
551 short interiorDofIdx = 0;
552 short exteriorDofIdx = 1;
553 EnergyModule::ExtensiveQuantities::updateEnergy(heatFlux,
563 res_nbinfo.faceArea);
569 const unsigned int numCells =
simulator_.model().numTotalDof();
570 for (
unsigned globI = 0; globI < numCells; ++globI) {
576 #pragma omp parallel for
578 for (
unsigned globI = 0; globI < numCells; ++globI) {
579 MatrixBlockTemp bMat;
580 Scalar volume =
simulator_.model().dofTotalVolume(globI);
581 Scalar storefac = volume / dt;
582 Evaluation storage = 0.0;
585 bMat[0][0] = storefac * storage.derivative(temperatureIdx);
589 const auto& floresInfo = this->
simulator_.problem().model().linearizer().getFloresInfo();
590 const bool enableDriftCompensation = Parameters::Get<Parameters::EnableDriftCompensationTemp>();
594 #pragma omp parallel for
596 for (
unsigned globI = 0; globI < numCells; ++globI) {
598 const auto& floresInfos = floresInfo[globI];
601 MatrixBlockTemp bMat;
602 for (
const auto& nbInfo : nbInfos) {
603 unsigned globJ = nbInfo.neighbor;
605 assert(globJ != globI);
606 const auto& darcyflux = floresInfos[loc].flow;
608 Evaluation flux = 0.0;
609 computeFluxTerm(intQuantsIn.fluidStateTemp(), intQuantsEx.fluidStateTemp(), darcyflux, flux);
611 Evaluation heatFlux = 0.0;
615 bMat[0][0] = heatFlux.derivative(temperatureIdx);
619 *nbInfo.matBlockAddress += bMat;
623 if (enableDriftCompensation) {
624 auto dofDriftRate = problem.drift()[globI]/dt;
625 const auto& fs = intQuantsIn.fluidStateTemp();
626 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++ phaseIdx) {
627 const unsigned activeCompIdx =
628 FluidSystem::canonicalToActiveCompIdx(FluidSystem::solventComponentIndex(phaseIdx));
629 auto drift_hrate = dofDriftRate[activeCompIdx]*getValue(fs.enthalpy(phaseIdx)) * getValue(fs.density(phaseIdx)) / getValue(fs.invB(phaseIdx));
636 const auto& wellPtrs =
simulator_.problem().wellModel().localNonshutWells();
637 for (
const auto& wellPtr : wellPtrs) {
659 const auto& eclWell = well.wellEcl();
660 std::size_t well_index =
simulator_.problem().wellModel().wellState().index(well.name()).value();
661 const auto& ws =
simulator_.problem().wellModel().wellState().well(well_index);
663 MatrixBlockTemp bMat;
664 for (std::size_t i = 0; i < ws.perf_data.size(); ++i) {
665 const auto globI = ws.perf_data.cell_index[i];
667 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++ phaseIdx) {
668 if (!FluidSystem::phaseIsActive(phaseIdx)) {
672 Evaluation rate = well.volumetricSurfaceRateForConnection(globI, phaseIdx);
673 if (rate > 0 && eclWell.isInjector()) {
674 fs.setTemperature(eclWell.inj_temperature());
675 const auto& rho = FluidSystem::density(fs, phaseIdx, fs.pvtRegionIndex());
676 fs.setDensity(phaseIdx, rho);
677 const auto& h = FluidSystem::enthalpy(fs, phaseIdx, fs.pvtRegionIndex());
678 fs.setEnthalpy(phaseIdx, h);
679 rate *= getValue(fs.enthalpy(phaseIdx)) * getValue(fs.density(phaseIdx)) / getValue(fs.invB(phaseIdx));
681 const Evaluation d = 1.0 - fs.Rv() * fs.Rs();
682 if (phaseIdx == gasPhaseIdx && d > 0) {
683 const auto& oilrate = well.volumetricSurfaceRateForConnection(globI, oilPhaseIdx);
684 rate -= oilrate * getValue(fs.Rs());
687 if (phaseIdx == oilPhaseIdx && d > 0) {
688 const auto& gasrate = well.volumetricSurfaceRateForConnection(globI, gasPhaseIdx);
689 rate -= gasrate * getValue(fs.Rv());
692 rate *= fs.enthalpy(phaseIdx) * getValue(fs.density(phaseIdx)) / getValue(fs.invB(phaseIdx));
697 bMat[0][0] = -rate.derivative(temperatureIdx);
703 template<
class Well,
class SingleWellState>
706 if (well.isInjector()) {
707 if (ws.
status != WellStatus::STOP) {
712 const int np =
simulator_.problem().wellModel().wellState().numPhases();
713 std::array<Scalar,2> weighted{0.0,0.0};
714 auto& [weighted_temperature, total_weight] = weighted;
715 for (std::size_t i = 0; i < ws.
perf_data.size(); ++i) {
716 const auto globI = ws.
perf_data.cell_index[i];
717 const auto& fs =
intQuants_[globI].fluidStateTemp();
718 Scalar weight_factor =
simulator_.problem().wellModel().computeTemperatureWeightFactor(i, np, fs, ws);
719 total_weight += weight_factor;
720 weighted_temperature += weight_factor * fs.temperature(0).value();
723 ws.
temperature = weighted_temperature / total_weight;
738template <
class TypeTag>
753 template <
class Restarter>
763 template <
class Restarter>
Definition: TemperatureModel.hpp:66
FluidStateTemp fluidState_
Definition: TemperatureModel.hpp:163
BlackOilFluidState< EvaluationTemp, FluidSystem, true, true, compositionSwitchEnabled, enableVapwat, enableBrine, enableSaltPrecipitation, enableDisgasInWater, Indices::numPhases > FluidStateTemp
Definition: TemperatureModel.hpp:94
void updateTemperature_(const Problem &problem, unsigned globalDofIdx, unsigned timeIdx)
Definition: TemperatureModel.hpp:98
void updateEnergyQuantities_(const Problem &problem, const unsigned globalSpaceIdx, const unsigned timeIdx)
Definition: TemperatureModel.hpp:104
void setFluidState(const FluidState &fs)
Definition: TemperatureModel.hpp:146
const FluidStateTemp & fluidStateTemp() const
Definition: TemperatureModel.hpp:142
Scalar rockFraction_
Definition: TemperatureModel.hpp:164
EvaluationTemp rockInternalEnergy_
Definition: TemperatureModel.hpp:161
EvaluationTemp totalThermalConductivity_
Definition: TemperatureModel.hpp:162
const Scalar & rockFraction() const
Definition: TemperatureModel.hpp:139
const EvaluationTemp & totalThermalConductivity() const
Definition: TemperatureModel.hpp:136
const EvaluationTemp & rockInternalEnergy() const
Definition: TemperatureModel.hpp:133
DenseAd::Evaluation< Scalar, 1 > EvaluationTemp
Definition: TemperatureModel.hpp:84
Contains the high level supplements required to extend the black oil model by energy.
Definition: blackoilenergymodules.hh:61
Definition: GenericTemperatureModel.hpp:52
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: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 > >::dofMapper_ const GetPropType< TypeTag, Properties::DofMapper > & dofMapper_
Definition: GenericTemperatureModel.hpp:87
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: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 > >::MatrixBlockTemp MatrixBlock< GetPropType< TypeTag, Properties::Scalar >, 1, 1 > MatrixBlockTemp
Definition: GenericTemperatureModel.hpp:55
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: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 > >::linearSolve_ bool linearSolve_(const EnergyMatrix &M, EnergyVector &x, EnergyVector &b)
Definition: GenericTemperatureModel_impl.hpp:129
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: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 > >::doTemp bool doTemp()
Definition: GenericTemperatureModel.hpp:61
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: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 > >::temperature_ std::vector< GetPropType< TypeTag, Properties::Scalar > > temperature_
Definition: GenericTemperatureModel.hpp:90
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:84
A sparse matrix interface backend for BCRSMatrix from dune-istl.
Definition: istlsparsematrixadapter.hh:43
Definition: SingleWellState.hpp:43
WellStatus status
Definition: SingleWellState.hpp:97
Scalar temperature
Definition: SingleWellState.hpp:105
PerfData< Scalar > perf_data
Definition: SingleWellState.hpp:143
Definition: BlackoilWellModel.hpp:87
void beginTimeStep()
Definition: TemperatureModel.hpp:768
void init()
Definition: TemperatureModel.hpp:767
const Scalar temperature(size_t) const
Definition: TemperatureModel.hpp:769
void deserialize(Restarter &)
This method restores the complete state of the temperature from disk.
Definition: TemperatureModel.hpp:764
void serialize(Restarter &)
This method writes the complete state of all temperature to the hard disk.
Definition: TemperatureModel.hpp:754
TemperatureModel(Simulator &)
Definition: TemperatureModel.hpp:746
A class which handles sequential implicit solution of the energy equation as specified in by TEMP.
Definition: TemperatureModel.hpp:179
void serialize(Restarter &)
This method writes the complete state of all temperature to the hard disk.
Definition: TemperatureModel.hpp:364
void assembleEquationWell(const Well &well)
Definition: TemperatureModel.hpp:657
void computeFluxTerm(const FluidStateTemp &fsIn, const FluidStateTemp &fsEx, const RateVector &darcyFlux, Evaluation &flux)
Definition: TemperatureModel.hpp:517
std::vector< int > interiorRows_
Definition: TemperatureModel.hpp:733
EnergyVector storage1_
Definition: TemperatureModel.hpp:727
bool converged(const int iter)
Definition: TemperatureModel.hpp:437
void computeStorageTerm(unsigned globI, LhsEval &storage)
Definition: TemperatureModel.hpp:491
void deserialize(Restarter &)
This method restores the complete state of the temperature from disk.
Definition: TemperatureModel.hpp:374
void endTimeStep(WellStateType &wellState)
Informs the temperature model that a time step has just been finished.
Definition: TemperatureModel.hpp:326
const Simulator & simulator_
Definition: TemperatureModel.hpp:726
void computeHeatFluxTerm(const IntensiveQuantitiesTemp &intQuantsIn, const IntensiveQuantitiesTemp &intQuantsEx, const ResidualNBInfo &res_nbinfo, Evaluation &heatFlux)
Definition: TemperatureModel.hpp:546
void beginTimeStep()
Definition: TemperatureModel.hpp:301
void computeWellTemperature(const Well &well, SingleWellState &ws)
Definition: TemperatureModel.hpp:704
void advanceTemperatureFields()
Definition: TemperatureModel.hpp:392
void updateStorageCache()
Definition: TemperatureModel.hpp:378
Scalar scalingFactor_
Definition: TemperatureModel.hpp:734
std::unique_ptr< SpareMatrixEnergyAdapter > energyMatrix_
Definition: TemperatureModel.hpp:731
void assembleEquations()
Definition: TemperatureModel.hpp:567
void solveAndUpdate()
Definition: TemperatureModel.hpp:415
std::vector< MatrixBlockTemp * > diagMatAddress_
Definition: TemperatureModel.hpp:730
std::vector< IntensiveQuantitiesTemp > intQuants_
Definition: TemperatureModel.hpp:728
std::vector< int > overlapRows_
Definition: TemperatureModel.hpp:732
SparseTable< NeighborInfoCPU > neighborInfo_
Definition: TemperatureModel.hpp:729
TemperatureModel(Simulator &simulator)
Definition: TemperatureModel.hpp:227
void init()
Definition: TemperatureModel.hpp:235
Definition: WellState.hpp:66
int numWells() const
Definition: WellState.hpp:99
const SingleWellState< Scalar, IndexTraits > & well(std::size_t well_index) const
Definition: WellState.hpp:290
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:43
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:53
a tag to mark properties as undefined
Definition: propertysystem.hh:38