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 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;
82
83 public:
84 using EvaluationTemp = DenseAd::Evaluation<Scalar, 1>;
85 using FluidStateTemp = BlackOilFluidState<EvaluationTemp,
86 FluidSystem,
87 true, //store temperature
88 true, // store enthalpy
89 compositionSwitchEnabled,
90 enableVapwat,
91 enableBrine,
92 enableSaltPrecipitation,
93 enableDisgasInWater,
94 enableSolvent,
95 Indices::numPhases>;
96
97
98
99 OPM_HOST_DEVICE void updateTemperature_(const Problem& problem,
100 unsigned globalDofIdx,
101 unsigned timeIdx)
102 {
103 const EvaluationTemp T = EvaluationTemp::createVariable(problem.temperature(globalDofIdx, timeIdx), 0);
104 fluidState_.setTemperature(T);
105 }
106
107 OPM_HOST_DEVICE void updateEnergyQuantities_(const Problem& problem,
108 const unsigned globalSpaceIdx,
109 const unsigned timeIdx)
110 {
111 // compute the specific enthalpy of the fluids, the specific enthalpy of the rock
112 // and the thermal conductivity coefficients
113 for (int phaseIdx = 0; phaseIdx < numPhases; ++ phaseIdx) {
114 if (!FluidSystem::phaseIsActive(phaseIdx)) {
115 continue;
116 }
117
118 const auto& h = FluidSystem::enthalpy(fluidState_, phaseIdx, problem.pvtRegionIndex(globalSpaceIdx));
119 fluidState_.setEnthalpy(phaseIdx, h);
120 }
121
122 const auto& solidEnergyLawParams = problem.solidEnergyLawParams(globalSpaceIdx, timeIdx);
123 rockInternalEnergy_ = SolidEnergyLaw::solidInternalEnergy(solidEnergyLawParams, fluidState_);
124
125 const auto& thermalConductionLawParams = problem.thermalConductionLawParams(globalSpaceIdx, timeIdx);
126 totalThermalConductivity_ = ThermalConductionLaw::thermalConductivity(thermalConductionLawParams, fluidState_);
127
128 // Retrieve the rock fraction from the problem
129 // Usually 1 - porosity, but if pvmult is used to modify porosity
130 // we will apply the same multiplier to the rock fraction
131 // i.e. pvmult*(1 - porosity) and thus interpret multpv as a volume
132 // multiplier. This is to avoid negative rock volume for pvmult*porosity > 1
133 rockFraction_ = problem.rockFraction(globalSpaceIdx, timeIdx);
134 }
135
136 OPM_HOST_DEVICE const EvaluationTemp& rockInternalEnergy() const
137 { return rockInternalEnergy_; }
138
139 OPM_HOST_DEVICE const EvaluationTemp& totalThermalConductivity() const
140 { return totalThermalConductivity_; }
141
142 OPM_HOST_DEVICE const Scalar& rockFraction() const
143 { return rockFraction_; }
144
145 OPM_HOST_DEVICE const FluidStateTemp& fluidStateTemp() const
146 { return fluidState_; }
147
148 template <class FluidState>
149 OPM_HOST_DEVICE void setFluidState(const FluidState& fs)
150 {
151 // copy the needed part of the fluid state
152 fluidState_.setPvtRegionIndex(fs.pvtRegionIndex());
153 fluidState_.setRs(getValue(fs.Rs()));
154 fluidState_.setRv(getValue(fs.Rv()));
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)));
160 }
161 }
162
163protected:
168};
169
175template <class TypeTag, bool enableTempV = getPropValue<TypeTag, Properties::EnergyModuleType>() == EnergyModules::SequentialImplicitThermal >
176class TemperatureModel : public GenericTemperatureModel<GetPropType<TypeTag, Properties::Grid>,
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>>
182{
200 using IndexTraits = typename FluidSystem::IndexTraitsType;
202
204 using FluidStateTemp = typename IntensiveQuantitiesTemp::FluidStateTemp;
205 using Evaluation = typename IntensiveQuantitiesTemp::EvaluationTemp;
206 using EnergyMatrix = typename BaseType::EnergyMatrix;
207 using EnergyVector = typename BaseType::EnergyVector;
208 using MatrixBlockTemp = typename BaseType::MatrixBlockTemp;
209 using SpareMatrixEnergyAdapter = typename Linear::IstlSparseMatrixAdapter<MatrixBlockTemp>;
210
211
212 //using ResidualNBInfo = typename LocalResidual::ResidualNBInfo;
213 struct ResidualNBInfo
214 {
215 double faceArea;
216 double inAlpha;
217 double outAlpha;
218 };
219
221
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;
228
229public:
230 explicit TemperatureModel(Simulator& simulator)
231 : BaseType(simulator.vanguard().gridView(),
232 simulator.vanguard().eclState(),
233 simulator.vanguard().cartesianIndexMapper(),
234 simulator.model().dofMapper())
235 , simulator_(simulator)
236 {}
237
238 void init()
239 {
240 const unsigned int numCells = simulator_.model().numTotalDof();
241 this->doInit(numCells);
242
243 if (!this->doTemp())
244 return;
245
246 // we need the storage term at start of the iteration (timeIdx = 1)
247 storage1_.resize(numCells);
248
249 // set the initial temperature
250 for (unsigned globI = 0; globI < numCells; ++globI) {
251 this->temperature_[globI] = simulator_.problem().initialFluidState(globI).temperature(0);
252 }
253 // keep a copy of the intensive quantities to simplify the update during
254 // the newton iterations
255 intQuants_.resize(numCells);
256
257 // find and store the overlap cells
258 const auto& elemMapper = simulator_.model().elementMapper();
260
261 // set the scaling factor
262 scalingFactor_ = getPropValue<TypeTag, Properties::BlackOilEnergyScalingFactor>();
263
264 // find the sparsity pattern of the temperature matrix
265 using NeighborSet = std::set<unsigned>;
266 std::vector<NeighborSet> neighbors(numCells);
267 Stencil stencil(this->gridView_, this->dofMapper_);
268 neighborInfo_.reserve(numCells, 6 * 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); // Do not include the primary dof in neighborInfo_
275 for (unsigned dofIdx = 0; dofIdx < stencil.numDof(); ++dofIdx) {
276 const unsigned neighborIdx = stencil.globalSpaceIndex(dofIdx);
277 neighbors[myIdx].insert(neighborIdx);
278 if (dofIdx > 0) {
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};
286 }
287 }
288 neighborInfo_.appendRow(loc_nbinfo.begin(), loc_nbinfo.end());
289 }
290 }
291 // allocate matrix for storing the Jacobian of the temperature residual
292 energyMatrix_ = std::make_unique<SpareMatrixEnergyAdapter>(simulator_);
293 diagMatAddress_.resize(numCells);
294 energyMatrix_->reserve(neighbors);
295 for (unsigned globI = 0; globI < numCells; globI++) {
296 const auto& nbInfos = neighborInfo_[globI];
297 diagMatAddress_[globI] = energyMatrix_->blockAddress(globI, globI);
298 for (auto& nbInfo : nbInfos) {
299 nbInfo.matBlockAddress = energyMatrix_->blockAddress(nbInfo.neighbor, globI);
300 }
301 }
302 }
303
305 {
306 if (!this->doTemp()) {
307 return;
308 }
309
310 // We use the specialized intensive quantities here with only the temperature derivative
311 const unsigned int numCells = simulator_.model().numTotalDof();
312 #ifdef _OPENMP
313 #pragma omp parallel for
314 #endif
315 for (unsigned globI = 0; globI < numCells; ++globI) {
316 intQuants_[globI].setFluidState(simulator_.model().intensiveQuantities(globI, /*timeIdx*/ 0).fluidState());
317 intQuants_[globI].updateTemperature_(simulator_.problem(), globI, /*timeIdx*/ 0);
318 intQuants_[globI].updateEnergyQuantities_(simulator_.problem(), globI, /*timeIdx*/ 0);
319 }
321
322 const int nw = simulator_.problem().wellModel().wellState().numWells();
323 this->energy_rates_.resize(nw, 0.0);
324 }
325
329 void endTimeStep(WellStateType& wellState)
330 {
331 if (!this->doTemp()) {
332 return;
333 }
334 OPM_TIMEBLOCK(TemperatureModel_endTimeStep);
335
336 // We use the specialized intensive quantities here with only the temperature derivative
337 const unsigned int numCells = simulator_.model().numTotalDof();
338 #ifdef _OPENMP
339 #pragma omp parallel for
340 #endif
341 for (unsigned globI = 0; globI < numCells; ++globI) {
342 intQuants_[globI].setFluidState(simulator_.model().intensiveQuantities(globI, /*timeIdx*/ 0).fluidState());
343 intQuants_[globI].updateTemperature_(simulator_.problem(), globI, /*timeIdx*/ 0);
344 intQuants_[globI].updateEnergyQuantities_(simulator_.problem(), globI, /*timeIdx*/ 0);
345 }
347
348 // update energy_rates
349 const int nw = wellState.numWells();
350 for (auto wellID = 0*nw; wellID < nw; ++wellID) {
351 auto& ws = wellState.well(wellID);
352 ws.energy_rate = this->energy_rates_[wellID];
353 }
354
355 // Update well temperature
356 const auto& wellPtrs = simulator_.problem().wellModel().localNonshutWells();
357 for (const auto& wellPtr : wellPtrs) {
358 auto& ws = wellState.well(wellPtr->name());
359 this->computeWellTemperature(*wellPtr, ws);
360 }
361 }
362
367 template <class Restarter>
368 void serialize(Restarter&)
369 { /* not implemented */ }
370
377 template <class Restarter>
378 void deserialize(Restarter&)
379 { /* not implemented */ }
380
381protected:
383 {
384 // we need the storage term at start of the iteration (timeIdx = 1)
385 const unsigned int numCells = simulator_.model().numTotalDof();
386 #ifdef _OPENMP
387 #pragma omp parallel for
388 #endif
389 for (unsigned globI = 0; globI < numCells; ++globI) {
390 Scalar storage = 0.0;
391 computeStorageTerm(globI, storage);
392 storage1_[globI] = storage;
393 }
394 }
395
397 {
398 OPM_TIMEBLOCK(TemperatureModel_advanceTemperatureFields);
399 const int max_iter = 20;
400 const int min_iter = 1;
401 bool is_converged = false;
402 // solve using Newton
403 for (int iter = 0; iter < max_iter; ++iter) {
405 if (iter >= min_iter && converged(iter)) {
406 is_converged = true;
407 break;
408 }
410 }
411 if (!is_converged) {
412 const auto msg =
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."),
415 max_iter);
416 OpmLog::debug(msg);
417 }
418 }
419
421 {
422 OPM_TIMEBLOCK(TemperatureModel_solveAndUpdate);
423 const unsigned int numCells = simulator_.model().numTotalDof();
424 EnergyVector dx(numCells);
425 bool conv = this->linearSolve_(this->energyMatrix_->istlMatrix(), dx, this->energyVector_);
426 if (!conv) {
427 if (simulator_.gridView().comm().rank() == 0) {
428 OpmLog::warning("Temp model: Linear solver did not converge. Temperature values not updated.");
429 }
430 }
431 else {
432 OPM_TIMEBLOCK(TemperatureModel_solveAndUpdate_update);
433 #ifdef _OPENMP
434 #pragma omp parallel for
435 #endif
436 for (unsigned globI = 0; globI < numCells; ++globI) {
437 this->temperature_[globI] -= std::clamp(dx[globI][0], -this->maxTempChange_, this->maxTempChange_);
438 intQuants_[globI].updateTemperature_(simulator_.problem(), globI, /*timeIdx*/ 0);
439 intQuants_[globI].updateEnergyQuantities_(simulator_.problem(), globI, /*timeIdx*/ 0);
440 }
441 }
442 }
443
444 bool converged(const int iter)
445 {
446 OPM_TIMEBLOCK(TemperatureModel_converged);
447 Scalar dt = simulator_.timeStepSize();
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();
452 const IsNumericalAquiferCell isNumericalAquiferCell(simulator_.gridView().grid());
453 Scalar sum_pv = 0.0;
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, /*timeIdx=*/0)
458 * simulator_.model().dofTotalVolume(globI);
459
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;
466 }
467 sum_pv += pvValue;
468 }
469 }
470 {
471 OPM_TIMEBLOCK(TemperatureModel_converged_communicate);
472
473 maxNorm = simulator_.gridView().comm().max(maxNorm);
474 sumNorm = simulator_.gridView().comm().sum(sumNorm);
475 sum_pv = simulator_.gridView().comm().sum(sum_pv);
476 sumNorm /= sum_pv;
477
478 // Use relaxed tolerance if the fraction of unconverged cells porevolume is less than relaxed_max_pv_fraction
479 sum_pv_not_converged = simulator_.gridView().comm().sum(sum_pv_not_converged);
480 }
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;
487
488 const auto msg = fmt::format(fmt::runtime("Temperature model (TEMP): Newton iter {}: "
489 "CNV(E): {:.1e}, EB: {:.1e}"),
490 iter, maxNorm, sumNorm);
491 OpmLog::debug(msg);
492 if (maxNorm < tolerance_cnv_energy && sumNorm < tolerance_energy_balance) {
493 const auto msg2 =
494 fmt::format(fmt::runtime("Temperature model (TEMP): Newton converged after {} iterations"),
495 iter);
496 OpmLog::debug(msg2);
497 return true;
498 }
499 return false;
500 }
501
502 template< class LhsEval>
503 void computeStorageTerm(unsigned globI, LhsEval& storage)
504 {
505 const auto& intQuants = intQuants_[globI];
506 const auto& poro = getValue(simulator_.model().intensiveQuantities(globI, /*timeIdx*/ 0).porosity());
507 // accumulate the internal energy of the fluids
508 const auto& fs = intQuants.fluidStateTemp();
509 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++ phaseIdx) {
510 if (!FluidSystem::phaseIsActive(phaseIdx)) {
511 continue;
512 }
513
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));
517
518 storage += poro*S*u*rho;
519 }
520
521 // add the internal energy of the rock
522 const Scalar rockFraction = intQuants.rockFraction();
523 const auto& uRock = decay<LhsEval>(intQuants.rockInternalEnergy());
524 storage += rockFraction*uRock;
525 storage*= scalingFactor_;
526 }
527
528 template <class RateVector>
529 void computeFluxTerm(const FluidStateTemp& fsIn,
530 const FluidStateTemp& fsEx,
531 const RateVector& darcyFlux,
532 Evaluation& flux)
533 {
534 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++ phaseIdx) {
535 if (!FluidSystem::phaseIsActive(phaseIdx)) {
536 continue;
537 }
538
539 const unsigned activeCompIdx =
540 FluidSystem::canonicalToActiveCompIdx(FluidSystem::solventComponentIndex(phaseIdx));
541 bool inIsUp = darcyFlux[activeCompIdx] > 0;
542 const auto& fs = inIsUp ? fsIn : fsEx;
543 if (inIsUp) {
544 flux += fs.enthalpy(phaseIdx)
545 * fs.density(phaseIdx)
546 * getValue(darcyFlux[activeCompIdx]);
547 }
548 else {
549 flux += getValue(fs.enthalpy(phaseIdx))
550 * getValue(fs.density(phaseIdx))
551 * getValue(darcyFlux[activeCompIdx]);
552 }
553 }
554 flux *= scalingFactor_;
555 }
556
557 template < class ResidualNBInfo>
559 const IntensiveQuantitiesTemp& intQuantsEx,
560 const ResidualNBInfo& res_nbinfo,
561 Evaluation& heatFlux)
562 {
563 short interiorDofIdx = 0; // NB
564 short exteriorDofIdx = 1; // NB
565 EnergyModule::ExtensiveQuantities::updateEnergy(heatFlux,
566 interiorDofIdx, // focusDofIndex,
567 interiorDofIdx,
568 exteriorDofIdx,
569 intQuantsIn,
570 intQuantsEx,
571 intQuantsIn.fluidStateTemp(),
572 intQuantsEx.fluidStateTemp(),
573 res_nbinfo.inAlpha,
574 res_nbinfo.outAlpha,
575 res_nbinfo.faceArea);
576 heatFlux *= scalingFactor_*res_nbinfo.faceArea;
577 }
578
580 {
581 OPM_TIMEBLOCK(TemperatureModel_assembleEquations);
582
583 const unsigned int numCells = simulator_.model().numTotalDof();
584 for (unsigned globI = 0; globI < numCells; ++globI) {
585 this->energyVector_[globI] = 0.0;
586 energyMatrix_->clearRow(globI, 0.0);
587 }
588 const Scalar dt = simulator_.timeStepSize();
589
590 // Storage term
591 {
592 OPM_TIMEBLOCK(TemperatureModel_assembleEquations_storage);
593 #ifdef _OPENMP
594 #pragma omp parallel for
595 #endif
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;
601 computeStorageTerm(globI, storage);
602 this->energyVector_[globI] += storefac * ( getValue(storage) - storage1_[globI] );
603 bMat[0][0] = storefac * storage.derivative(temperatureIdx);
604 *diagMatAddress_[globI] += bMat;
605 }
606 }
607
608 // Flux term
609 {
610 OPM_TIMEBLOCK(TemperatureModel_assembleEquations_flux);
611 const auto& floresInfo = this->simulator_.problem().model().linearizer().getFloresInfo();
612 const bool enableDriftCompensation = Parameters::Get<Parameters::EnableDriftCompensationTemp>();
613 const auto& problem = simulator_.problem();
614
615 #ifdef _OPENMP
616 #pragma omp parallel for
617 #endif
618 for (unsigned globI = 0; globI < numCells; ++globI) {
619 const auto& nbInfos = neighborInfo_[globI];
620 const auto& floresInfos = floresInfo[globI];
621 int loc = 0;
622 const auto& intQuantsIn = intQuants_[globI];
623 MatrixBlockTemp bMat;
624 for (const auto& nbInfo : nbInfos) {
625 unsigned globJ = nbInfo.neighbor;
626 const auto& intQuantsEx = intQuants_[globJ];
627 assert(globJ != globI);
628 const auto& darcyflux = floresInfos[loc].flow;
629 // compute convective flux
630 Evaluation flux = 0.0;
631 computeFluxTerm(intQuantsIn.fluidStateTemp(), intQuantsEx.fluidStateTemp(), darcyflux, flux);
632 // compute conductive flux
633 Evaluation heatFlux = 0.0;
634 computeHeatFluxTerm(intQuantsIn, intQuantsEx, nbInfo.res_nbinfo, heatFlux);
635 heatFlux += flux;
636 this->energyVector_[globI] += getValue(heatFlux);
637 bMat[0][0] = heatFlux.derivative(temperatureIdx);
638 *diagMatAddress_[globI] += bMat;
639 bMat *= -1.0;
640 //SparseAdapter syntax: jacobian_->addToBlock(globJ, globI, bMat);
641 *nbInfo.matBlockAddress += bMat;
642 loc++;
643 }
644
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));
652 this->energyVector_[globI] -= drift_hrate*scalingFactor_;
653 }
654 }
655 }
656 }
657
658 // Well terms
659 {
660 OPM_TIMEBLOCK(TemperatureModel_assembleEquations_wells);
661 const auto& wellPtrs = simulator_.problem().wellModel().localNonshutWells();
662 for (const auto& wellPtr : wellPtrs) {
663 this->assembleEquationWell(*wellPtr);
664 }
665 }
666
667 // For standard ilu0 we need to set the overlapping cells to identity. But since we use "paroverilu0"
668 // here this is not needed.
669 // if (simulator_.gridView().comm().size() > 1) {
670 // // Set dirichlet conditions for overlapping cells
671 // // loop over precalculated overlap rows and columns
672 // for (const auto row : overlapRows_) {
673 // // Zero out row.
674 // (*this->energyMatrix_)[row] = 0.0;
675 //
676 // //diagonal block set to diag(1.0).
677 // (*this->energyMatrix_)[row][row][0][0] = 1.0;
678 // }
679 //}
680 }
681
682 template<class Well>
683 void assembleEquationWell(const Well& well)
684 {
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);
688 this->energy_rates_[well_index] = 0.0;
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];
692 auto fs = intQuants_[globI].fluidStateTemp(); //copy to make it possible to change the temp in the injector
693 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++ phaseIdx) {
694 if (!FluidSystem::phaseIsActive(phaseIdx)) {
695 continue;
696 }
697
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));
706 } else {
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());
711 rate /= d;
712 }
713 if (phaseIdx == oilPhaseIdx && d > 0) {
714 const auto& gasrate = well.volumetricSurfaceRateForConnection(globI, gasPhaseIdx);
715 rate -= gasrate * getValue(fs.Rv());
716 rate /= d;
717 }
718 rate *= fs.enthalpy(phaseIdx) * getValue(fs.density(phaseIdx)) / getValue(fs.invB(phaseIdx));
719 }
720 this->energy_rates_[well_index] += getValue(rate);
721 rate *= scalingFactor_;
722 this->energyVector_[globI] -= getValue(rate);
723 bMat[0][0] = -rate.derivative(temperatureIdx);
724 *diagMatAddress_[globI] += bMat;
725 }
726 }
727 }
728
729 template<class Well, class SingleWellState>
730 void computeWellTemperature(const Well& well, SingleWellState& ws)
731 {
732 if (well.isInjector()) {
733 if (ws.status != WellStatus::STOP) {
734 ws.temperature = well.wellEcl().inj_temperature();
735 return;
736 }
737 }
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(/*phaseIdx*/0).value();
747 }
748 //simulator_.gridView().comm().sum(weighted.data(), 2);
749 ws.temperature = weighted_temperature / total_weight;
750 }
751
752 const Simulator& simulator_;
753 EnergyVector storage1_;
754 std::vector<IntensiveQuantitiesTemp> intQuants_;
756 std::vector<MatrixBlockTemp*> diagMatAddress_{};
757 std::unique_ptr<SpareMatrixEnergyAdapter> energyMatrix_;
758 std::vector<int> overlapRows_;
759 std::vector<int> interiorRows_;
760 Scalar scalingFactor_{1.0};
761};
762
763// need for the old linearizer
764template <class TypeTag>
765class TemperatureModel<TypeTag, false>
766{
770
771public:
773 { }
774
779 template <class Restarter>
780 void serialize(Restarter&)
781 { /* not implemented */ }
782
789 template <class Restarter>
790 void deserialize(Restarter&)
791 { /* not implemented */ }
792
793 void init() {}
795 const Scalar temperature(size_t /*globalIdx*/) const
796 {
797 return 273.15; // return 0C to make the compiler happy
798 }
799};
800
801} // namespace Opm
802
803#endif // OPM_TEMPERATURE_MODEL_HPP
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
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