28#ifndef OPM_TRACER_MODEL_HPP
29#define OPM_TRACER_MODEL_HPP
31#include <opm/common/OpmLog/OpmLog.hpp>
32#include <opm/common/TimingMacros.hpp>
34#include <opm/input/eclipse/Schedule/Well/Well.hpp>
35#include <opm/input/eclipse/Schedule/Well/WellConnections.hpp>
37#include <opm/grid/utility/ElementChunks.hpp>
55template<
class TypeTag,
class MyTypeTag>
69template <
class TypeTag>
71 GetPropType<TypeTag, Properties::GridView>,
72 GetPropType<TypeTag, Properties::DofMapper>,
73 GetPropType<TypeTag, Properties::Stencil>,
74 GetPropType<TypeTag, Properties::FluidSystem>,
75 GetPropType<TypeTag, Properties::Scalar>>
93 using TracerEvaluation = DenseAd::Evaluation<Scalar,1>;
99 enum { numEq = getPropValue<TypeTag, Properties::NumEq>() };
100 enum { numPhases = FluidSystem::numPhases };
101 enum { waterPhaseIdx = FluidSystem::waterPhaseIdx };
102 enum { oilPhaseIdx = FluidSystem::oilPhaseIdx };
103 enum { gasPhaseIdx = FluidSystem::gasPhaseIdx };
107 :
BaseType(simulator.vanguard().gridView(),
108 simulator.vanguard().eclState(),
109 simulator.vanguard().cartesianIndexMapper(),
110 simulator.model().dofMapper(),
111 simulator.vanguard().cellCentroids())
113 ,
tbatch({waterPhaseIdx, oilPhaseIdx, gasPhaseIdx})
142 gasPhaseIdx, oilPhaseIdx, waterPhaseIdx);
147 for (std::size_t tracerIdx = 0; tracerIdx < this->
tracerPhaseIdx_.size(); ++tracerIdx) {
149 if (! FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx)){
150 throw std::runtime_error(
"Water tracer specified for non-water fluid system: " +
151 this->
name(tracerIdx));
156 else if (this->
tracerPhaseIdx_[tracerIdx] == FluidSystem::oilPhaseIdx) {
157 if (! FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx)){
158 throw std::runtime_error(
"Oil tracer specified for non-oil fluid system: " +
159 this->
name(tracerIdx));
164 else if (this->
tracerPhaseIdx_[tracerIdx] == FluidSystem::gasPhaseIdx) {
165 if (! FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx)){
166 throw std::runtime_error(
"Gas tracer specified for non-gas fluid system: " +
167 this->
name(tracerIdx));
186 for (
auto& tr : this->
tbatch) {
187 if (tr.numTracer() != 0) {
192 tr.mat = std::make_unique<TracerMatrix>(*base);
204 OPM_TIMEBLOCK(tracerUpdateCache);
217 OPM_TIMEBLOCK(tracerAdvance);
225 template <
class Restarter>
235 template <
class Restarter>
239 template<
class Serializer>
242 serializer(
static_cast<BaseType&
>(*
this));
252 template<TracerTypeIdx Index>
254 const unsigned globalDofIdx,
255 const unsigned timeIdx)
const
257 const auto& intQuants =
simulator_.model().intensiveQuantities(globalDofIdx, timeIdx);
258 const auto& fs = intQuants.fluidState();
259 constexpr Scalar min_volume = 1e-10;
262 return std::max(decay<Scalar>(fs.saturation(tracerPhaseIdx)) *
263 decay<Scalar>(fs.invB(tracerPhaseIdx)) *
264 decay<Scalar>(intQuants.porosity()),
268 if (tracerPhaseIdx == FluidSystem::oilPhaseIdx && FluidSystem::enableVaporizedOil()) {
269 return std::max(decay<Scalar>(fs.saturation(FluidSystem::gasPhaseIdx)) *
270 decay<Scalar>(fs.invB(FluidSystem::gasPhaseIdx)) *
271 decay<Scalar>(fs.Rv()) *
272 decay<Scalar>(intQuants.porosity()),
277 else if (tracerPhaseIdx == FluidSystem::gasPhaseIdx && FluidSystem::enableDissolvedGas()) {
278 return std::max(decay<Scalar>(fs.saturation(FluidSystem::oilPhaseIdx)) *
279 decay<Scalar>(fs.invB(FluidSystem::oilPhaseIdx)) *
280 decay<Scalar>(fs.Rs()) *
281 decay<Scalar>(intQuants.porosity()),
289 template<TracerTypeIdx Index>
290 std::pair<TracerEvaluation, bool>
292 const ElementContext& elemCtx,
293 const unsigned scvfIdx,
294 const unsigned timeIdx)
const
296 const auto& stencil = elemCtx.stencil(timeIdx);
297 const auto& scvf = stencil.interiorFace(scvfIdx);
299 const auto& extQuants = elemCtx.extensiveQuantities(scvfIdx, timeIdx);
300 const unsigned inIdx = extQuants.interiorIndex();
306 upIdx = extQuants.upstreamIndex(tracerPhaseIdx);
307 const auto& intQuants = elemCtx.intensiveQuantities(upIdx, timeIdx);
308 const auto& fs = intQuants.fluidState();
309 v = decay<Scalar>(extQuants.volumeFlux(tracerPhaseIdx)) *
310 decay<Scalar>(fs.invB(tracerPhaseIdx));
312 if (tracerPhaseIdx == FluidSystem::oilPhaseIdx && FluidSystem::enableVaporizedOil()) {
313 upIdx = extQuants.upstreamIndex(FluidSystem::gasPhaseIdx);
315 const auto& intQuants = elemCtx.intensiveQuantities(upIdx, timeIdx);
316 const auto& fs = intQuants.fluidState();
317 v = decay<Scalar>(fs.invB(FluidSystem::gasPhaseIdx)) *
318 decay<Scalar>(extQuants.volumeFlux(FluidSystem::gasPhaseIdx)) *
319 decay<Scalar>(fs.Rv());
322 else if (tracerPhaseIdx == FluidSystem::gasPhaseIdx && FluidSystem::enableDissolvedGas()) {
323 upIdx = extQuants.upstreamIndex(FluidSystem::oilPhaseIdx);
325 const auto& intQuants = elemCtx.intensiveQuantities(upIdx, timeIdx);
326 const auto& fs = intQuants.fluidState();
327 v = decay<Scalar>(fs.invB(FluidSystem::oilPhaseIdx)) *
328 decay<Scalar>(extQuants.volumeFlux(FluidSystem::oilPhaseIdx)) *
329 decay<Scalar>(fs.Rs());
337 const Scalar A = scvf.area();
338 return inIdx == upIdx
339 ? std::pair{A * v * variable<TracerEvaluation>(1.0, 0),
true}
340 : std::pair{A * v,
false};
343 template<TracerTypeIdx Index,
class TrRe>
351 return tr.storageOfTimeIndex1_[tIdx][I][
Index];
353 return computeVolume_<Index>(tr.phaseIdx_, I1, 1) *
354 tr.concentration_[tIdx][I1][
Index];
360 const ElementContext& elemCtx,
361 const Scalar scvVolume,
367 if (tr.numTracer() == 0) {
371 const TracerEvaluation fVol = computeVolume_<Free>(tr.phaseIdx_, I, 0) * variable<TracerEvaluation>(1.0, 0);
372 const TracerEvaluation sVol = computeVolume_<Solution>(tr.phaseIdx_, I, 0) * variable<TracerEvaluation>(1.0, 0);
373 dVol_[
Solution][tr.phaseIdx_][I] += sVol.value() * scvVolume -
vol1_[1][tr.phaseIdx_][I];
374 dVol_[
Free][tr.phaseIdx_][I] += fVol.value() * scvVolume -
vol1_[0][tr.phaseIdx_][I];
375 for (
int tIdx = 0; tIdx < tr.numTracer(); ++tIdx) {
377 const Scalar fStorageOfTimeIndex0 = fVol.value() * tr.concentration_[tIdx][I][
Free];
378 const Scalar fLocalStorage = (fStorageOfTimeIndex0 - storage1_<Free>(tr, tIdx, I, I1,
379 elemCtx.enableStorageCache())) * scvVolume / dt;
380 tr.residual_[tIdx][I][
Free] += fLocalStorage;
383 const Scalar sStorageOfTimeIndex0 = sVol.value() * tr.concentration_[tIdx][I][
Solution];
384 const Scalar sLocalStorage = (sStorageOfTimeIndex0 - storage1_<Solution>(tr, tIdx, I, I1,
385 elemCtx.enableStorageCache())) * scvVolume / dt;
386 tr.residual_[tIdx][I][
Solution] += sLocalStorage;
390 (*tr.mat)[I][I][
Free][
Free] += fVol.derivative(0) * scvVolume/dt;
396 const ElementContext& elemCtx,
402 if (tr.numTracer() == 0) {
406 const auto& [fFlux, isUpF] = computeFlux_<Free>(tr.phaseIdx_, elemCtx, scvfIdx, 0);
407 const auto& [sFlux, isUpS] = computeFlux_<Solution>(tr.phaseIdx_, elemCtx, scvfIdx, 0);
409 dVol_[
Free][tr.phaseIdx_][I] += fFlux.value() * dt;
410 const int fGlobalUpIdx = isUpF ? I : J;
411 const int sGlobalUpIdx = isUpS ? I : J;
412 for (
int tIdx = 0; tIdx < tr.numTracer(); ++tIdx) {
414 tr.residual_[tIdx][I][
Free] += fFlux.value()*tr.concentration_[tIdx][fGlobalUpIdx][
Free];
415 tr.residual_[tIdx][I][
Solution] += sFlux.value()*tr.concentration_[tIdx][sGlobalUpIdx][
Solution];
420 (*tr.mat)[J][I][
Free][
Free] = -fFlux.derivative(0);
421 (*tr.mat)[I][I][
Free][
Free] += fFlux.derivative(0);
429 template<
class TrRe,
class Well>
433 if (tr.numTracer() == 0) {
437 const auto& eclWell = well.wellEcl();
443 auto* mswTracerRate = eclWell.isMultiSegment()
446 for (
int tIdx = 0; tIdx < tr.numTracer(); ++tIdx) {
447 tracerRate[tr.idx_[tIdx]] = {this->
name(tr.idx_[tIdx]), 0.0};
448 freeTracerRate[tr.idx_[tIdx]] = {this->
wellfname(tr.idx_[tIdx]), 0.0};
449 solTracerRate[tr.idx_[tIdx]] = {this->
wellsname(tr.idx_[tIdx]), 0.0};
450 if (eclWell.isMultiSegment()) {
451 auto& wtr = mswTracerRate->at(tr.idx_[tIdx]) = {this->
name(tr.idx_[tIdx])};;
452 wtr.rate.reserve(eclWell.getConnections().size());
453 for (std::size_t i = 0; i < eclWell.getConnections().size(); ++i) {
454 wtr.rate.emplace(eclWell.getConnections().get(i).segment(), 0.0);
459 std::vector<Scalar> wtracer(tr.numTracer());
460 for (
int tIdx = 0; tIdx < tr.numTracer(); ++tIdx) {
462 simulator_.problem().wellModel().summaryState());
466 const auto& ws =
simulator_.problem().wellModel().wellState().well(well.name());
467 for (std::size_t i = 0; i < ws.perf_data.size(); ++i) {
468 const auto I = ws.perf_data.cell_index[i];
469 const Scalar rate = well.volumetricSurfaceRateForConnection(I, tr.phaseIdx_);
471 if (tr.phaseIdx_ == FluidSystem::oilPhaseIdx && FluidSystem::enableVaporizedOil()) {
472 rate_s = ws.perf_data.phase_mixing_rates[i][ws.vaporized_oil];
474 else if (tr.phaseIdx_ == FluidSystem::gasPhaseIdx && FluidSystem::enableDissolvedGas()) {
475 rate_s = ws.perf_data.phase_mixing_rates[i][ws.dissolved_gas];
481 const Scalar rate_f = rate - rate_s;
483 for (
int tIdx = 0; tIdx < tr.numTracer(); ++tIdx) {
484 const Scalar delta = rate_f * wtracer[tIdx];
486 tr.residual_[tIdx][I][
Free] -= delta;
490 tracerRate[tr.idx_[tIdx]].rate += delta;
491 freeTracerRate[tr.idx_[tIdx]].rate += delta;
492 if (eclWell.isMultiSegment()) {
493 (*mswTracerRate)[tr.idx_[tIdx]].rate[eclWell.getConnections().get(i).segment()] += delta;
496 dVol_[
Free][tr.phaseIdx_][I] -= rate_f * dt;
498 else if (rate_f < 0) {
499 for (
int tIdx = 0; tIdx < tr.numTracer(); ++tIdx) {
500 const Scalar delta = rate_f * wtracer[tIdx];
503 tracerRate[tr.idx_[tIdx]].rate += delta;
504 freeTracerRate[tr.idx_[tIdx]].rate += delta;
507 tr.residual_[tIdx][I][
Free] -= rate_f * tr.concentration_[tIdx][I][
Free];
509 dVol_[
Free][tr.phaseIdx_][I] -= rate_f * dt;
512 (*tr.mat)[I][I][
Free][
Free] -= rate_f * variable<TracerEvaluation>(1.0, 0).derivative(0);
515 for (
int tIdx = 0; tIdx < tr.numTracer(); ++tIdx) {
517 tr.residual_[tIdx][I][
Solution] -= rate_s * tr.concentration_[tIdx][I][
Solution];
522 (*tr.mat)[I][I][
Solution][
Solution] -= rate_s * variable<TracerEvaluation>(1.0, 0).derivative(0);
532 if (tr.numTracer() == 0) {
537 if (tr.phaseIdx_ == FluidSystem::waterPhaseIdx ||
538 (tr.phaseIdx_ == FluidSystem::gasPhaseIdx && !FluidSystem::enableDissolvedGas()) ||
539 (tr.phaseIdx_ == FluidSystem::oilPhaseIdx && !FluidSystem::enableVaporizedOil()))
545 const Scalar& dfVol =
dVol_[
Free][tr.phaseIdx_][I];
548 for (
int tIdx = 0; tIdx < tr.numTracer(); ++tIdx) {
550 const auto delta = (dfVol / dt) * tr.concentration_[tIdx][I][
Free];
551 tr.residual_[tIdx][I][
Free] -= delta;
552 tr.residual_[tIdx][I][
Solution] += delta;
555 const auto delta = (dsVol / dt) * tr.concentration_[tIdx][I][
Solution];
556 tr.residual_[tIdx][I][
Free] += delta;
557 tr.residual_[tIdx][I][
Solution] -= delta;
563 const auto delta = (dfVol / dt) * variable<TracerEvaluation>(1.0, 0).derivative(0);
564 (*tr.mat)[I][I][
Free][
Free] -= delta;
568 const auto delta = (dsVol / dt) * variable<TracerEvaluation>(1.0, 0).derivative(0);
585 OPM_TIMEBLOCK(tracerAssemble);
587 if (tr.numTracer() != 0) {
589 for (
int tIdx = 0; tIdx < tr.numTracer(); ++tIdx) {
590 tr.residual_[tIdx] = 0.0;
604 const auto& wellPtrs =
simulator_.problem().wellModel().localNonshutWells();
609 for (
const auto& wellPtr : wellPtrs) {
611 const auto& eclWell = wellPtr->wellEcl();
615 auto* mswTracerRate = eclWell.isMultiSegment()
628 #pragma omp parallel for
632 const Scalar dt = elemCtx.simulator().timeStepSize();
634 for (
const auto& elem : chunk) {
635 elemCtx.updateStencil(elem);
636 const std::size_t I = elemCtx.globalSpaceIndex( 0, 0);
638 if (elem.partitionType() != Dune::InteriorEntity) {
642 for (
const auto& tr :
tbatch) {
643 if (tr.numTracer() != 0) {
644 (*tr.mat)[I][I][0][0] = 1.;
645 (*tr.mat)[I][I][1][1] = 1.;
650 elemCtx.updateAllIntensiveQuantities();
651 elemCtx.updateAllExtensiveQuantities();
653 const Scalar extrusionFactor =
654 elemCtx.intensiveQuantities( 0, 0).extrusionFactor();
655 Valgrind::CheckDefined(extrusionFactor);
656 assert(isfinite(extrusionFactor));
657 assert(extrusionFactor > 0.0);
658 const Scalar scvVolume =
659 elemCtx.stencil(0).subControlVolume( 0).volume() * extrusionFactor;
660 const std::size_t I1 = elemCtx.globalSpaceIndex( 0, 1);
665 if (tr.numTracer() == 0) {
671 const std::size_t numInteriorFaces = elemCtx.numInteriorFaces(0);
672 for (
unsigned scvfIdx = 0; scvfIdx < numInteriorFaces; scvfIdx++) {
673 const auto& face = elemCtx.stencil(0).interiorFace(scvfIdx);
674 const unsigned j = face.exteriorIndex();
675 const unsigned J = elemCtx.globalSpaceIndex( j, 0);
677 if (tr.numTracer() == 0) {
686 if (tr.numTracer() == 0) {
695 "assembleTracerEquations() failed: ",
700 if (tr.numTracer() == 0) {
705 simulator_.gridView().communicate(handle, Dune::InteriorBorder_All_Interface,
706 Dune::ForwardCommunication);
710 template<TracerTypeIdx Index,
class TrRe>
712 const Scalar scvVolume,
713 const unsigned globalDofIdx)
715 const Scalar vol1 = computeVolume_<Index>(tr.phaseIdx_, globalDofIdx, 0);
716 vol1_[
Index][tr.phaseIdx_][globalDofIdx] = vol1 * scvVolume;
717 dVol_[
Index][tr.phaseIdx_][globalDofIdx] = 0.0;
718 for (
int tIdx = 0; tIdx < tr.numTracer(); ++tIdx) {
719 tr.storageOfTimeIndex1_[tIdx][globalDofIdx][
Index] =
720 vol1 * tr.concentrationInitial_[tIdx][globalDofIdx][
Index];
727 if (tr.numTracer() != 0) {
728 tr.concentrationInitial_ = tr.concentration_;
734 #pragma omp parallel for
739 for (
const auto& elem : chunk) {
740 elemCtx.updatePrimaryStencil(elem);
741 elemCtx.updatePrimaryIntensiveQuantities(0);
742 const Scalar extrusionFactor = elemCtx.intensiveQuantities( 0, 0).extrusionFactor();
743 const Scalar scvVolume = elemCtx.stencil(0).subControlVolume( 0).volume() * extrusionFactor;
744 const unsigned globalDofIdx = elemCtx.globalSpaceIndex(0, 0);
747 if (tr.numTracer() == 0) {
752 updateElem<Free>(tr, scvVolume, globalDofIdx);
753 updateElem<Solution>(tr, scvVolume, globalDofIdx);
759 template<TracerTypeIdx Index,
class TrRe>
761 const std::vector<TracerVector>& dx,
764 const unsigned globalDofIdx,
765 std::vector<TracerVectorSingle>& sc)
767 constexpr Scalar tol_gas_sat = 1e-6;
768 tr.concentration_[tIdx][globalDofIdx][
Index] -= dx[tIdx][globalDofIdx][
Index];
769 if (tr.concentration_[tIdx][globalDofIdx][
Index] < 0.0 || S < tol_gas_sat) {
770 tr.concentration_[tIdx][globalDofIdx][
Index] = 0.0;
772 sc[tr.idx_[tIdx]][globalDofIdx] = tr.concentration_[tIdx][globalDofIdx][
Index];
775 template<TracerTypeIdx Index,
class TrRe>
786 for (
int tIdx = 0; tIdx < tr.numTracer(); ++tIdx) {
788 const Scalar delta = rate * tr.concentration_[tIdx][I][
Index];
789 tracerRate[tr.idx_[tIdx]].rate += delta;
790 splitRate[tr.idx_[tIdx]].rate += delta;
791 if (eclWell.isMultiSegment()) {
792 (*mswTracerRate)[tr.idx_[tIdx]].rate[eclWell.getConnections().get(i).segment()] += delta;
803 if (tr.numTracer() == 0) {
809 std::vector<TracerVector> dx(tr.concentration_);
810 for (
int tIdx = 0; tIdx < tr.numTracer(); ++tIdx) {
816 OpmLog::warning(
"### Tracer model: Linear solver did not converge. ###");
819 OPM_TIMEBLOCK(tracerPost);
821 for (
int tIdx = 0; tIdx < tr.numTracer(); ++tIdx) {
822 for (std::size_t globalDofIdx = 0; globalDofIdx < tr.concentration_[tIdx].size(); ++globalDofIdx) {
825 const auto& intQuants =
simulator_.model().intensiveQuantities(globalDofIdx, 0);
826 const auto& fs = intQuants.fluidState();
827 const Scalar Sf = decay<Scalar>(fs.saturation(tr.phaseIdx_));
830 if (tr.phaseIdx_ == FluidSystem::gasPhaseIdx && FluidSystem::enableDissolvedGas()) {
831 Ss = decay<Scalar>(fs.saturation(FluidSystem::oilPhaseIdx));
833 else if (tr.phaseIdx_ == FluidSystem::oilPhaseIdx && FluidSystem::enableVaporizedOil()) {
834 Ss = decay<Scalar>(fs.saturation(FluidSystem::gasPhaseIdx));
843 const auto& wellPtrs =
simulator_.problem().wellModel().localNonshutWells();
844 for (
const auto& wellPtr : wellPtrs) {
845 const auto& eclWell = wellPtr->wellEcl();
848 if (!eclWell.isProducer()) {
852 Scalar rateWellPos = 0.0;
853 Scalar rateWellNeg = 0.0;
854 const std::size_t well_index =
simulator_.problem().wellModel().wellState().index(eclWell.name()).value();
855 const auto& ws =
simulator_.problem().wellModel().wellState().well(well_index);
859 auto* mswTracerRate = eclWell.isMultiSegment() ? &this->
mSwTracerRate_[eclWell.seqIndex()] :
nullptr;
861 for (std::size_t i = 0; i < ws.perf_data.size(); ++i) {
862 const auto I = ws.perf_data.cell_index[i];
863 const Scalar rate = wellPtr->volumetricSurfaceRateForConnection(I, tr.phaseIdx_);
866 if (tr.phaseIdx_ == FluidSystem::oilPhaseIdx && FluidSystem::enableVaporizedOil()) {
867 rate_s = ws.perf_data.phase_mixing_rates[i][ws.vaporized_oil];
869 else if (tr.phaseIdx_ == FluidSystem::gasPhaseIdx && FluidSystem::enableDissolvedGas()) {
870 rate_s = ws.perf_data.phase_mixing_rates[i][ws.dissolved_gas];
876 const Scalar rate_f = rate - rate_s;
877 assignRates<Free>(tr, eclWell, i, I, rate_f,
878 tracerRate, mswTracerRate, freeTracerRate);
879 assignRates<Solution>(tr, eclWell, i, I, rate_s,
880 tracerRate, mswTracerRate, solTracerRate);
895 const Scalar official_well_rate_total =
896 simulator_.problem().wellModel().wellState().well(well_index).surface_rates[tr.phaseIdx_];
898 const Scalar rateWellTotal = official_well_rate_total;
900 if (rateWellTotal > rateWellNeg) {
901 constexpr Scalar bucketPrDay = 10.0 / (1000. * 3600. * 24.);
902 const Scalar factor = (rateWellTotal < -bucketPrDay) ? rateWellTotal / rateWellNeg : 0.0;
903 for (
int tIdx = 0; tIdx < tr.numTracer(); ++tIdx) {
904 tracerRate[tIdx].rate *= factor;
921 template <
typename TV>
930 std::unique_ptr<TracerMatrix>
mat;
941 result.
idx_ = {1,2,3};
950 template<
class Serializer>
960 {
return idx_.size(); }
964 const int numGridDof = concentration.size();
965 idx_.emplace_back(idx);
973 std::array<TracerBatch<TracerVector>,numPhases>
tbatch;
977 std::array<std::array<std::vector<Scalar>,numPhases>,2>
vol1_;
978 std::array<std::array<std::vector<Scalar>,numPhases>,2>
dVol_;
#define OPM_END_PARALLEL_TRY_CATCH_LOG(obptc_logger, obptc_prefix, obptc_output, comm)
Catch exception, log, and throw in a parallel try-catch clause.
Definition: DeferredLoggingErrorHelpers.hpp:202
#define OPM_BEGIN_PARALLEL_TRY_CATCH()
Macro to setup the try of a parallel try-catch.
Definition: DeferredLoggingErrorHelpers.hpp:158
A datahandle sending data located in multiple vectors.
Definition: DeferredLogger.hpp:57
Definition: GenericTracerModel.hpp:56
Opm::GenericTracerModel< GetPropType< TypeTag, Properties::Grid >, GetPropType< TypeTag, Properties::GridView >, GetPropType< TypeTag, Properties::DofMapper >, GetPropType< TypeTag, Properties::Stencil >, GetPropType< TypeTag, Properties::FluidSystem >, GetPropType< TypeTag, Properties::Scalar > >::numTracers int numTracers() const
Return the number of tracers considered by the tracerModel.
Definition: GenericTracerModel_impl.hpp:162
Opm::GenericTracerModel< GetPropType< TypeTag, Properties::Grid >, GetPropType< TypeTag, Properties::GridView >, GetPropType< TypeTag, Properties::DofMapper >, GetPropType< TypeTag, Properties::Stencil >, GetPropType< TypeTag, Properties::FluidSystem >, GetPropType< TypeTag, Properties::Scalar > >::doInit void doInit(bool rst, std::size_t numGridDof, std::size_t gasPhaseIdx, std::size_t oilPhaseIdx, std::size_t waterPhaseIdx)
Initialize all internal data structures needed by the tracer module.
Definition: GenericTracerModel_impl.hpp:227
Opm::GenericTracerModel< GetPropType< TypeTag, Properties::Grid >, GetPropType< TypeTag, Properties::GridView >, GetPropType< TypeTag, Properties::DofMapper >, GetPropType< TypeTag, Properties::Stencil >, GetPropType< TypeTag, Properties::FluidSystem >, GetPropType< TypeTag, Properties::Scalar > >::tracerConcentration_ std::vector< TracerVector > tracerConcentration_
Definition: GenericTracerModel.hpp:158
Opm::GenericTracerModel< GetPropType< TypeTag, Properties::Grid >, GetPropType< TypeTag, Properties::GridView >, GetPropType< TypeTag, Properties::DofMapper >, GetPropType< TypeTag, Properties::Stencil >, GetPropType< TypeTag, Properties::FluidSystem >, GetPropType< TypeTag, Properties::Scalar > >::tracerPhaseIdx_ std::vector< int > tracerPhaseIdx_
Definition: GenericTracerModel.hpp:156
Opm::GenericTracerModel< GetPropType< TypeTag, Properties::Grid >, GetPropType< TypeTag, Properties::GridView >, GetPropType< TypeTag, Properties::DofMapper >, GetPropType< TypeTag, Properties::Stencil >, GetPropType< TypeTag, Properties::FluidSystem >, GetPropType< TypeTag, Properties::Scalar > >::wellFreeTracerRate_ std::unordered_map< int, std::vector< WellTracerRate< GetPropType< TypeTag, Properties::Scalar > > > > wellFreeTracerRate_
Definition: GenericTracerModel.hpp:165
Opm::GenericTracerModel< GetPropType< TypeTag, Properties::Grid >, GetPropType< TypeTag, Properties::GridView >, GetPropType< TypeTag, Properties::DofMapper >, GetPropType< TypeTag, Properties::Stencil >, GetPropType< TypeTag, Properties::FluidSystem >, GetPropType< TypeTag, Properties::Scalar > >::wellTracerRate_ std::unordered_map< int, std::vector< WellTracerRate< GetPropType< TypeTag, Properties::Scalar > > > > wellTracerRate_
Definition: GenericTracerModel.hpp:164
Opm::GenericTracerModel< GetPropType< TypeTag, Properties::Grid >, GetPropType< TypeTag, Properties::GridView >, GetPropType< TypeTag, Properties::DofMapper >, GetPropType< TypeTag, Properties::Stencil >, GetPropType< TypeTag, Properties::FluidSystem >, GetPropType< TypeTag, Properties::Scalar > >::mSwTracerRate_ std::unordered_map< int, std::vector< MSWellTracerRate< GetPropType< TypeTag, Properties::Scalar > > > > mSwTracerRate_
Definition: GenericTracerModel.hpp:168
Opm::GenericTracerModel< GetPropType< TypeTag, Properties::Grid >, GetPropType< TypeTag, Properties::GridView >, GetPropType< TypeTag, Properties::DofMapper >, GetPropType< TypeTag, Properties::Stencil >, GetPropType< TypeTag, Properties::FluidSystem >, GetPropType< TypeTag, Properties::Scalar > >::wellfname std::string wellfname(int tracerIdx) const
Definition: GenericTracerModel_impl.hpp:183
Opm::GenericTracerModel< GetPropType< TypeTag, Properties::Grid >, GetPropType< TypeTag, Properties::GridView >, GetPropType< TypeTag, Properties::DofMapper >, GetPropType< TypeTag, Properties::Stencil >, GetPropType< TypeTag, Properties::FluidSystem >, GetPropType< TypeTag, Properties::Scalar > >::name const std::string & name(int tracerIdx) const
Return the tracer name.
Definition: GenericTracerModel_impl.hpp:220
Opm::GenericTracerModel< GetPropType< TypeTag, Properties::Grid >, GetPropType< TypeTag, Properties::GridView >, GetPropType< TypeTag, Properties::DofMapper >, GetPropType< TypeTag, Properties::Stencil >, GetPropType< TypeTag, Properties::FluidSystem >, GetPropType< TypeTag, Properties::Scalar > >::wellsname std::string wellsname(int tracerIdx) const
Definition: GenericTracerModel_impl.hpp:190
Opm::GenericTracerModel< GetPropType< TypeTag, Properties::Grid >, GetPropType< TypeTag, Properties::GridView >, GetPropType< TypeTag, Properties::DofMapper >, GetPropType< TypeTag, Properties::Stencil >, GetPropType< TypeTag, Properties::FluidSystem >, GetPropType< TypeTag, Properties::Scalar > >::freeTracerConcentration_ std::vector< TracerVectorSingle > freeTracerConcentration_
Definition: GenericTracerModel.hpp:160
Opm::GenericTracerModel< GetPropType< TypeTag, Properties::Grid >, GetPropType< TypeTag, Properties::GridView >, GetPropType< TypeTag, Properties::DofMapper >, GetPropType< TypeTag, Properties::Stencil >, GetPropType< TypeTag, Properties::FluidSystem >, GetPropType< TypeTag, Properties::Scalar > >::TracerVector Dune::BlockVector< Dune::FieldVector< GetPropType< TypeTag, Properties::Scalar >, 2 > > TracerVector
Definition: GenericTracerModel.hpp:60
Opm::GenericTracerModel< GetPropType< TypeTag, Properties::Grid >, GetPropType< TypeTag, Properties::GridView >, GetPropType< TypeTag, Properties::DofMapper >, GetPropType< TypeTag, Properties::Stencil >, GetPropType< TypeTag, Properties::FluidSystem >, GetPropType< TypeTag, Properties::Scalar > >::TracerTypeIdx TracerTypeIdx
Tracer type index.
Definition: GenericTracerModel.hpp:146
Opm::GenericTracerModel< GetPropType< TypeTag, Properties::Grid >, GetPropType< TypeTag, Properties::GridView >, GetPropType< TypeTag, Properties::DofMapper >, GetPropType< TypeTag, Properties::Stencil >, GetPropType< TypeTag, Properties::FluidSystem >, GetPropType< TypeTag, Properties::Scalar > >::Free @ Free
Definition: GenericTracerModel.hpp:147
Opm::GenericTracerModel< GetPropType< TypeTag, Properties::Grid >, GetPropType< TypeTag, Properties::GridView >, GetPropType< TypeTag, Properties::DofMapper >, GetPropType< TypeTag, Properties::Stencil >, GetPropType< TypeTag, Properties::FluidSystem >, GetPropType< TypeTag, Properties::Scalar > >::Solution @ Solution
Definition: GenericTracerModel.hpp:148
Opm::GenericTracerModel< GetPropType< TypeTag, Properties::Grid >, GetPropType< TypeTag, Properties::GridView >, GetPropType< TypeTag, Properties::DofMapper >, GetPropType< TypeTag, Properties::Stencil >, GetPropType< TypeTag, Properties::FluidSystem >, GetPropType< TypeTag, Properties::Scalar > >::linearSolveBatchwise_ bool linearSolveBatchwise_(const TracerMatrix &M, std::vector< TracerVector > &x, std::vector< TracerVector > &b)
Definition: GenericTracerModel_impl.hpp:444
Opm::GenericTracerModel< GetPropType< TypeTag, Properties::Grid >, GetPropType< TypeTag, Properties::GridView >, GetPropType< TypeTag, Properties::DofMapper >, GetPropType< TypeTag, Properties::Stencil >, GetPropType< TypeTag, Properties::FluidSystem >, GetPropType< TypeTag, Properties::Scalar > >::TracerMatrix Dune::BCRSMatrix< Opm::MatrixBlock< GetPropType< TypeTag, Properties::Scalar >, 2, 2 > > TracerMatrix
Definition: GenericTracerModel.hpp:59
Opm::GenericTracerModel< GetPropType< TypeTag, Properties::Grid >, GetPropType< TypeTag, Properties::GridView >, GetPropType< TypeTag, Properties::DofMapper >, GetPropType< TypeTag, Properties::Stencil >, GetPropType< TypeTag, Properties::FluidSystem >, GetPropType< TypeTag, Properties::Scalar > >::solTracerConcentration_ std::vector< TracerVectorSingle > solTracerConcentration_
Definition: GenericTracerModel.hpp:161
Opm::GenericTracerModel< GetPropType< TypeTag, Properties::Grid >, GetPropType< TypeTag, Properties::GridView >, GetPropType< TypeTag, Properties::DofMapper >, GetPropType< TypeTag, Properties::Stencil >, GetPropType< TypeTag, Properties::FluidSystem >, GetPropType< TypeTag, Properties::Scalar > >::currentConcentration_ GetPropType< TypeTag, Properties::Scalar > currentConcentration_(const Well &eclWell, const std::string &trName, const SummaryState &summaryState) const
Definition: GenericTracerModel_impl.hpp:211
Opm::GenericTracerModel< GetPropType< TypeTag, Properties::Grid >, GetPropType< TypeTag, Properties::GridView >, GetPropType< TypeTag, Properties::DofMapper >, GetPropType< TypeTag, Properties::Stencil >, GetPropType< TypeTag, Properties::FluidSystem >, GetPropType< TypeTag, Properties::Scalar > >::TracerVectorSingle Dune::BlockVector< Dune::FieldVector< GetPropType< TypeTag, Properties::Scalar >, 1 > > TracerVectorSingle
Definition: GenericTracerModel.hpp:58
Opm::GenericTracerModel< GetPropType< TypeTag, Properties::Grid >, GetPropType< TypeTag, Properties::GridView >, GetPropType< TypeTag, Properties::DofMapper >, GetPropType< TypeTag, Properties::Stencil >, GetPropType< TypeTag, Properties::FluidSystem >, GetPropType< TypeTag, Properties::Scalar > >::tracerMatrix_ std::unique_ptr< TracerMatrix > tracerMatrix_
Definition: GenericTracerModel.hpp:159
Opm::GenericTracerModel< GetPropType< TypeTag, Properties::Grid >, GetPropType< TypeTag, Properties::GridView >, GetPropType< TypeTag, Properties::DofMapper >, GetPropType< TypeTag, Properties::Stencil >, GetPropType< TypeTag, Properties::FluidSystem >, GetPropType< TypeTag, Properties::Scalar > >::wellSolTracerRate_ std::unordered_map< int, std::vector< WellTracerRate< GetPropType< TypeTag, Properties::Scalar > > > > wellSolTracerRate_
Definition: GenericTracerModel.hpp:166
static unsigned maxThreads()
Return the maximum number of threads of the current process.
Definition: threadmanager.hpp:66
A class which handles tracers as specified in by ECL.
Definition: TracerModel.hpp:76
void updateStorageCache()
Definition: TracerModel.hpp:724
std::array< std::array< std::vector< Scalar >, numPhases >, 2 > vol1_
Definition: TracerModel.hpp:977
TracerBatch< TracerVector > & oil_
Definition: TracerModel.hpp:975
void advanceTracerFields()
Definition: TracerModel.hpp:798
void init(bool rst)
Definition: TracerModel.hpp:139
void assembleTracerEquationSource(TrRe &tr, const Scalar dt, unsigned I)
Definition: TracerModel.hpp:528
TracerModel(Simulator &simulator)
Definition: TracerModel.hpp:106
void assignRates(const TrRe &tr, const Well &eclWell, const std::size_t i, const std::size_t I, const Scalar rate, std::vector< WellTracerRate< Scalar > > &tracerRate, std::vector< MSWellTracerRate< Scalar > > *mswTracerRate, std::vector< WellTracerRate< Scalar > > &splitRate)
Definition: TracerModel.hpp:776
void copyForOutput(TrRe &tr, const std::vector< TracerVector > &dx, const Scalar S, const unsigned tIdx, const unsigned globalDofIdx, std::vector< TracerVectorSingle > &sc)
Definition: TracerModel.hpp:760
void assembleTracerEquationFlux(TrRe &tr, const ElementContext &elemCtx, unsigned scvfIdx, unsigned I, unsigned J, const Scalar dt)
Definition: TracerModel.hpp:395
void assembleTracerEquationWell(TrRe &tr, const Well &well)
Definition: TracerModel.hpp:430
void updateElem(TrRe &tr, const Scalar scvVolume, const unsigned globalDofIdx)
Definition: TracerModel.hpp:711
Scalar computeVolume_(const int tracerPhaseIdx, const unsigned globalDofIdx, const unsigned timeIdx) const
Definition: TracerModel.hpp:253
typename BaseType::TracerTypeIdx TracerTypeIdx
Definition: TracerModel.hpp:247
void assembleTracerEquationVolume(TrRe &tr, const ElementContext &elemCtx, const Scalar scvVolume, const Scalar dt, unsigned I, unsigned I1)
Definition: TracerModel.hpp:359
Scalar storage1_(const TrRe &tr, const unsigned tIdx, const unsigned I, const unsigned I1, const bool cache)
Definition: TracerModel.hpp:344
void serialize(Restarter &)
This method writes the complete state of all tracer to the hard disk.
Definition: TracerModel.hpp:226
TracerBatch< TracerVector > & gas_
Definition: TracerModel.hpp:976
void assembleTracerEquations_()
Definition: TracerModel.hpp:574
std::array< TracerBatch< TracerVector >, numPhases > tbatch
Definition: TracerModel.hpp:973
std::pair< TracerEvaluation, bool > computeFlux_(const int tracerPhaseIdx, const ElementContext &elemCtx, const unsigned scvfIdx, const unsigned timeIdx) const
Definition: TracerModel.hpp:291
void beginTimeStep()
Definition: TracerModel.hpp:198
void serializeOp(Serializer &serializer)
Definition: TracerModel.hpp:240
void prepareTracerBatches()
Definition: TracerModel.hpp:145
std::array< std::array< std::vector< Scalar >, numPhases >, 2 > dVol_
Definition: TracerModel.hpp:978
ElementChunks< GridView, Dune::Partitions::All > element_chunks_
Definition: TracerModel.hpp:979
TracerBatch< TracerVector > & wat_
Definition: TracerModel.hpp:974
void endTimeStep()
Informs the tracer model that a time step has just been finished.
Definition: TracerModel.hpp:211
Simulator & simulator_
Definition: TracerModel.hpp:911
void deserialize(Restarter &)
This method restores the complete state of the tracer from disk.
Definition: TracerModel.hpp:236
A data handle sending multiple data store in vectors attached to cells.
Definition: VectorVectorDataHandle.hpp:50
int Index
The type of an index of a degree of freedom.
Definition: overlaptypes.hh:44
Definition: blackoilmodel.hh:77
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: WellTracerRate.hpp:60
Definition: TracerModel.hpp:56
a tag to mark properties as undefined
Definition: propertysystem.hh:38
Definition: TracerModel.hpp:923
const int phaseIdx_
Definition: TracerModel.hpp:925
void addTracer(const int idx, const TV &concentration)
Definition: TracerModel.hpp:962
bool operator==(const TracerBatch &rhs) const
Definition: TracerModel.hpp:932
static TracerBatch serializationTestObject()
Definition: TracerModel.hpp:938
std::vector< TV > concentrationInitial_
Definition: TracerModel.hpp:926
void serializeOp(Serializer &serializer)
Definition: TracerModel.hpp:951
TracerBatch(int phaseIdx=0)
Definition: TracerModel.hpp:957
int numTracer() const
Definition: TracerModel.hpp:959
std::vector< int > idx_
Definition: TracerModel.hpp:924
std::unique_ptr< TracerMatrix > mat
Definition: TracerModel.hpp:930
std::vector< TV > storageOfTimeIndex1_
Definition: TracerModel.hpp:928
std::vector< TV > residual_
Definition: TracerModel.hpp:929
std::vector< TV > concentration_
Definition: TracerModel.hpp:927
Definition: WellTracerRate.hpp:33