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
35
42
43#include <algorithm>
44#include <cassert>
45#include <cmath>
46#include <cstddef>
47#include <string>
48#include <vector>
49
50namespace Opm::Properties {
51
52template<class TypeTag, class MyTypeTag>
55};
56
57} // namespace Opm::Properties
58
59namespace Opm {
60
61template<typename Scalar, typename IndexTraits> class WellState;
62
63
64template <class TypeTag>
66{
73
74
75
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;
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 Indices::numPhases>;
95
96
97
98 void updateTemperature_(const Problem& problem, unsigned globalDofIdx, unsigned timeIdx)
99 {
100 const EvaluationTemp T = EvaluationTemp::createVariable(problem.temperature(globalDofIdx, timeIdx), 0);
101 fluidState_.setTemperature(T);
102 }
103
104 void updateEnergyQuantities_(const Problem& problem,
105 const unsigned globalSpaceIdx,
106 const unsigned timeIdx)
107 {
108 // compute the specific enthalpy of the fluids, the specific enthalpy of the rock
109 // and the thermal conductivity coefficients
110 for (int phaseIdx = 0; phaseIdx < numPhases; ++ phaseIdx) {
111 if (!FluidSystem::phaseIsActive(phaseIdx)) {
112 continue;
113 }
114
115 const auto& h = FluidSystem::enthalpy(fluidState_, phaseIdx, problem.pvtRegionIndex(globalSpaceIdx));
116 fluidState_.setEnthalpy(phaseIdx, h);
117 }
118
119 const auto& solidEnergyLawParams = problem.solidEnergyLawParams(globalSpaceIdx, timeIdx);
120 rockInternalEnergy_ = SolidEnergyLaw::solidInternalEnergy(solidEnergyLawParams, fluidState_);
121
122 const auto& thermalConductionLawParams = problem.thermalConductionLawParams(globalSpaceIdx, timeIdx);
123 totalThermalConductivity_ = ThermalConductionLaw::thermalConductivity(thermalConductionLawParams, fluidState_);
124
125 // Retrieve the rock fraction from the problem
126 // Usually 1 - porosity, but if pvmult is used to modify porosity
127 // we will apply the same multiplier to the rock fraction
128 // i.e. pvmult*(1 - porosity) and thus interpret multpv as a volume
129 // multiplier. This is to avoid negative rock volume for pvmult*porosity > 1
130 rockFraction_ = problem.rockFraction(globalSpaceIdx, timeIdx);
131 }
132
134 { return rockInternalEnergy_; }
135
137 { return totalThermalConductivity_; }
138
139 const Scalar& rockFraction() const
140 { return rockFraction_; }
141
143 { return fluidState_; }
144
145 template <class FluidState>
146 void setFluidState(const FluidState& fs)
147 {
148 // copy the needed part of the fluid state
149 fluidState_.setPvtRegionIndex(fs.pvtRegionIndex());
150 fluidState_.setRs(getValue(fs.Rs()));
151 fluidState_.setRv(getValue(fs.Rv()));
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)));
157 }
158 }
159
160protected:
165};
166
172template <class TypeTag, bool enableTempV = getPropValue<TypeTag, Properties::EnergyModuleType>() == EnergyModules::SequentialImplicitThermal >
173class TemperatureModel : public GenericTemperatureModel<GetPropType<TypeTag, Properties::Grid>,
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>>
179{
197 using IndexTraits = typename FluidSystem::IndexTraitsType;
199
201 using FluidStateTemp = typename IntensiveQuantitiesTemp::FluidStateTemp;
202 using Evaluation = typename IntensiveQuantitiesTemp::EvaluationTemp;
203 using EnergyMatrix = typename BaseType::EnergyMatrix;
204 using EnergyVector = typename BaseType::EnergyVector;
205 using MatrixBlockTemp = typename BaseType::MatrixBlockTemp;
206 using SpareMatrixEnergyAdapter = typename Linear::IstlSparseMatrixAdapter<MatrixBlockTemp>;
207
208
209 //using ResidualNBInfo = typename LocalResidual::ResidualNBInfo;
210 struct ResidualNBInfo
211 {
212 double faceArea;
213 double inAlpha;
214 double outAlpha;
215 };
216
218
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;
225
226public:
227 explicit TemperatureModel(Simulator& simulator)
228 : BaseType(simulator.vanguard().gridView(),
229 simulator.vanguard().eclState(),
230 simulator.vanguard().cartesianIndexMapper(),
231 simulator.model().dofMapper())
232 , simulator_(simulator)
233 {}
234
235 void init()
236 {
237 const unsigned int numCells = simulator_.model().numTotalDof();
238 this->doInit(numCells);
239
240 if (!this->doTemp())
241 return;
242
243 // we need the storage term at start of the iteration (timeIdx = 1)
244 storage1_.resize(numCells);
245
246 // set the initial temperature
247 for (unsigned globI = 0; globI < numCells; ++globI) {
248 this->temperature_[globI] = simulator_.problem().initialFluidState(globI).temperature(0);
249 }
250 // keep a copy of the intensive quantities to simplify the update during
251 // the newton iterations
252 intQuants_.resize(numCells);
253
254 // find and store the overlap cells
255 const auto& elemMapper = simulator_.model().elementMapper();
257
258 // set the scaling factor
259 scalingFactor_ = getPropValue<TypeTag, Properties::BlackOilEnergyScalingFactor>();
260
261 // find the sparsity pattern of the temperature matrix
262 using NeighborSet = std::set<unsigned>;
263 std::vector<NeighborSet> neighbors(numCells);
264 Stencil stencil(this->gridView_, this->dofMapper_);
265 neighborInfo_.reserve(numCells, 6 * 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); // Do not include the primary dof in neighborInfo_
272 for (unsigned dofIdx = 0; dofIdx < stencil.numDof(); ++dofIdx) {
273 const unsigned neighborIdx = stencil.globalSpaceIndex(dofIdx);
274 neighbors[myIdx].insert(neighborIdx);
275 if (dofIdx > 0) {
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};
283 }
284 }
285 neighborInfo_.appendRow(loc_nbinfo.begin(), loc_nbinfo.end());
286 }
287 }
288 // allocate matrix for storing the Jacobian of the temperature residual
289 energyMatrix_ = std::make_unique<SpareMatrixEnergyAdapter>(simulator_);
290 diagMatAddress_.resize(numCells);
291 energyMatrix_->reserve(neighbors);
292 for (unsigned globI = 0; globI < numCells; globI++) {
293 const auto& nbInfos = neighborInfo_[globI];
294 diagMatAddress_[globI] = energyMatrix_->blockAddress(globI, globI);
295 for (auto& nbInfo : nbInfos) {
296 nbInfo.matBlockAddress = energyMatrix_->blockAddress(nbInfo.neighbor, globI);
297 }
298 }
299 }
300
302 {
303 if (!this->doTemp()) {
304 return;
305 }
306
307 // We use the specialized intensive quantities here with only the temperature derivative
308 const unsigned int numCells = simulator_.model().numTotalDof();
309 #ifdef _OPENMP
310 #pragma omp parallel for
311 #endif
312 for (unsigned globI = 0; globI < numCells; ++globI) {
313 intQuants_[globI].setFluidState(simulator_.model().intensiveQuantities(globI, /*timeIdx*/ 0).fluidState());
314 intQuants_[globI].updateTemperature_(simulator_.problem(), globI, /*timeIdx*/ 0);
315 intQuants_[globI].updateEnergyQuantities_(simulator_.problem(), globI, /*timeIdx*/ 0);
316 }
318
319 const int nw = simulator_.problem().wellModel().wellState().numWells();
320 this->energy_rates_.resize(nw, 0.0);
321 }
322
326 void endTimeStep(WellStateType& wellState)
327 {
328 if (!this->doTemp()) {
329 return;
330 }
331
332 // We use the specialized intensive quantities here with only the temperature derivative
333 const unsigned int numCells = simulator_.model().numTotalDof();
334 #ifdef _OPENMP
335 #pragma omp parallel for
336 #endif
337 for (unsigned globI = 0; globI < numCells; ++globI) {
338 intQuants_[globI].setFluidState(simulator_.model().intensiveQuantities(globI, /*timeIdx*/ 0).fluidState());
339 intQuants_[globI].updateTemperature_(simulator_.problem(), globI, /*timeIdx*/ 0);
340 intQuants_[globI].updateEnergyQuantities_(simulator_.problem(), globI, /*timeIdx*/ 0);
341 }
343
344 // update energy_rates
345 const int nw = wellState.numWells();
346 for (auto wellID = 0*nw; wellID < nw; ++wellID) {
347 auto& ws = wellState.well(wellID);
348 ws.energy_rate = this->energy_rates_[wellID];
349 }
350
351 // Update well temperature
352 const auto& wellPtrs = simulator_.problem().wellModel().localNonshutWells();
353 for (const auto& wellPtr : wellPtrs) {
354 auto& ws = wellState.well(wellPtr->name());
355 this->computeWellTemperature(*wellPtr, ws);
356 }
357 }
358
363 template <class Restarter>
364 void serialize(Restarter&)
365 { /* not implemented */ }
366
373 template <class Restarter>
374 void deserialize(Restarter&)
375 { /* not implemented */ }
376
377protected:
379 {
380 // we need the storage term at start of the iteration (timeIdx = 1)
381 const unsigned int numCells = simulator_.model().numTotalDof();
382 #ifdef _OPENMP
383 #pragma omp parallel for
384 #endif
385 for (unsigned globI = 0; globI < numCells; ++globI) {
386 Scalar storage = 0.0;
387 computeStorageTerm(globI, storage);
388 storage1_[globI] = storage;
389 }
390 }
391
393 {
394 const int max_iter = 20;
395 const int min_iter = 1;
396 bool is_converged = false;
397 // solve using Newton
398 for (int iter = 0; iter < max_iter; ++iter) {
400 if (iter >= min_iter && converged(iter)) {
401 is_converged = true;
402 break;
403 }
405 }
406 if (!is_converged) {
407 const auto msg =
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."),
410 max_iter);
411 OpmLog::debug(msg);
412 }
413 }
414
416 {
417 const unsigned int numCells = simulator_.model().numTotalDof();
418 EnergyVector dx(numCells);
419 bool conv = this->linearSolve_(this->energyMatrix_->istlMatrix(), dx, this->energyVector_);
420 if (!conv) {
421 if (simulator_.gridView().comm().rank() == 0) {
422 OpmLog::warning("Temp model: Linear solver did not converge. Temperature values not updated.");
423 }
424 }
425 else {
426 #ifdef _OPENMP
427 #pragma omp parallel for
428 #endif
429 for (unsigned globI = 0; globI < numCells; ++globI) {
430 this->temperature_[globI] -= std::clamp(dx[globI][0], -this->maxTempChange_, this->maxTempChange_);
431 intQuants_[globI].updateTemperature_(simulator_.problem(), globI, /*timeIdx*/ 0);
432 intQuants_[globI].updateEnergyQuantities_(simulator_.problem(), globI, /*timeIdx*/ 0);
433 }
434 }
435 }
436
437 bool converged(const int iter)
438 {
439 Scalar dt = simulator_.timeStepSize();
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();
444 const IsNumericalAquiferCell isNumericalAquiferCell(simulator_.gridView().grid());
445 Scalar sum_pv = 0.0;
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, /*timeIdx=*/0)
450 * simulator_.model().dofTotalVolume(globI);
451
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;
458 }
459 sum_pv += pvValue;
460 }
461 }
462 maxNorm = simulator_.gridView().comm().max(maxNorm);
463 sumNorm = simulator_.gridView().comm().sum(sumNorm);
464 sum_pv = simulator_.gridView().comm().sum(sum_pv);
465 sumNorm /= sum_pv;
466
467 // Use relaxed tolerance if the fraction of unconverged cells porevolume is less than relaxed_max_pv_fraction
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;
475
476 const auto msg = fmt::format(fmt::runtime("Temperature model (TEMP): Newton iter {}: "
477 "CNV(E): {:.1e}, EB: {:.1e}"),
478 iter, maxNorm, sumNorm);
479 OpmLog::debug(msg);
480 if (maxNorm < tolerance_cnv_energy && sumNorm < tolerance_energy_balance) {
481 const auto msg2 =
482 fmt::format(fmt::runtime("Temperature model (TEMP): Newton converged after {} iterations"),
483 iter);
484 OpmLog::debug(msg2);
485 return true;
486 }
487 return false;
488 }
489
490 template< class LhsEval>
491 void computeStorageTerm(unsigned globI, LhsEval& storage)
492 {
493 const auto& intQuants = intQuants_[globI];
494 const auto& poro = getValue(simulator_.model().intensiveQuantities(globI, /*timeIdx*/ 0).porosity());
495 // accumulate the internal energy of the fluids
496 const auto& fs = intQuants.fluidStateTemp();
497 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++ phaseIdx) {
498 if (!FluidSystem::phaseIsActive(phaseIdx)) {
499 continue;
500 }
501
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));
505
506 storage += poro*S*u*rho;
507 }
508
509 // add the internal energy of the rock
510 const Scalar rockFraction = intQuants.rockFraction();
511 const auto& uRock = decay<LhsEval>(intQuants.rockInternalEnergy());
512 storage += rockFraction*uRock;
513 storage*= scalingFactor_;
514 }
515
516 template <class RateVector>
517 void computeFluxTerm(const FluidStateTemp& fsIn,
518 const FluidStateTemp& fsEx,
519 const RateVector& darcyFlux,
520 Evaluation& flux)
521 {
522 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++ phaseIdx) {
523 if (!FluidSystem::phaseIsActive(phaseIdx)) {
524 continue;
525 }
526
527 const unsigned activeCompIdx =
528 FluidSystem::canonicalToActiveCompIdx(FluidSystem::solventComponentIndex(phaseIdx));
529 bool inIsUp = darcyFlux[activeCompIdx] > 0;
530 const auto& fs = inIsUp ? fsIn : fsEx;
531 if (inIsUp) {
532 flux += fs.enthalpy(phaseIdx)
533 * fs.density(phaseIdx)
534 * getValue(darcyFlux[activeCompIdx]);
535 }
536 else {
537 flux += getValue(fs.enthalpy(phaseIdx))
538 * getValue(fs.density(phaseIdx))
539 * getValue(darcyFlux[activeCompIdx]);
540 }
541 }
542 flux *= scalingFactor_;
543 }
544
545 template < class ResidualNBInfo>
547 const IntensiveQuantitiesTemp& intQuantsEx,
548 const ResidualNBInfo& res_nbinfo,
549 Evaluation& heatFlux)
550 {
551 short interiorDofIdx = 0; // NB
552 short exteriorDofIdx = 1; // NB
553 EnergyModule::ExtensiveQuantities::updateEnergy(heatFlux,
554 interiorDofIdx, // focusDofIndex,
555 interiorDofIdx,
556 exteriorDofIdx,
557 intQuantsIn,
558 intQuantsEx,
559 intQuantsIn.fluidStateTemp(),
560 intQuantsEx.fluidStateTemp(),
561 res_nbinfo.inAlpha,
562 res_nbinfo.outAlpha,
563 res_nbinfo.faceArea);
564 heatFlux *= scalingFactor_*res_nbinfo.faceArea;
565 }
566
568 {
569 const unsigned int numCells = simulator_.model().numTotalDof();
570 for (unsigned globI = 0; globI < numCells; ++globI) {
571 this->energyVector_[globI] = 0.0;
572 energyMatrix_->clearRow(globI, 0.0);
573 }
574 Scalar dt = simulator_.timeStepSize();
575 #ifdef _OPENMP
576 #pragma omp parallel for
577 #endif
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;
583 computeStorageTerm(globI, storage);
584 this->energyVector_[globI] += storefac * ( getValue(storage) - storage1_[globI] );
585 bMat[0][0] = storefac * storage.derivative(temperatureIdx);
586 *diagMatAddress_[globI] += bMat;
587 }
588
589 const auto& floresInfo = this->simulator_.problem().model().linearizer().getFloresInfo();
590 const bool enableDriftCompensation = Parameters::Get<Parameters::EnableDriftCompensationTemp>();
591 const auto& problem = simulator_.problem();
592
593 #ifdef _OPENMP
594 #pragma omp parallel for
595 #endif
596 for (unsigned globI = 0; globI < numCells; ++globI) {
597 const auto& nbInfos = neighborInfo_[globI];
598 const auto& floresInfos = floresInfo[globI];
599 int loc = 0;
600 const auto& intQuantsIn = intQuants_[globI];
601 MatrixBlockTemp bMat;
602 for (const auto& nbInfo : nbInfos) {
603 unsigned globJ = nbInfo.neighbor;
604 const auto& intQuantsEx = intQuants_[globJ];
605 assert(globJ != globI);
606 const auto& darcyflux = floresInfos[loc].flow;
607 // compute convective flux
608 Evaluation flux = 0.0;
609 computeFluxTerm(intQuantsIn.fluidStateTemp(), intQuantsEx.fluidStateTemp(), darcyflux, flux);
610 // compute conductive flux
611 Evaluation heatFlux = 0.0;
612 computeHeatFluxTerm(intQuantsIn, intQuantsEx, nbInfo.res_nbinfo, heatFlux);
613 heatFlux += flux;
614 this->energyVector_[globI] += getValue(heatFlux);
615 bMat[0][0] = heatFlux.derivative(temperatureIdx);
616 *diagMatAddress_[globI] += bMat;
617 bMat *= -1.0;
618 //SparseAdapter syntax: jacobian_->addToBlock(globJ, globI, bMat);
619 *nbInfo.matBlockAddress += bMat;
620 loc++;
621 }
622
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));
630 this->energyVector_[globI] -= drift_hrate*scalingFactor_;
631 }
632 }
633 }
634
635 // Well terms
636 const auto& wellPtrs = simulator_.problem().wellModel().localNonshutWells();
637 for (const auto& wellPtr : wellPtrs) {
638 this->assembleEquationWell(*wellPtr);
639 }
640
641 // For standard ilu0 we need to set the overlapping cells to identity. But since we use "paroverilu0"
642 // here this is not needed.
643 // if (simulator_.gridView().comm().size() > 1) {
644 // // Set dirichlet conditions for overlapping cells
645 // // loop over precalculated overlap rows and columns
646 // for (const auto row : overlapRows_) {
647 // // Zero out row.
648 // (*this->energyMatrix_)[row] = 0.0;
649 //
650 // //diagonal block set to diag(1.0).
651 // (*this->energyMatrix_)[row][row][0][0] = 1.0;
652 // }
653 //}
654 }
655
656 template<class Well>
657 void assembleEquationWell(const Well& well)
658 {
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);
662 this->energy_rates_[well_index] = 0.0;
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];
666 auto fs = intQuants_[globI].fluidStateTemp(); //copy to make it possible to change the temp in the injector
667 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++ phaseIdx) {
668 if (!FluidSystem::phaseIsActive(phaseIdx)) {
669 continue;
670 }
671
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));
680 } else {
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());
685 rate /= d;
686 }
687 if (phaseIdx == oilPhaseIdx && d > 0) {
688 const auto& gasrate = well.volumetricSurfaceRateForConnection(globI, gasPhaseIdx);
689 rate -= gasrate * getValue(fs.Rv());
690 rate /= d;
691 }
692 rate *= fs.enthalpy(phaseIdx) * getValue(fs.density(phaseIdx)) / getValue(fs.invB(phaseIdx));
693 }
694 this->energy_rates_[well_index] += getValue(rate);
695 rate *= scalingFactor_;
696 this->energyVector_[globI] -= getValue(rate);
697 bMat[0][0] = -rate.derivative(temperatureIdx);
698 *diagMatAddress_[globI] += bMat;
699 }
700 }
701 }
702
703 template<class Well, class SingleWellState>
704 void computeWellTemperature(const Well& well, SingleWellState& ws)
705 {
706 if (well.isInjector()) {
707 if (ws.status != WellStatus::STOP) {
708 ws.temperature = well.wellEcl().inj_temperature();
709 return;
710 }
711 }
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(/*phaseIdx*/0).value();
721 }
722 //simulator_.gridView().comm().sum(weighted.data(), 2);
723 ws.temperature = weighted_temperature / total_weight;
724 }
725
726 const Simulator& simulator_;
727 EnergyVector storage1_;
728 std::vector<IntensiveQuantitiesTemp> intQuants_;
730 std::vector<MatrixBlockTemp*> diagMatAddress_{};
731 std::unique_ptr<SpareMatrixEnergyAdapter> energyMatrix_;
732 std::vector<int> overlapRows_;
733 std::vector<int> interiorRows_;
734 Scalar scalingFactor_{1.0};
735};
736
737// need for the old linearizer
738template <class TypeTag>
739class TemperatureModel<TypeTag, false>
740{
744
745public:
747 { }
748
753 template <class Restarter>
754 void serialize(Restarter&)
755 { /* not implemented */ }
756
763 template <class Restarter>
764 void deserialize(Restarter&)
765 { /* not implemented */ }
766
767 void init() {}
769 const Scalar temperature(size_t /*globalIdx*/) const
770 {
771 return 273.15; // return 0C to make the compiler happy
772 }
773};
774
775} // namespace Opm
776
777#endif // OPM_TEMPERATURE_MODEL_HPP
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
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