175class TemperatureModel :
public GenericTemperatureModel<GetPropType<TypeTag, Properties::Grid>,
176 GetPropType<TypeTag, Properties::GridView>,
177 GetPropType<TypeTag, Properties::DofMapper>,
178 GetPropType<TypeTag, Properties::Stencil>,
179 GetPropType<TypeTag, Properties::FluidSystem>,
180 GetPropType<TypeTag, Properties::Scalar>>
182 using BaseType = GenericTemperatureModel<GetPropType<TypeTag, Properties::Grid>,
199 using IndexTraits =
typename FluidSystem::IndexTraitsType;
203 using FluidStateTemp =
typename IntensiveQuantitiesTemp::FluidStateTemp;
204 using Evaluation =
typename IntensiveQuantitiesTemp::EvaluationTemp;
205 using EnergyMatrix =
typename BaseType::EnergyMatrix;
206 using EnergyVector =
typename BaseType::EnergyVector;
207 using MatrixBlockTemp =
typename BaseType::MatrixBlockTemp;
212 struct ResidualNBInfo
222 enum { numPhases = FluidSystem::numPhases };
223 enum { waterPhaseIdx = FluidSystem::waterPhaseIdx };
224 enum { oilPhaseIdx = FluidSystem::oilPhaseIdx };
225 enum { gasPhaseIdx = FluidSystem::gasPhaseIdx };
226 static constexpr int temperatureIdx = 0;
229 explicit TemperatureModel(Simulator& simulator)
230 : BaseType(simulator.vanguard().gridView(),
231 simulator.vanguard().eclState(),
232 simulator.vanguard().cartesianIndexMapper(),
233 simulator.model().dofMapper())
234 , simulator_(simulator)
239 const unsigned int numCells = simulator_.model().numTotalDof();
246 storage1_.resize(numCells);
249 for (
unsigned globI = 0; globI < numCells; ++globI) {
250 this->temperature_[globI] = simulator_.problem().initialFluidState(globI).temperature(0);
254 intQuants_.resize(numCells);
257 const auto& elemMapper = simulator_.model().elementMapper();
258 detail::findOverlapAndInterior(simulator_.vanguard().grid(), elemMapper, overlapRows_, interiorRows_);
264 using NeighborSet = std::set<unsigned>;
265 std::vector<NeighborSet> neighbors(numCells);
266 Stencil stencil(this->gridView_, this->dofMapper_);
267 neighborInfo_.reserve(numCells, 6 * numCells);
268 std::vector<NeighborInfoCPU> loc_nbinfo;
269 for (
const auto& elem : elements(this->gridView_)) {
270 stencil.update(elem);
271 for (
unsigned primaryDofIdx = 0; primaryDofIdx < stencil.numPrimaryDof(); ++primaryDofIdx) {
272 const unsigned myIdx = stencil.globalSpaceIndex(primaryDofIdx);
273 loc_nbinfo.resize(stencil.numDof() - 1);
274 for (
unsigned dofIdx = 0; dofIdx < stencil.numDof(); ++dofIdx) {
275 const unsigned neighborIdx = stencil.globalSpaceIndex(dofIdx);
276 neighbors[myIdx].insert(neighborIdx);
278 const auto scvfIdx = dofIdx - 1;
279 const auto& scvf = stencil.interiorFace(scvfIdx);
280 const Scalar area = scvf.area();
281 Scalar inAlpha = simulator_.problem().thermalHalfTransmissibility(myIdx, neighborIdx);
282 Scalar outAlpha = simulator_.problem().thermalHalfTransmissibility(neighborIdx, myIdx);
283 ResidualNBInfo nbinfo{area, inAlpha, outAlpha};
284 loc_nbinfo[dofIdx - 1] = NeighborInfoCPU{neighborIdx, nbinfo,
nullptr};
287 neighborInfo_.appendRow(loc_nbinfo.begin(), loc_nbinfo.end());
291 energyMatrix_ = std::make_unique<SpareMatrixEnergyAdapter>(simulator_);
292 diagMatAddress_.resize(numCells);
293 energyMatrix_->reserve(neighbors);
294 for (
unsigned globI = 0; globI < numCells; globI++) {
295 const auto& nbInfos = neighborInfo_[globI];
296 diagMatAddress_[globI] = energyMatrix_->blockAddress(globI, globI);
297 for (
auto& nbInfo : nbInfos) {
298 nbInfo.matBlockAddress = energyMatrix_->blockAddress(nbInfo.neighbor, globI);
305 if (!this->doTemp()) {
310 const unsigned int numCells = simulator_.model().numTotalDof();
312 #pragma omp parallel for
314 for (
unsigned globI = 0; globI < numCells; ++globI) {
315 intQuants_[globI].setFluidState(simulator_.model().intensiveQuantities(globI, 0).fluidState());
316 intQuants_[globI].updateTemperature_(simulator_.problem(), globI, 0);
317 intQuants_[globI].updateEnergyQuantities_(simulator_.problem(), globI, 0);
319 updateStorageCache();
321 const int nw = simulator_.problem().wellModel().wellState().numWells();
322 this->energy_rates_.resize(nw, 0.0);
330 if (!this->doTemp()) {
335 const unsigned int numCells = simulator_.model().numTotalDof();
337 #pragma omp parallel for
339 for (
unsigned globI = 0; globI < numCells; ++globI) {
340 intQuants_[globI].setFluidState(simulator_.model().intensiveQuantities(globI, 0).fluidState());
341 intQuants_[globI].updateTemperature_(simulator_.problem(), globI, 0);
342 intQuants_[globI].updateEnergyQuantities_(simulator_.problem(), globI, 0);
344 advanceTemperatureFields();
347 const int nw = wellState.numWells();
348 for (
auto wellID = 0*nw; wellID < nw; ++wellID) {
349 auto& ws = wellState.well(wellID);
350 ws.energy_rate = this->energy_rates_[wellID];
354 const auto& wellPtrs = simulator_.problem().wellModel().localNonshutWells();
355 for (
const auto& wellPtr : wellPtrs) {
356 auto& ws = wellState.well(wellPtr->name());
357 this->computeWellTemperature(*wellPtr, ws);
365 template <
class Restarter>
375 template <
class Restarter>
380 void updateStorageCache()
383 const unsigned int numCells = simulator_.model().numTotalDof();
385 #pragma omp parallel for
387 for (
unsigned globI = 0; globI < numCells; ++globI) {
388 Scalar storage = 0.0;
389 computeStorageTerm(globI, storage);
390 storage1_[globI] = storage;
394 void advanceTemperatureFields()
396 const int max_iter = 20;
397 const int min_iter = 1;
398 bool is_converged =
false;
400 for (
int iter = 0; iter < max_iter; ++iter) {
402 if (iter >= min_iter && converged(iter)) {
410 fmt::format(fmt::runtime(
"Temperature model (TEMP): Newton did not converge after {} iterations. \n"
411 "The Simulator will continue to the next step with an unconverged solution."),
417 void solveAndUpdate()
419 const unsigned int numCells = simulator_.model().numTotalDof();
420 EnergyVector dx(numCells);
421 bool conv = this->linearSolve_(this->energyMatrix_->istlMatrix(), dx, this->energyVector_);
423 if (simulator_.gridView().comm().rank() == 0) {
424 OpmLog::warning(
"Temp model: Linear solver did not converge. Temperature values not updated.");
429 #pragma omp parallel for
431 for (
unsigned globI = 0; globI < numCells; ++globI) {
432 this->temperature_[globI] -= std::clamp(dx[globI][0], -this->maxTempChange_, this->maxTempChange_);
433 intQuants_[globI].updateTemperature_(simulator_.problem(), globI, 0);
434 intQuants_[globI].updateEnergyQuantities_(simulator_.problem(), globI, 0);
439 bool converged(
const int iter)
441 Scalar dt = simulator_.timeStepSize();
442 Scalar maxNorm = 0.0;
443 Scalar sumNorm = 0.0;
444 const auto tolerance_cnv_energy_strict = Parameters::Get<Parameters::ToleranceCnvEnergy<Scalar>>();
445 const auto& elemMapper = simulator_.model().elementMapper();
446 const IsNumericalAquiferCell isNumericalAquiferCell(simulator_.gridView().grid());
448 Scalar sum_pv_not_converged = 0.0;
449 for (
const auto& elem : elements(simulator_.gridView(), Dune::Partitions::interior)) {
450 unsigned globI = elemMapper.index(elem);
451 const auto pvValue = simulator_.problem().referencePorosity(globI, 0)
452 * simulator_.model().dofTotalVolume(globI);
454 const Scalar scaled_norm = dt * std::abs(this->energyVector_[globI])/ pvValue;
455 maxNorm = max(maxNorm, scaled_norm);
456 sumNorm += scaled_norm;
457 if (!isNumericalAquiferCell(elem)) {
458 if (scaled_norm > tolerance_cnv_energy_strict) {
459 sum_pv_not_converged += pvValue;
464 maxNorm = simulator_.gridView().comm().max(maxNorm);
465 sumNorm = simulator_.gridView().comm().sum(sumNorm);
466 sum_pv = simulator_.gridView().comm().sum(sum_pv);
470 sum_pv_not_converged = simulator_.gridView().comm().sum(sum_pv_not_converged);
471 Scalar relaxed_max_pv_fraction = Parameters::Get<Parameters::RelaxedMaxPvFraction<Scalar>>();
472 const bool relax = (sum_pv_not_converged / sum_pv) < relaxed_max_pv_fraction;
473 const auto tolerance_energy_balance = relax? Parameters::Get<Parameters::ToleranceEnergyBalanceRelaxed<Scalar>>():
474 Parameters::
Get<Parameters::ToleranceEnergyBalance<Scalar>>();
475 const bool tolerance_cnv_energy = relax? Parameters::Get<Parameters::ToleranceCnvEnergyRelaxed<Scalar>>():
476 tolerance_cnv_energy_strict;
478 const auto msg = fmt::format(fmt::runtime(
"Temperature model (TEMP): Newton iter {}: "
479 "CNV(E): {:.1e}, EB: {:.1e}"),
480 iter, maxNorm, sumNorm);
482 if (maxNorm < tolerance_cnv_energy && sumNorm < tolerance_energy_balance) {
484 fmt::format(fmt::runtime(
"Temperature model (TEMP): Newton converged after {} iterations"),
492 template<
class LhsEval>
493 void computeStorageTerm(
unsigned globI, LhsEval& storage)
495 const auto& intQuants = intQuants_[globI];
496 const auto& poro = getValue(simulator_.model().intensiveQuantities(globI, 0).porosity());
498 const auto& fs = intQuants.fluidStateTemp();
499 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++ phaseIdx) {
500 if (!FluidSystem::phaseIsActive(phaseIdx)) {
504 const auto& u = decay<LhsEval>(fs.internalEnergy(phaseIdx));
505 const auto& S = decay<LhsEval>(fs.saturation(phaseIdx));
506 const auto& rho = decay<LhsEval>(fs.density(phaseIdx));
508 storage += poro*S*u*rho;
512 const Scalar rockFraction = intQuants.rockFraction();
513 const auto& uRock = decay<LhsEval>(intQuants.rockInternalEnergy());
514 storage += rockFraction*uRock;
515 storage*= scalingFactor_;
518 template <
class RateVector>
519 void computeFluxTerm(
const FluidStateTemp& fsIn,
520 const FluidStateTemp& fsEx,
521 const RateVector& darcyFlux,
524 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++ phaseIdx) {
525 if (!FluidSystem::phaseIsActive(phaseIdx)) {
529 const unsigned activeCompIdx =
530 FluidSystem::canonicalToActiveCompIdx(FluidSystem::solventComponentIndex(phaseIdx));
531 bool inIsUp = darcyFlux[activeCompIdx] > 0;
532 const auto& fs = inIsUp ? fsIn : fsEx;
534 flux += fs.enthalpy(phaseIdx)
535 * fs.density(phaseIdx)
536 * getValue(darcyFlux[activeCompIdx]);
539 flux += getValue(fs.enthalpy(phaseIdx))
540 * getValue(fs.density(phaseIdx))
541 * getValue(darcyFlux[activeCompIdx]);
544 flux *= scalingFactor_;
547 template <
class Res
idualNBInfo>
548 void computeHeatFluxTerm(
const IntensiveQuantitiesTemp& intQuantsIn,
549 const IntensiveQuantitiesTemp& intQuantsEx,
550 const ResidualNBInfo& res_nbinfo,
551 Evaluation& heatFlux)
553 short interiorDofIdx = 0;
554 short exteriorDofIdx = 1;
555 EnergyModule::ExtensiveQuantities::updateEnergy(heatFlux,
561 intQuantsIn.fluidStateTemp(),
562 intQuantsEx.fluidStateTemp(),
565 res_nbinfo.faceArea);
566 heatFlux *= scalingFactor_*res_nbinfo.faceArea;
569 void assembleEquations()
571 const unsigned int numCells = simulator_.model().numTotalDof();
572 for (
unsigned globI = 0; globI < numCells; ++globI) {
573 this->energyVector_[globI] = 0.0;
574 energyMatrix_->clearRow(globI, 0.0);
576 Scalar dt = simulator_.timeStepSize();
578 #pragma omp parallel for
580 for (
unsigned globI = 0; globI < numCells; ++globI) {
581 MatrixBlockTemp bMat;
582 Scalar volume = simulator_.model().dofTotalVolume(globI);
583 Scalar storefac = volume / dt;
584 Evaluation storage = 0.0;
585 computeStorageTerm(globI, storage);
586 this->energyVector_[globI] += storefac * ( getValue(storage) - storage1_[globI] );
587 bMat[0][0] = storefac * storage.derivative(temperatureIdx);
588 *diagMatAddress_[globI] += bMat;
591 const auto& floresInfo = this->simulator_.problem().model().linearizer().getFloresInfo();
592 const bool enableDriftCompensation = Parameters::Get<Parameters::EnableDriftCompensationTemp>();
593 const auto& problem = simulator_.problem();
596 #pragma omp parallel for
598 for (
unsigned globI = 0; globI < numCells; ++globI) {
599 const auto& nbInfos = neighborInfo_[globI];
600 const auto& floresInfos = floresInfo[globI];
602 const auto& intQuantsIn = intQuants_[globI];
603 MatrixBlockTemp bMat;
604 for (
const auto& nbInfo : nbInfos) {
605 unsigned globJ = nbInfo.neighbor;
606 const auto& intQuantsEx = intQuants_[globJ];
607 assert(globJ != globI);
608 const auto& darcyflux = floresInfos[loc].flow;
610 Evaluation flux = 0.0;
611 computeFluxTerm(intQuantsIn.fluidStateTemp(), intQuantsEx.fluidStateTemp(), darcyflux, flux);
613 Evaluation heatFlux = 0.0;
614 computeHeatFluxTerm(intQuantsIn, intQuantsEx, nbInfo.res_nbinfo, heatFlux);
616 this->energyVector_[globI] += getValue(heatFlux);
617 bMat[0][0] = heatFlux.derivative(temperatureIdx);
618 *diagMatAddress_[globI] += bMat;
621 *nbInfo.matBlockAddress += bMat;
625 if (enableDriftCompensation) {
626 auto dofDriftRate = problem.drift()[globI]/dt;
627 const auto& fs = intQuantsIn.fluidStateTemp();
628 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++ phaseIdx) {
629 const unsigned activeCompIdx =
630 FluidSystem::canonicalToActiveCompIdx(FluidSystem::solventComponentIndex(phaseIdx));
631 auto drift_hrate = dofDriftRate[activeCompIdx]*getValue(fs.enthalpy(phaseIdx)) * getValue(fs.density(phaseIdx)) / getValue(fs.invB(phaseIdx));
632 this->energyVector_[globI] -= drift_hrate*scalingFactor_;
638 const auto& wellPtrs = simulator_.problem().wellModel().localNonshutWells();
639 for (
const auto& wellPtr : wellPtrs) {
640 this->assembleEquationWell(*wellPtr);
659 void assembleEquationWell(
const Well& well)
661 const auto& eclWell = well.wellEcl();
662 std::size_t well_index = simulator_.problem().wellModel().wellState().index(well.name()).value();
663 const auto& ws = simulator_.problem().wellModel().wellState().well(well_index);
664 this->energy_rates_[well_index] = 0.0;
665 MatrixBlockTemp bMat;
666 for (std::size_t i = 0; i < ws.perf_data.size(); ++i) {
667 const auto globI = ws.perf_data.cell_index[i];
668 auto fs = intQuants_[globI].fluidStateTemp();
669 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++ phaseIdx) {
670 if (!FluidSystem::phaseIsActive(phaseIdx)) {
674 Evaluation rate = well.volumetricSurfaceRateForConnection(globI, phaseIdx);
675 if (rate > 0 && eclWell.isInjector()) {
676 fs.setTemperature(eclWell.inj_temperature());
677 const auto& rho = FluidSystem::density(fs, phaseIdx, fs.pvtRegionIndex());
678 fs.setDensity(phaseIdx, rho);
679 const auto& h = FluidSystem::enthalpy(fs, phaseIdx, fs.pvtRegionIndex());
680 fs.setEnthalpy(phaseIdx, h);
681 rate *= getValue(fs.enthalpy(phaseIdx)) * getValue(fs.density(phaseIdx)) / getValue(fs.invB(phaseIdx));
683 const Evaluation d = 1.0 - fs.Rv() * fs.Rs();
684 if (phaseIdx == gasPhaseIdx && d > 0) {
685 const auto& oilrate = well.volumetricSurfaceRateForConnection(globI, oilPhaseIdx);
686 rate -= oilrate * getValue(fs.Rs());
689 if (phaseIdx == oilPhaseIdx && d > 0) {
690 const auto& gasrate = well.volumetricSurfaceRateForConnection(globI, gasPhaseIdx);
691 rate -= gasrate * getValue(fs.Rv());
694 rate *= fs.enthalpy(phaseIdx) * getValue(fs.density(phaseIdx)) / getValue(fs.invB(phaseIdx));
696 this->energy_rates_[well_index] += getValue(rate);
697 rate *= scalingFactor_;
698 this->energyVector_[globI] -= getValue(rate);
699 bMat[0][0] = -rate.derivative(temperatureIdx);
700 *diagMatAddress_[globI] += bMat;
705 template<
class Well,
class SingleWellState>
706 void computeWellTemperature(
const Well& well, SingleWellState& ws)
708 if (well.isInjector()) {
709 if (ws.status != WellStatus::STOP) {
710 ws.temperature = well.wellEcl().inj_temperature();
714 const int np = simulator_.problem().wellModel().wellState().numPhases();
715 std::array<Scalar,2> weighted{0.0,0.0};
716 auto& [weighted_temperature, total_weight] = weighted;
717 for (std::size_t i = 0; i < ws.perf_data.size(); ++i) {
718 const auto globI = ws.perf_data.cell_index[i];
719 const auto& fs = intQuants_[globI].fluidStateTemp();
720 Scalar weight_factor = simulator_.problem().wellModel().computeTemperatureWeightFactor(i, np, fs, ws);
721 total_weight += weight_factor;
722 weighted_temperature += weight_factor * fs.temperature(0).value();
725 ws.temperature = weighted_temperature / total_weight;
728 const Simulator& simulator_;
729 EnergyVector storage1_;
730 std::vector<IntensiveQuantitiesTemp> intQuants_;
731 SparseTable<NeighborInfoCPU> neighborInfo_{};
732 std::vector<MatrixBlockTemp*> diagMatAddress_{};
733 std::unique_ptr<SpareMatrixEnergyAdapter> energyMatrix_;
734 std::vector<int> overlapRows_;
735 std::vector<int> interiorRows_;
736 Scalar scalingFactor_{1.0};