69class TracerModel :
public GenericTracerModel<GetPropType<TypeTag, Properties::Grid>,
70 GetPropType<TypeTag, Properties::GridView>,
71 GetPropType<TypeTag, Properties::DofMapper>,
72 GetPropType<TypeTag, Properties::Stencil>,
73 GetPropType<TypeTag, Properties::FluidSystem>,
74 GetPropType<TypeTag, Properties::Scalar>>
76 using BaseType = GenericTracerModel<GetPropType<TypeTag, Properties::Grid>,
92 using TracerEvaluation = DenseAd::Evaluation<Scalar,1>;
94 using TracerMatrix =
typename BaseType::TracerMatrix;
95 using TracerVector =
typename BaseType::TracerVector;
96 using TracerVectorSingle =
typename BaseType::TracerVectorSingle;
99 enum { numPhases = FluidSystem::numPhases };
100 enum { waterPhaseIdx = FluidSystem::waterPhaseIdx };
101 enum { oilPhaseIdx = FluidSystem::oilPhaseIdx };
102 enum { gasPhaseIdx = FluidSystem::gasPhaseIdx };
105 explicit TracerModel(Simulator& simulator)
106 : BaseType(simulator.vanguard().gridView(),
107 simulator.vanguard().eclState(),
108 simulator.vanguard().cartesianIndexMapper(),
109 simulator.model().dofMapper(),
110 simulator.vanguard().cellCentroids())
111 , simulator_(simulator)
112 , tbatch({waterPhaseIdx, oilPhaseIdx, gasPhaseIdx})
140 this->
doInit(rst, simulator_.model().numGridDof(),
141 gasPhaseIdx, oilPhaseIdx, waterPhaseIdx);
144 void prepareTracerBatches()
146 for (std::size_t tracerIdx = 0; tracerIdx < this->tracerPhaseIdx_.size(); ++tracerIdx) {
147 if (this->tracerPhaseIdx_[tracerIdx] == FluidSystem::waterPhaseIdx) {
148 if (! FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx)){
149 throw std::runtime_error(
"Water tracer specified for non-water fluid system: " +
150 this->
name(tracerIdx));
153 wat_.addTracer(tracerIdx, this->tracerConcentration_[tracerIdx]);
155 else if (this->tracerPhaseIdx_[tracerIdx] == FluidSystem::oilPhaseIdx) {
156 if (! FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx)){
157 throw std::runtime_error(
"Oil tracer specified for non-oil fluid system: " +
158 this->
name(tracerIdx));
161 oil_.addTracer(tracerIdx, this->tracerConcentration_[tracerIdx]);
163 else if (this->tracerPhaseIdx_[tracerIdx] == FluidSystem::gasPhaseIdx) {
164 if (! FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx)){
165 throw std::runtime_error(
"Gas tracer specified for non-gas fluid system: " +
166 this->
name(tracerIdx));
169 gas_.addTracer(tracerIdx, this->tracerConcentration_[tracerIdx]);
173 vol1_[0][this->tracerPhaseIdx_[tracerIdx]].
174 resize(this->freeTracerConcentration_[tracerIdx].size());
175 vol1_[1][this->tracerPhaseIdx_[tracerIdx]].
176 resize(this->freeTracerConcentration_[tracerIdx].size());
177 dVol_[0][this->tracerPhaseIdx_[tracerIdx]].
178 resize(this->solTracerConcentration_[tracerIdx].size());
179 dVol_[1][this->tracerPhaseIdx_[tracerIdx]].
180 resize(this->solTracerConcentration_[tracerIdx].size());
184 TracerMatrix* base = this->tracerMatrix_.get();
185 for (
auto& tr : this->tbatch) {
186 if (tr.numTracer() != 0) {
187 if (this->tracerMatrix_) {
188 tr.mat = std::move(this->tracerMatrix_);
191 tr.mat = std::make_unique<TracerMatrix>(*base);
203 OPM_TIMEBLOCK(tracerUpdateCache);
204 updateStorageCache();
216 OPM_TIMEBLOCK(tracerAdvance);
217 advanceTracerFields();
224 template <
class Restarter>
234 template <
class Restarter>
238 template<
class Serializer>
239 void serializeOp(Serializer& serializer)
241 serializer(
static_cast<BaseType&
>(*
this));
247 using BaseType::Free;
248 using BaseType::Solution;
251 template<TracerTypeIdx Index>
252 Scalar computeVolume_(
const int tracerPhaseIdx,
253 const unsigned globalDofIdx,
254 const unsigned timeIdx)
const
256 const auto& intQuants = simulator_.model().intensiveQuantities(globalDofIdx, timeIdx);
257 const auto& fs = intQuants.fluidState();
258 constexpr Scalar min_volume = 1e-10;
260 if constexpr (Index == Free) {
261 return std::max(decay<Scalar>(fs.saturation(tracerPhaseIdx)) *
262 decay<Scalar>(fs.invB(tracerPhaseIdx)) *
263 decay<Scalar>(intQuants.porosity()),
267 if (tracerPhaseIdx == FluidSystem::oilPhaseIdx && FluidSystem::enableVaporizedOil()) {
268 return std::max(decay<Scalar>(fs.saturation(FluidSystem::gasPhaseIdx)) *
269 decay<Scalar>(fs.invB(FluidSystem::gasPhaseIdx)) *
270 decay<Scalar>(fs.Rv()) *
271 decay<Scalar>(intQuants.porosity()),
276 else if (tracerPhaseIdx == FluidSystem::gasPhaseIdx && FluidSystem::enableDissolvedGas()) {
277 return std::max(decay<Scalar>(fs.saturation(FluidSystem::oilPhaseIdx)) *
278 decay<Scalar>(fs.invB(FluidSystem::oilPhaseIdx)) *
279 decay<Scalar>(fs.Rs()) *
280 decay<Scalar>(intQuants.porosity()),
288 template<TracerTypeIdx Index>
289 std::pair<TracerEvaluation, bool>
290 computeFlux_(
const int tracerPhaseIdx,
291 const ElementContext& elemCtx,
292 const unsigned scvfIdx,
293 const unsigned timeIdx)
const
295 const auto& stencil = elemCtx.stencil(timeIdx);
296 const auto& scvf = stencil.interiorFace(scvfIdx);
298 const auto& extQuants = elemCtx.extensiveQuantities(scvfIdx, timeIdx);
299 const unsigned inIdx = extQuants.interiorIndex();
304 if constexpr (
Index == Free) {
305 upIdx = extQuants.upstreamIndex(tracerPhaseIdx);
306 const auto& intQuants = elemCtx.intensiveQuantities(upIdx, timeIdx);
307 const auto& fs = intQuants.fluidState();
308 v = decay<Scalar>(extQuants.volumeFlux(tracerPhaseIdx)) *
309 decay<Scalar>(fs.invB(tracerPhaseIdx));
311 if (tracerPhaseIdx == FluidSystem::oilPhaseIdx && FluidSystem::enableVaporizedOil()) {
312 upIdx = extQuants.upstreamIndex(FluidSystem::gasPhaseIdx);
314 const auto& intQuants = elemCtx.intensiveQuantities(upIdx, timeIdx);
315 const auto& fs = intQuants.fluidState();
316 v = decay<Scalar>(fs.invB(FluidSystem::gasPhaseIdx)) *
317 decay<Scalar>(extQuants.volumeFlux(FluidSystem::gasPhaseIdx)) *
318 decay<Scalar>(fs.Rv());
321 else if (tracerPhaseIdx == FluidSystem::gasPhaseIdx && FluidSystem::enableDissolvedGas()) {
322 upIdx = extQuants.upstreamIndex(FluidSystem::oilPhaseIdx);
324 const auto& intQuants = elemCtx.intensiveQuantities(upIdx, timeIdx);
325 const auto& fs = intQuants.fluidState();
326 v = decay<Scalar>(fs.invB(FluidSystem::oilPhaseIdx)) *
327 decay<Scalar>(extQuants.volumeFlux(FluidSystem::oilPhaseIdx)) *
328 decay<Scalar>(fs.Rs());
336 const Scalar A = scvf.area();
337 return inIdx == upIdx
338 ? std::pair{A * v * variable<TracerEvaluation>(1.0, 0),
true}
339 : std::pair{A * v,
false};
342 template<TracerTypeIdx Index,
class TrRe>
343 Scalar storage1_(
const TrRe& tr,
350 return tr.storageOfTimeIndex1_[tIdx][I][
Index];
352 return computeVolume_<Index>(tr.phaseIdx_, I1, 1) *
353 tr.concentration_[tIdx][I1][
Index];
358 void assembleTracerEquationVolume(TrRe& tr,
359 const ElementContext& elemCtx,
360 const Scalar scvVolume,
366 if (tr.numTracer() == 0) {
370 const TracerEvaluation fVol = computeVolume_<Free>(tr.phaseIdx_, I, 0) * variable<TracerEvaluation>(1.0, 0);
371 const TracerEvaluation sVol = computeVolume_<Solution>(tr.phaseIdx_, I, 0) * variable<TracerEvaluation>(1.0, 0);
372 dVol_[Solution][tr.phaseIdx_][I] += sVol.value() * scvVolume - vol1_[1][tr.phaseIdx_][I];
373 dVol_[Free][tr.phaseIdx_][I] += fVol.value() * scvVolume - vol1_[0][tr.phaseIdx_][I];
374 for (
int tIdx = 0; tIdx < tr.numTracer(); ++tIdx) {
376 const Scalar fStorageOfTimeIndex0 = fVol.value() * tr.concentration_[tIdx][I][Free];
377 const Scalar fLocalStorage = (fStorageOfTimeIndex0 - storage1_<Free>(tr, tIdx, I, I1,
378 elemCtx.enableStorageCache())) * scvVolume / dt;
379 tr.residual_[tIdx][I][Free] += fLocalStorage;
382 const Scalar sStorageOfTimeIndex0 = sVol.value() * tr.concentration_[tIdx][I][Solution];
383 const Scalar sLocalStorage = (sStorageOfTimeIndex0 - storage1_<Solution>(tr, tIdx, I, I1,
384 elemCtx.enableStorageCache())) * scvVolume / dt;
385 tr.residual_[tIdx][I][Solution] += sLocalStorage;
389 (*tr.mat)[I][I][Free][Free] += fVol.derivative(0) * scvVolume/dt;
390 (*tr.mat)[I][I][Solution][Solution] += sVol.derivative(0) * scvVolume/dt;
394 void assembleTracerEquationFlux(TrRe& tr,
395 const ElementContext& elemCtx,
401 if (tr.numTracer() == 0) {
405 const auto& [fFlux, isUpF] = computeFlux_<Free>(tr.phaseIdx_, elemCtx, scvfIdx, 0);
406 const auto& [sFlux, isUpS] = computeFlux_<Solution>(tr.phaseIdx_, elemCtx, scvfIdx, 0);
407 dVol_[Solution][tr.phaseIdx_][I] += sFlux.value() * dt;
408 dVol_[Free][tr.phaseIdx_][I] += fFlux.value() * dt;
409 const int fGlobalUpIdx = isUpF ? I : J;
410 const int sGlobalUpIdx = isUpS ? I : J;
411 for (
int tIdx = 0; tIdx < tr.numTracer(); ++tIdx) {
413 tr.residual_[tIdx][I][Free] += fFlux.value()*tr.concentration_[tIdx][fGlobalUpIdx][Free];
414 tr.residual_[tIdx][I][Solution] += sFlux.value()*tr.concentration_[tIdx][sGlobalUpIdx][Solution];
419 (*tr.mat)[J][I][Free][Free] = -fFlux.derivative(0);
420 (*tr.mat)[I][I][Free][Free] += fFlux.derivative(0);
423 (*tr.mat)[J][I][Solution][Solution] = -sFlux.derivative(0);
424 (*tr.mat)[I][I][Solution][Solution] += sFlux.derivative(0);
428 template<
class TrRe,
class Well>
429 void assembleTracerEquationWell(TrRe& tr,
432 if (tr.numTracer() == 0) {
436 const auto& eclWell = well.wellEcl();
439 auto& tracerRate = this->wellTracerRate_[eclWell.seqIndex()];
440 tracerRate.reserve(tr.numTracer());
441 auto& solTracerRate = this->wellSolTracerRate_[eclWell.seqIndex()];
442 solTracerRate.reserve(tr.numTracer());
443 auto& freeTracerRate = this->wellFreeTracerRate_[eclWell.seqIndex()];
444 freeTracerRate.reserve(tr.numTracer());
445 auto* mswTracerRate = eclWell.isMultiSegment()
446 ? &this->mSwTracerRate_[eclWell.seqIndex()]
449 mswTracerRate->reserve(tr.numTracer());
451 for (
int tIdx = 0; tIdx < tr.numTracer(); ++tIdx) {
452 tracerRate.emplace_back(this->
name(tr.idx_[tIdx]), 0.0);
453 freeTracerRate.emplace_back(this->wellfname(tr.idx_[tIdx]), 0.0);
454 solTracerRate.emplace_back(this->wellsname(tr.idx_[tIdx]), 0.0);
455 if (eclWell.isMultiSegment()) {
456 auto& wtr = mswTracerRate->emplace_back(this->
name(tr.idx_[tIdx]));
457 wtr.rate.reserve(eclWell.getConnections().size());
458 for (std::size_t i = 0; i < eclWell.getConnections().size(); ++i) {
459 wtr.rate.emplace(eclWell.getConnections().get(i).segment(), 0.0);
464 std::vector<Scalar> wtracer(tr.numTracer());
465 for (
int tIdx = 0; tIdx < tr.numTracer(); ++tIdx) {
466 wtracer[tIdx] = this->currentConcentration_(eclWell, this->
name(tr.idx_[tIdx]),
467 simulator_.problem().wellModel().summaryState());
470 const Scalar dt = simulator_.timeStepSize();
471 const auto& ws = simulator_.problem().wellModel().wellState().well(well.name());
472 for (std::size_t i = 0; i < ws.perf_data.size(); ++i) {
473 const auto I = ws.perf_data.cell_index[i];
474 const Scalar rate = well.volumetricSurfaceRateForConnection(I, tr.phaseIdx_);
476 if (tr.phaseIdx_ == FluidSystem::oilPhaseIdx && FluidSystem::enableVaporizedOil()) {
477 rate_s = ws.perf_data.phase_mixing_rates[i][ws.vaporized_oil];
479 else if (tr.phaseIdx_ == FluidSystem::gasPhaseIdx && FluidSystem::enableDissolvedGas()) {
480 rate_s = ws.perf_data.phase_mixing_rates[i][ws.dissolved_gas];
486 const Scalar rate_f = rate - rate_s;
488 for (
int tIdx = 0; tIdx < tr.numTracer(); ++tIdx) {
489 const Scalar delta = rate_f * wtracer[tIdx];
491 tr.residual_[tIdx][I][Free] -= delta;
495 tracerRate[tIdx].rate += delta;
496 freeTracerRate[tIdx].rate += delta;
497 if (eclWell.isMultiSegment()) {
498 (*mswTracerRate)[tIdx].rate[eclWell.getConnections().get(i).segment()] += delta;
501 dVol_[Free][tr.phaseIdx_][I] -= rate_f * dt;
503 else if (rate_f < 0) {
504 for (
int tIdx = 0; tIdx < tr.numTracer(); ++tIdx) {
505 const Scalar delta = rate_f * wtracer[tIdx];
508 tracerRate[tIdx].rate += delta;
509 freeTracerRate[tIdx].rate += delta;
512 tr.residual_[tIdx][I][Free] -= rate_f * tr.concentration_[tIdx][I][Free];
514 dVol_[Free][tr.phaseIdx_][I] -= rate_f * dt;
517 (*tr.mat)[I][I][Free][Free] -= rate_f * variable<TracerEvaluation>(1.0, 0).derivative(0);
520 for (
int tIdx = 0; tIdx < tr.numTracer(); ++tIdx) {
522 tr.residual_[tIdx][I][Solution] -= rate_s * tr.concentration_[tIdx][I][Solution];
524 dVol_[Solution][tr.phaseIdx_][I] -= rate_s * dt;
527 (*tr.mat)[I][I][Solution][Solution] -= rate_s * variable<TracerEvaluation>(1.0, 0).derivative(0);
533 void assembleTracerEquationSource(TrRe& tr,
537 if (tr.numTracer() == 0) {
542 if (tr.phaseIdx_ == FluidSystem::waterPhaseIdx ||
543 (tr.phaseIdx_ == FluidSystem::gasPhaseIdx && !FluidSystem::enableDissolvedGas()) ||
544 (tr.phaseIdx_ == FluidSystem::oilPhaseIdx && !FluidSystem::enableVaporizedOil()))
549 const Scalar& dsVol = dVol_[Solution][tr.phaseIdx_][I];
550 const Scalar& dfVol = dVol_[Free][tr.phaseIdx_][I];
553 for (
int tIdx = 0; tIdx < tr.numTracer(); ++tIdx) {
555 const auto delta = (dfVol / dt) * tr.concentration_[tIdx][I][Free];
556 tr.residual_[tIdx][I][Free] -= delta;
557 tr.residual_[tIdx][I][Solution] += delta;
560 const auto delta = (dsVol / dt) * tr.concentration_[tIdx][I][Solution];
561 tr.residual_[tIdx][I][Free] += delta;
562 tr.residual_[tIdx][I][Solution] -= delta;
568 const auto delta = (dfVol / dt) * variable<TracerEvaluation>(1.0, 0).derivative(0);
569 (*tr.mat)[I][I][Free][Free] -= delta;
570 (*tr.mat)[I][I][Solution][Free] += delta;
573 const auto delta = (dsVol / dt) * variable<TracerEvaluation>(1.0, 0).derivative(0);
574 (*tr.mat)[I][I][Free][Solution] += delta;
575 (*tr.mat)[I][I][Solution][Solution] -= delta;
579 void assembleTracerEquations_()
587 DeferredLogger local_deferredLogger{};
588 OPM_BEGIN_PARALLEL_TRY_CATCH()
590 OPM_TIMEBLOCK(tracerAssemble);
591 for (
auto& tr : tbatch) {
592 if (tr.numTracer() != 0) {
594 for (
int tIdx = 0; tIdx < tr.numTracer(); ++tIdx) {
595 tr.residual_[tIdx] = 0.0;
600 this->wellTracerRate_.clear();
601 this->wellFreeTracerRate_.clear();
602 this->wellSolTracerRate_.clear();
605 const auto num_msw = this->mSwTracerRate_.size();
606 this->mSwTracerRate_.clear();
609 const auto& wellPtrs = simulator_.problem().wellModel().localNonshutWells();
610 this->wellTracerRate_.reserve(wellPtrs.size());
611 this->wellFreeTracerRate_.reserve(wellPtrs.size());
612 this->wellSolTracerRate_.reserve(wellPtrs.size());
613 this->mSwTracerRate_.reserve(num_msw);
614 for (
const auto& wellPtr : wellPtrs) {
615 for (
auto& tr : tbatch) {
616 this->assembleTracerEquationWell(tr, *wellPtr);
622 #pragma omp parallel for
624 for (
const auto& chunk : element_chunks_) {
625 ElementContext elemCtx(simulator_);
626 const Scalar dt = elemCtx.simulator().timeStepSize();
628 for (
const auto& elem : chunk) {
629 elemCtx.updateStencil(elem);
630 const std::size_t I = elemCtx.globalSpaceIndex( 0, 0);
632 if (elem.partitionType() != Dune::InteriorEntity) {
636 for (
const auto& tr : tbatch) {
637 if (tr.numTracer() != 0) {
638 (*tr.mat)[I][I][0][0] = 1.;
639 (*tr.mat)[I][I][1][1] = 1.;
644 elemCtx.updateAllIntensiveQuantities();
645 elemCtx.updateAllExtensiveQuantities();
647 const Scalar extrusionFactor =
648 elemCtx.intensiveQuantities( 0, 0).extrusionFactor();
649 Valgrind::CheckDefined(extrusionFactor);
650 assert(isfinite(extrusionFactor));
651 assert(extrusionFactor > 0.0);
652 const Scalar scvVolume =
653 elemCtx.stencil(0).subControlVolume( 0).volume() * extrusionFactor;
654 const std::size_t I1 = elemCtx.globalSpaceIndex( 0, 1);
658 for (
auto& tr : tbatch) {
659 if (tr.numTracer() == 0) {
662 this->assembleTracerEquationVolume(tr, elemCtx, scvVolume, dt, I, I1);
665 const std::size_t numInteriorFaces = elemCtx.numInteriorFaces(0);
666 for (
unsigned scvfIdx = 0; scvfIdx < numInteriorFaces; scvfIdx++) {
667 const auto& face = elemCtx.stencil(0).interiorFace(scvfIdx);
668 const unsigned j = face.exteriorIndex();
669 const unsigned J = elemCtx.globalSpaceIndex( j, 0);
670 for (
auto& tr : tbatch) {
671 if (tr.numTracer() == 0) {
674 this->assembleTracerEquationFlux(tr, elemCtx, scvfIdx, I, J, dt);
679 for (
auto& tr : tbatch) {
680 if (tr.numTracer() == 0) {
683 this->assembleTracerEquationSource(tr, dt, I);
688 OPM_END_PARALLEL_TRY_CATCH_LOG(local_deferredLogger,
689 "assembleTracerEquations() failed: ",
690 true, simulator_.gridView().comm())
693 for (auto& tr : tbatch) {
694 if (tr.numTracer() == 0) {
697 auto handle = VectorVectorDataHandle<GridView, std::vector<TracerVector>>(tr.residual_,
698 simulator_.gridView());
699 simulator_.gridView().communicate(handle, Dune::InteriorBorder_All_Interface,
700 Dune::ForwardCommunication);
704 template<TracerTypeIdx Index,
class TrRe>
705 void updateElem(TrRe& tr,
706 const Scalar scvVolume,
707 const unsigned globalDofIdx)
709 const Scalar vol1 = computeVolume_<Index>(tr.phaseIdx_, globalDofIdx, 0);
710 vol1_[
Index][tr.phaseIdx_][globalDofIdx] = vol1 * scvVolume;
711 dVol_[
Index][tr.phaseIdx_][globalDofIdx] = 0.0;
712 for (
int tIdx = 0; tIdx < tr.numTracer(); ++tIdx) {
713 tr.storageOfTimeIndex1_[tIdx][globalDofIdx][
Index] =
714 vol1 * tr.concentrationInitial_[tIdx][globalDofIdx][
Index];
718 void updateStorageCache()
720 for (
auto& tr : tbatch) {
721 if (tr.numTracer() != 0) {
722 tr.concentrationInitial_ = tr.concentration_;
728 #pragma omp parallel for
730 for (
const auto& chunk : element_chunks_) {
731 ElementContext elemCtx(simulator_);
733 for (
const auto& elem : chunk) {
734 elemCtx.updatePrimaryStencil(elem);
735 elemCtx.updatePrimaryIntensiveQuantities(0);
736 const Scalar extrusionFactor = elemCtx.intensiveQuantities( 0, 0).extrusionFactor();
737 const Scalar scvVolume = elemCtx.stencil(0).subControlVolume( 0).volume() * extrusionFactor;
738 const unsigned globalDofIdx = elemCtx.globalSpaceIndex(0, 0);
740 for (
auto& tr : tbatch) {
741 if (tr.numTracer() == 0) {
746 updateElem<Free>(tr, scvVolume, globalDofIdx);
747 updateElem<Solution>(tr, scvVolume, globalDofIdx);
753 template<TracerTypeIdx Index,
class TrRe>
754 void copyForOutput(TrRe& tr,
755 const std::vector<TracerVector>& dx,
758 const unsigned globalDofIdx,
759 std::vector<TracerVectorSingle>& sc)
761 constexpr Scalar tol_gas_sat = 1e-6;
762 tr.concentration_[tIdx][globalDofIdx][
Index] -= dx[tIdx][globalDofIdx][
Index];
763 if (tr.concentration_[tIdx][globalDofIdx][Index] < 0.0 || S < tol_gas_sat) {
764 tr.concentration_[tIdx][globalDofIdx][
Index] = 0.0;
766 sc[tr.idx_[tIdx]][globalDofIdx] = tr.concentration_[tIdx][globalDofIdx][
Index];
769 template<TracerTypeIdx Index,
class TrRe>
770 void assignRates(
const TrRe& tr,
775 std::vector<WellTracerRate<Scalar>>& tracerRate,
776 std::vector<MSWellTracerRate<Scalar>>* mswTracerRate,
777 std::vector<WellTracerRate<Scalar>>& splitRate)
780 for (
int tIdx = 0; tIdx < tr.numTracer(); ++tIdx) {
782 const Scalar delta = rate * tr.concentration_[tIdx][I][
Index];
783 tracerRate[tIdx].rate += delta;
784 splitRate[tIdx].rate += delta;
785 if (eclWell.isMultiSegment()) {
786 (*mswTracerRate)[tIdx].rate[eclWell.getConnections().get(i).segment()] += delta;
792 void advanceTracerFields()
794 assembleTracerEquations_();
796 for (
auto& tr : tbatch) {
797 if (tr.numTracer() == 0) {
803 std::vector<TracerVector> dx(tr.concentration_);
804 for (
int tIdx = 0; tIdx < tr.numTracer(); ++tIdx) {
808 const bool converged = this->linearSolveBatchwise_(*tr.mat, dx, tr.residual_);
810 OpmLog::warning(
"### Tracer model: Linear solver did not converge. ###");
813 OPM_TIMEBLOCK(tracerPost);
815 for (
int tIdx = 0; tIdx < tr.numTracer(); ++tIdx) {
816 for (std::size_t globalDofIdx = 0; globalDofIdx < tr.concentration_[tIdx].size(); ++globalDofIdx) {
819 const auto& intQuants = simulator_.model().intensiveQuantities(globalDofIdx, 0);
820 const auto& fs = intQuants.fluidState();
821 const Scalar Sf = decay<Scalar>(fs.saturation(tr.phaseIdx_));
824 if (tr.phaseIdx_ == FluidSystem::gasPhaseIdx && FluidSystem::enableDissolvedGas()) {
825 Ss = decay<Scalar>(fs.saturation(FluidSystem::oilPhaseIdx));
827 else if (tr.phaseIdx_ == FluidSystem::oilPhaseIdx && FluidSystem::enableVaporizedOil()) {
828 Ss = decay<Scalar>(fs.saturation(FluidSystem::gasPhaseIdx));
831 copyForOutput<Free>(tr, dx, Sf, tIdx, globalDofIdx, this->freeTracerConcentration_);
832 copyForOutput<Solution>(tr, dx, Ss, tIdx, globalDofIdx, this->solTracerConcentration_);
837 const auto& wellPtrs = simulator_.problem().wellModel().localNonshutWells();
838 for (
const auto& wellPtr : wellPtrs) {
839 const auto& eclWell = wellPtr->wellEcl();
842 if (!eclWell.isProducer()) {
846 Scalar rateWellPos = 0.0;
847 Scalar rateWellNeg = 0.0;
848 const std::size_t well_index = simulator_.problem().wellModel().wellState().index(eclWell.name()).value();
849 const auto& ws = simulator_.problem().wellModel().wellState().well(well_index);
850 auto& tracerRate = this->wellTracerRate_[eclWell.seqIndex()];
851 auto& freeTracerRate = this->wellFreeTracerRate_[eclWell.seqIndex()];
852 auto& solTracerRate = this->wellSolTracerRate_[eclWell.seqIndex()];
853 auto* mswTracerRate = eclWell.isMultiSegment() ? &this->mSwTracerRate_[eclWell.seqIndex()] :
nullptr;
855 for (std::size_t i = 0; i < ws.perf_data.size(); ++i) {
856 const auto I = ws.perf_data.cell_index[i];
857 const Scalar rate = wellPtr->volumetricSurfaceRateForConnection(I, tr.phaseIdx_);
860 if (tr.phaseIdx_ == FluidSystem::oilPhaseIdx && FluidSystem::enableVaporizedOil()) {
861 rate_s = ws.perf_data.phase_mixing_rates[i][ws.vaporized_oil];
863 else if (tr.phaseIdx_ == FluidSystem::gasPhaseIdx && FluidSystem::enableDissolvedGas()) {
864 rate_s = ws.perf_data.phase_mixing_rates[i][ws.dissolved_gas];
870 const Scalar rate_f = rate - rate_s;
871 assignRates<Free>(tr, eclWell, i, I, rate_f,
872 tracerRate, mswTracerRate, freeTracerRate);
873 assignRates<Solution>(tr, eclWell, i, I, rate_s,
874 tracerRate, mswTracerRate, solTracerRate);
889 const Scalar official_well_rate_total =
890 simulator_.problem().wellModel().wellState().well(well_index).surface_rates[tr.phaseIdx_];
892 const Scalar rateWellTotal = official_well_rate_total;
894 if (rateWellTotal > rateWellNeg) {
895 constexpr Scalar bucketPrDay = 10.0 / (1000. * 3600. * 24.);
896 const Scalar factor = (rateWellTotal < -bucketPrDay) ? rateWellTotal / rateWellNeg : 0.0;
897 for (
int tIdx = 0; tIdx < tr.numTracer(); ++tIdx) {
898 tracerRate[tIdx].rate *= factor;
905 Simulator& simulator_;
915 template <
typename TV>
918 std::vector<int> idx_;
920 std::vector<TV> concentrationInitial_;
921 std::vector<TV> concentration_;
922 std::vector<TV> storageOfTimeIndex1_;
923 std::vector<TV> residual_;
924 std::unique_ptr<TracerMatrix> mat;
926 bool operator==(
const TracerBatch& rhs)
const
928 return this->concentrationInitial_ == rhs.concentrationInitial_ &&
929 this->concentration_ == rhs.concentration_;
932 static TracerBatch serializationTestObject()
934 TracerBatch<TV> result(4);
935 result.idx_ = {1,2,3};
936 result.concentrationInitial_ = {5.0, 6.0};
937 result.concentration_ = {7.0, 8.0};
938 result.storageOfTimeIndex1_ = {9.0, 10.0, 11.0};
939 result.residual_ = {12.0, 13.0};
944 template<
class Serializer>
945 void serializeOp(Serializer& serializer)
947 serializer(concentrationInitial_);
948 serializer(concentration_);
951 TracerBatch(
int phaseIdx = 0) : phaseIdx_(phaseIdx) {}
953 int numTracer()
const
954 {
return idx_.size(); }
956 void addTracer(
const int idx,
const TV& concentration)
958 const int numGridDof = concentration.size();
959 idx_.emplace_back(idx);
960 concentrationInitial_.emplace_back(concentration);
961 concentration_.emplace_back(concentration);
962 residual_.emplace_back(numGridDof);
963 storageOfTimeIndex1_.emplace_back(numGridDof);
967 std::array<TracerBatch<TracerVector>,numPhases> tbatch;
971 std::array<std::array<std::vector<Scalar>,numPhases>,2> vol1_;
972 std::array<std::array<std::vector<Scalar>,numPhases>,2> dVol_;
973 ElementChunks<GridView, Dune::Partitions::All> element_chunks_;