TemperatureModel.hpp
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3/*
4 This file is part of the Open Porous Media project (OPM).
5
6 OPM is free software: you can redistribute it and/or modify
7 it under the terms of the GNU General Public License as published by
8 the Free Software Foundation, either version 2 of the License, or
9 (at your option) any later version.
10
11 OPM is distributed in the hope that it will be useful,
12 but WITHOUT ANY WARRANTY; without even the implied warranty of
13 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14 GNU General Public License for more details.
15
16 You should have received a copy of the GNU General Public License
17 along with OPM. If not, see <http://www.gnu.org/licenses/>.
18
19 Consult the COPYING file in the top-level source directory of this
20 module for the precise wording of the license and the list of
21 copyright holders.
22*/
28#ifndef OPM_TEMPERATURE_MODEL_HPP
29#define OPM_TEMPERATURE_MODEL_HPP
30
31#include <opm/common/OpmLog/OpmLog.hpp>
32#include <opm/common/utility/gpuDecorators.hpp>
33
36
43
44#include <algorithm>
45#include <cassert>
46#include <cmath>
47#include <cstddef>
48#include <string>
49#include <vector>
50
51namespace Opm::Properties {
52
53template<class TypeTag, class MyTypeTag>
56};
57
58} // namespace Opm::Properties
59
60namespace Opm {
61
62template<typename Scalar, typename IndexTraits> class WellState;
63
64
65template <class TypeTag>
67{
74
75
76
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;
84
85 public:
86 using EvaluationTemp = DenseAd::Evaluation<Scalar, 1>;
87 using FluidStateTemp = BlackOilFluidState<EvaluationTemp,
88 FluidSystem,
89 true, //store temperature
90 true, // store enthalpy
91 compositionSwitchEnabled,
92 enableVapwat,
93 enableBrine,
94 enableSaltPrecipitation,
95 enableDisgasInWater,
96 enableSolvent,
97 Indices::numPhases>;
98
99
100
101 OPM_HOST_DEVICE void updateTemperature_(const Problem& problem,
102 unsigned globalDofIdx,
103 unsigned timeIdx)
104 {
105 const EvaluationTemp T = EvaluationTemp::createVariable(problem.temperature(globalDofIdx, timeIdx), 0);
106 fluidState_.setTemperature(T);
107 }
108
109 OPM_HOST_DEVICE void updateEnergyQuantities_(const Problem& problem,
110 const unsigned globalSpaceIdx,
111 const unsigned timeIdx)
112 {
113 // compute the specific enthalpy of the fluids, the specific enthalpy of the rock
114 // and the thermal conductivity coefficients
115 for (int phaseIdx = 0; phaseIdx < numPhases; ++ phaseIdx) {
116 if (!FluidSystem::phaseIsActive(phaseIdx)) {
117 continue;
118 }
119
120 const auto& h = FluidSystem::enthalpy(fluidState_, phaseIdx, problem.pvtRegionIndex(globalSpaceIdx));
121 fluidState_.setEnthalpy(phaseIdx, h);
122 }
123
124 const auto& solidEnergyLawParams = problem.solidEnergyLawParams(globalSpaceIdx, timeIdx);
125 rockInternalEnergy_ = SolidEnergyLaw::solidInternalEnergy(solidEnergyLawParams, fluidState_);
126
127 const auto& thermalConductionLawParams = problem.thermalConductionLawParams(globalSpaceIdx, timeIdx);
128 totalThermalConductivity_ = ThermalConductionLaw::thermalConductivity(thermalConductionLawParams, fluidState_);
129
130 // Retrieve the rock fraction from the problem
131 // Usually 1 - porosity, but if pvmult is used to modify porosity
132 // we will apply the same multiplier to the rock fraction
133 // i.e. pvmult*(1 - porosity) and thus interpret multpv as a volume
134 // multiplier. This is to avoid negative rock volume for pvmult*porosity > 1
135 rockFraction_ = problem.rockFraction(globalSpaceIdx, timeIdx);
136 }
137
138 OPM_HOST_DEVICE const EvaluationTemp& rockInternalEnergy() const
139 { return rockInternalEnergy_; }
140
141 OPM_HOST_DEVICE const EvaluationTemp& totalThermalConductivity() const
142 { return totalThermalConductivity_; }
143
144 OPM_HOST_DEVICE const Scalar& rockFraction() const
145 { return rockFraction_; }
146
147 OPM_HOST_DEVICE const FluidStateTemp& fluidStateTemp() const
148 { return fluidState_; }
149
150 template <class FluidState>
151 OPM_HOST_DEVICE void setFluidState(const FluidState& fs)
152 {
153 // copy the needed part of the fluid state
154 fluidState_.setPvtRegionIndex(fs.pvtRegionIndex());
155 fluidState_.setRs(getValue(fs.Rs()));
156 fluidState_.setRv(getValue(fs.Rv()));
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)));
162 }
163 }
164
165protected:
170};
171
177template <class TypeTag, bool enableTempV = getPropValue<TypeTag, Properties::EnergyModuleType>() == EnergyModules::SequentialImplicitThermal >
178class TemperatureModel : public GenericTemperatureModel<GetPropType<TypeTag, Properties::Grid>,
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>>
184{
202 using IndexTraits = typename FluidSystem::IndexTraitsType;
204
206 using FluidStateTemp = typename IntensiveQuantitiesTemp::FluidStateTemp;
207 using Evaluation = typename IntensiveQuantitiesTemp::EvaluationTemp;
208 using EnergyMatrix = typename BaseType::EnergyMatrix;
209 using EnergyVector = typename BaseType::EnergyVector;
210 using MatrixBlockTemp = typename BaseType::MatrixBlockTemp;
211 using SpareMatrixEnergyAdapter = typename Linear::IstlSparseMatrixAdapter<MatrixBlockTemp>;
212
213
214 //using ResidualNBInfo = typename LocalResidual::ResidualNBInfo;
215 struct ResidualNBInfo
216 {
217 double faceArea;
218 double inAlpha;
219 double outAlpha;
220 };
221
223
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;
230
231public:
232 explicit TemperatureModel(Simulator& simulator)
233 : BaseType(simulator.vanguard().gridView(),
234 simulator.vanguard().eclState(),
235 simulator.vanguard().cartesianIndexMapper(),
236 simulator.model().dofMapper())
237 , simulator_(simulator)
238 {}
239
240 void init()
241 {
242 const unsigned int numCells = simulator_.model().numTotalDof();
243 this->doInit(numCells);
244
245 if (!this->doTemp())
246 return;
247
248 // we need the storage term at start of the iteration (timeIdx = 1)
249 storage1_.resize(numCells);
250
251 // set the initial temperature
252 for (unsigned globI = 0; globI < numCells; ++globI) {
253 this->temperature_[globI] = simulator_.problem().initialFluidState(globI).temperature(0);
254 }
255 // keep a copy of the intensive quantities to simplify the update during
256 // the newton iterations
257 intQuants_.resize(numCells);
258
259 // find and store the overlap cells
260 const auto& elemMapper = simulator_.model().elementMapper();
262
263 // set the scaling factor
264 scalingFactor_ = getPropValue<TypeTag, Properties::BlackOilEnergyScalingFactor>();
265
266 // find the sparsity pattern of the temperature matrix
267 using NeighborSet = std::set<unsigned>;
268 std::vector<NeighborSet> neighbors(numCells);
269 Stencil stencil(this->gridView_, this->dofMapper_);
270 neighborInfo_.reserve(numCells, 6 * 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); // Do not include the primary dof in neighborInfo_
277 for (unsigned dofIdx = 0; dofIdx < stencil.numDof(); ++dofIdx) {
278 const unsigned neighborIdx = stencil.globalSpaceIndex(dofIdx);
279 neighbors[myIdx].insert(neighborIdx);
280 if (dofIdx > 0) {
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};
288 }
289 }
290 neighborInfo_.appendRow(loc_nbinfo.begin(), loc_nbinfo.end());
291 }
292 }
293 // allocate matrix for storing the Jacobian of the temperature residual
294 energyMatrix_ = std::make_unique<SpareMatrixEnergyAdapter>(simulator_);
295 diagMatAddress_.resize(numCells);
296 energyMatrix_->reserve(neighbors);
297 for (unsigned globI = 0; globI < numCells; globI++) {
298 const auto& nbInfos = neighborInfo_[globI];
299 diagMatAddress_[globI] = energyMatrix_->blockAddress(globI, globI);
300 for (auto& nbInfo : nbInfos) {
301 nbInfo.matBlockAddress = energyMatrix_->blockAddress(nbInfo.neighbor, globI);
302 }
303 }
304 }
305
307 {
308 if (!this->doTemp()) {
309 return;
310 }
311
312 // We use the specialized intensive quantities here with only the temperature derivative
313 const unsigned int numCells = simulator_.model().numTotalDof();
314 #ifdef _OPENMP
315 #pragma omp parallel for
316 #endif
317 for (unsigned globI = 0; globI < numCells; ++globI) {
318 intQuants_[globI].setFluidState(simulator_.model().intensiveQuantities(globI, /*timeIdx*/ 0).fluidState());
319 intQuants_[globI].updateTemperature_(simulator_.problem(), globI, /*timeIdx*/ 0);
320 intQuants_[globI].updateEnergyQuantities_(simulator_.problem(), globI, /*timeIdx*/ 0);
321 }
323
324 const int nw = simulator_.problem().wellModel().wellState().numWells();
325 this->energy_rates_.resize(nw, 0.0);
326 }
327
331 void endTimeStep(WellStateType& wellState)
332 {
333 if (!this->doTemp()) {
334 return;
335 }
336 OPM_TIMEBLOCK(TemperatureModel_endTimeStep);
337
338 // We use the specialized intensive quantities here with only the temperature derivative
339 const unsigned int numCells = simulator_.model().numTotalDof();
340 #ifdef _OPENMP
341 #pragma omp parallel for
342 #endif
343 for (unsigned globI = 0; globI < numCells; ++globI) {
344 intQuants_[globI].setFluidState(simulator_.model().intensiveQuantities(globI, /*timeIdx*/ 0).fluidState());
345 intQuants_[globI].updateTemperature_(simulator_.problem(), globI, /*timeIdx*/ 0);
346 intQuants_[globI].updateEnergyQuantities_(simulator_.problem(), globI, /*timeIdx*/ 0);
347 }
349
350 // update energy_rates
351 const int nw = wellState.numWells();
352 for (auto wellID = 0*nw; wellID < nw; ++wellID) {
353 auto& ws = wellState.well(wellID);
354 ws.energy_rate = this->energy_rates_[wellID];
355 }
356
357 // Update well temperature
358 const auto& wellPtrs = simulator_.problem().wellModel().localNonshutWells();
359 for (const auto& wellPtr : wellPtrs) {
360 auto& ws = wellState.well(wellPtr->name());
361 this->computeWellTemperature(*wellPtr, ws);
362 }
363 }
364
369 template <class Restarter>
370 void serialize(Restarter&)
371 { /* not implemented */ }
372
379 template <class Restarter>
380 void deserialize(Restarter&)
381 { /* not implemented */ }
382
383protected:
385 {
386 // we need the storage term at start of the iteration (timeIdx = 1)
387 const unsigned int numCells = simulator_.model().numTotalDof();
388 #ifdef _OPENMP
389 #pragma omp parallel for
390 #endif
391 for (unsigned globI = 0; globI < numCells; ++globI) {
392 Scalar storage = 0.0;
393 computeStorageTerm(globI, storage);
394 storage1_[globI] = storage;
395 }
396 }
397
399 {
400 OPM_TIMEBLOCK(TemperatureModel_advanceTemperatureFields);
401 const int max_iter = 20;
402 const int min_iter = 1;
403 bool is_converged = false;
404 // solve using Newton
405 for (int iter = 0; iter < max_iter; ++iter) {
407 if (iter >= min_iter && converged(iter)) {
408 is_converged = true;
409 break;
410 }
412 }
413 if (!is_converged) {
414 const auto msg =
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."),
417 max_iter);
418 OpmLog::debug(msg);
419 }
420 }
421
423 {
424 OPM_TIMEBLOCK(TemperatureModel_solveAndUpdate);
425 const unsigned int numCells = simulator_.model().numTotalDof();
426 EnergyVector dx(numCells);
427 bool conv = this->linearSolve_(this->energyMatrix_->istlMatrix(), dx, this->energyVector_);
428 if (!conv) {
429 if (simulator_.gridView().comm().rank() == 0) {
430 OpmLog::warning("Temp model: Linear solver did not converge. Temperature values not updated.");
431 }
432 }
433 else {
434 OPM_TIMEBLOCK(TemperatureModel_solveAndUpdate_update);
435 #ifdef _OPENMP
436 #pragma omp parallel for
437 #endif
438 for (unsigned globI = 0; globI < numCells; ++globI) {
439 this->temperature_[globI] -= std::clamp(dx[globI][0], -this->maxTempChange_, this->maxTempChange_);
440 intQuants_[globI].updateTemperature_(simulator_.problem(), globI, /*timeIdx*/ 0);
441 intQuants_[globI].updateEnergyQuantities_(simulator_.problem(), globI, /*timeIdx*/ 0);
442 }
443 }
444 }
445
446 bool converged(const int iter)
447 {
448 OPM_TIMEBLOCK(TemperatureModel_converged);
449 Scalar dt = simulator_.timeStepSize();
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();
454 const IsNumericalAquiferCell isNumericalAquiferCell(simulator_.gridView().grid());
455 Scalar sum_pv = 0.0;
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, /*timeIdx=*/0)
460 * simulator_.model().dofTotalVolume(globI);
461
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;
468 }
469 sum_pv += pvValue;
470 }
471 }
472 {
473 OPM_TIMEBLOCK(TemperatureModel_converged_communicate);
474
475 maxNorm = simulator_.gridView().comm().max(maxNorm);
476 sumNorm = simulator_.gridView().comm().sum(sumNorm);
477 sum_pv = simulator_.gridView().comm().sum(sum_pv);
478 sumNorm /= sum_pv;
479
480 // Use relaxed tolerance if the fraction of unconverged cells porevolume is less than relaxed_max_pv_fraction
481 sum_pv_not_converged = simulator_.gridView().comm().sum(sum_pv_not_converged);
482 }
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;
489
490 const auto msg = fmt::format(fmt::runtime("Temperature model (TEMP): Newton iter {}: "
491 "CNV(E): {:.1e}, EB: {:.1e}"),
492 iter, maxNorm, sumNorm);
493 OpmLog::debug(msg);
494 if (maxNorm < tolerance_cnv_energy && sumNorm < tolerance_energy_balance) {
495 const auto msg2 =
496 fmt::format(fmt::runtime("Temperature model (TEMP): Newton converged after {} iterations"),
497 iter);
498 OpmLog::debug(msg2);
499 return true;
500 }
501 return false;
502 }
503
504 template< class LhsEval>
505 void computeStorageTerm(unsigned globI, LhsEval& storage)
506 {
507 const auto& intQuants = intQuants_[globI];
508 const auto& poro = getValue(simulator_.model().intensiveQuantities(globI, /*timeIdx*/ 0).porosity());
509 // accumulate the internal energy of the fluids
510 const auto& fs = intQuants.fluidStateTemp();
511 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++ phaseIdx) {
512 if (!FluidSystem::phaseIsActive(phaseIdx)) {
513 continue;
514 }
515
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));
519
520 storage += poro*S*u*rho;
521 }
522
523 // add the internal energy of the rock
524 const Scalar rockFraction = intQuants.rockFraction();
525 const auto& uRock = decay<LhsEval>(intQuants.rockInternalEnergy());
526 storage += rockFraction*uRock;
527 storage*= scalingFactor_;
528 }
529
530 template <class RateVector>
531 void computeFluxTerm(const FluidStateTemp& fsIn,
532 const FluidStateTemp& fsEx,
533 const RateVector& darcyFlux,
534 Evaluation& flux)
535 {
536 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++ phaseIdx) {
537 if (!FluidSystem::phaseIsActive(phaseIdx)) {
538 continue;
539 }
540
541 const unsigned activeCompIdx =
542 FluidSystem::canonicalToActiveCompIdx(FluidSystem::solventComponentIndex(phaseIdx));
543 bool inIsUp = darcyFlux[activeCompIdx] > 0;
544 const auto& fs = inIsUp ? fsIn : fsEx;
545 if (inIsUp) {
546 flux += fs.enthalpy(phaseIdx)
547 * fs.density(phaseIdx)
548 * getValue(darcyFlux[activeCompIdx]);
549 }
550 else {
551 flux += getValue(fs.enthalpy(phaseIdx))
552 * getValue(fs.density(phaseIdx))
553 * getValue(darcyFlux[activeCompIdx]);
554 }
555 }
556 flux *= scalingFactor_;
557 }
558
559 template < class ResidualNBInfo>
561 const IntensiveQuantitiesTemp& intQuantsEx,
562 const ResidualNBInfo& res_nbinfo,
563 Evaluation& heatFlux)
564 {
565 short interiorDofIdx = 0; // NB
566 short exteriorDofIdx = 1; // NB
567 EnergyModule::ExtensiveQuantities::updateEnergy(heatFlux,
568 interiorDofIdx, // focusDofIndex,
569 interiorDofIdx,
570 exteriorDofIdx,
571 intQuantsIn,
572 intQuantsEx,
573 intQuantsIn.fluidStateTemp(),
574 intQuantsEx.fluidStateTemp(),
575 res_nbinfo.inAlpha,
576 res_nbinfo.outAlpha,
577 res_nbinfo.faceArea);
578 heatFlux *= scalingFactor_*res_nbinfo.faceArea;
579 }
580
582 {
583 OPM_TIMEBLOCK(TemperatureModel_assembleEquations);
584
585 const unsigned int numCells = simulator_.model().numTotalDof();
586 for (unsigned globI = 0; globI < numCells; ++globI) {
587 this->energyVector_[globI] = 0.0;
588 energyMatrix_->clearRow(globI, 0.0);
589 }
590 const Scalar dt = simulator_.timeStepSize();
591
592 // Storage term
593 {
594 OPM_TIMEBLOCK(TemperatureModel_assembleEquations_storage);
595 #ifdef _OPENMP
596 #pragma omp parallel for
597 #endif
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;
603 computeStorageTerm(globI, storage);
604 this->energyVector_[globI] += storefac * ( getValue(storage) - storage1_[globI] );
605 bMat[0][0] = storefac * storage.derivative(temperatureIdx);
606 *diagMatAddress_[globI] += bMat;
607 }
608 }
609
610 // Flux term
611 {
612 OPM_TIMEBLOCK(TemperatureModel_assembleEquations_flux);
613 const auto& floresInfo = this->simulator_.problem().model().linearizer().getFloresInfo();
614 const bool enableDriftCompensation = Parameters::Get<Parameters::EnableDriftCompensationTemp>();
615 const auto& problem = simulator_.problem();
616
617 #ifdef _OPENMP
618 #pragma omp parallel for
619 #endif
620 for (unsigned globI = 0; globI < numCells; ++globI) {
621 const auto& nbInfos = neighborInfo_[globI];
622 const auto& floresInfos = floresInfo[globI];
623 int loc = 0;
624 const auto& intQuantsIn = intQuants_[globI];
625 MatrixBlockTemp bMat;
626 for (const auto& nbInfo : nbInfos) {
627 unsigned globJ = nbInfo.neighbor;
628 const auto& intQuantsEx = intQuants_[globJ];
629 assert(globJ != globI);
630 const auto& darcyflux = floresInfos[loc].flow;
631 // compute convective flux
632 Evaluation flux = 0.0;
633 computeFluxTerm(intQuantsIn.fluidStateTemp(), intQuantsEx.fluidStateTemp(), darcyflux, flux);
634 // compute conductive flux
635 Evaluation heatFlux = 0.0;
636 computeHeatFluxTerm(intQuantsIn, intQuantsEx, nbInfo.res_nbinfo, heatFlux);
637 heatFlux += flux;
638 this->energyVector_[globI] += getValue(heatFlux);
639 bMat[0][0] = heatFlux.derivative(temperatureIdx);
640 *diagMatAddress_[globI] += bMat;
641 bMat *= -1.0;
642 //SparseAdapter syntax: jacobian_->addToBlock(globJ, globI, bMat);
643 *nbInfo.matBlockAddress += bMat;
644 loc++;
645 }
646
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));
654 this->energyVector_[globI] -= drift_hrate*scalingFactor_;
655 }
656 }
657 }
658 }
659
660 // Well terms
661 {
662 OPM_TIMEBLOCK(TemperatureModel_assembleEquations_wells);
663 const auto& wellPtrs = simulator_.problem().wellModel().localNonshutWells();
664 for (const auto& wellPtr : wellPtrs) {
665 this->assembleEquationWell(*wellPtr);
666 }
667 }
668
669 // For standard ilu0 we need to set the overlapping cells to identity. But since we use "paroverilu0"
670 // here this is not needed.
671 // if (simulator_.gridView().comm().size() > 1) {
672 // // Set dirichlet conditions for overlapping cells
673 // // loop over precalculated overlap rows and columns
674 // for (const auto row : overlapRows_) {
675 // // Zero out row.
676 // (*this->energyMatrix_)[row] = 0.0;
677 //
678 // //diagonal block set to diag(1.0).
679 // (*this->energyMatrix_)[row][row][0][0] = 1.0;
680 // }
681 //}
682 }
683
684 template<class Well>
685 void assembleEquationWell(const Well& well)
686 {
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);
690 this->energy_rates_[well_index] = 0.0;
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];
694 auto fs = intQuants_[globI].fluidStateTemp(); //copy to make it possible to change the temp in the injector
695 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++ phaseIdx) {
696 if (!FluidSystem::phaseIsActive(phaseIdx)) {
697 continue;
698 }
699
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));
708 } else {
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());
713 rate /= d;
714 }
715 if (phaseIdx == oilPhaseIdx && d > 0) {
716 const auto& gasrate = well.volumetricSurfaceRateForConnection(globI, gasPhaseIdx);
717 rate -= gasrate * getValue(fs.Rv());
718 rate /= d;
719 }
720 rate *= fs.enthalpy(phaseIdx) * getValue(fs.density(phaseIdx)) / getValue(fs.invB(phaseIdx));
721 }
722 this->energy_rates_[well_index] += getValue(rate);
723 rate *= scalingFactor_;
724 this->energyVector_[globI] -= getValue(rate);
725 bMat[0][0] = -rate.derivative(temperatureIdx);
726 *diagMatAddress_[globI] += bMat;
727 }
728 }
729 }
730
731 template<class Well, class SingleWellState>
732 void computeWellTemperature(const Well& well, SingleWellState& ws)
733 {
734 if (well.isInjector()) {
735 if (ws.status != WellStatus::STOP) {
736 ws.temperature = well.wellEcl().inj_temperature();
737 return;
738 }
739 }
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(/*phaseIdx*/0).value();
749 }
750 //simulator_.gridView().comm().sum(weighted.data(), 2);
751 ws.temperature = weighted_temperature / total_weight;
752 }
753
754 const Simulator& simulator_;
755 EnergyVector storage1_;
756 std::vector<IntensiveQuantitiesTemp> intQuants_;
758 std::vector<MatrixBlockTemp*> diagMatAddress_{};
759 std::unique_ptr<SpareMatrixEnergyAdapter> energyMatrix_;
760 std::vector<int> overlapRows_;
761 std::vector<int> interiorRows_;
762 Scalar scalingFactor_{1.0};
763};
764
765// need for the old linearizer
766template <class TypeTag>
767class TemperatureModel<TypeTag, false>
768{
772
773public:
775 { }
776
781 template <class Restarter>
782 void serialize(Restarter&)
783 { /* not implemented */ }
784
791 template <class Restarter>
792 void deserialize(Restarter&)
793 { /* not implemented */ }
794
795 void init() {}
797 const Scalar temperature(size_t /*globalIdx*/) const
798 {
799 return 273.15; // return 0C to make the compiler happy
800 }
801};
802
803} // namespace Opm
804
805#endif // OPM_TEMPERATURE_MODEL_HPP
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
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