TracerModel.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_TRACER_MODEL_HPP
29#define OPM_TRACER_MODEL_HPP
30
31#include <opm/common/OpmLog/OpmLog.hpp>
32#include <opm/common/TimingMacros.hpp>
33
34#include <opm/input/eclipse/Schedule/Well/Well.hpp>
35#include <opm/input/eclipse/Schedule/Well/WellConnections.hpp>
36
37#include <opm/grid/utility/ElementChunks.hpp>
38
41
45
46#include <array>
47#include <cstddef>
48#include <memory>
49#include <stdexcept>
50#include <string>
51#include <vector>
52
53namespace Opm::Properties {
54
55template<class TypeTag, class MyTypeTag>
58};
59
60} // namespace Opm::Properties
61
62namespace Opm {
63
69template <class TypeTag>
70class TracerModel : public GenericTracerModel<GetPropType<TypeTag, Properties::Grid>,
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>>
76{
92
93 using TracerEvaluation = DenseAd::Evaluation<Scalar,1>;
94
95 using TracerMatrix = typename BaseType::TracerMatrix;
96 using TracerVector = typename BaseType::TracerVector;
97 using TracerVectorSingle = typename BaseType::TracerVectorSingle;
98
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 };
104
105public:
106 explicit TracerModel(Simulator& simulator)
107 : BaseType(simulator.vanguard().gridView(),
108 simulator.vanguard().eclState(),
109 simulator.vanguard().cartesianIndexMapper(),
110 simulator.model().dofMapper(),
111 simulator.vanguard().cellCentroids())
112 , simulator_(simulator)
113 , tbatch({waterPhaseIdx, oilPhaseIdx, gasPhaseIdx})
114 , wat_(tbatch[0])
115 , oil_(tbatch[1])
116 , gas_(tbatch[2])
117 , element_chunks_(simulator.gridView(), Dune::Partitions::all, ThreadManager::maxThreads())
118 { }
119
120
121 /*
122 The initialization of the tracer model is a three step process:
123
124 1. The init() method is called. This will allocate buffers and initialize
125 some phase index stuff. If this is a normal run the initial tracer
126 concentrations will be assigned from the TBLK or TVDPF keywords.
127
128 2. [Restart only:] The tracer concentration are read from the restart
129 file and the concentrations are applied with repeated calls to the
130 setTracerConcentration() method. This is currently done in the
131 eclwriter::beginRestart() method.
132
133 3. Internally the tracer model manages the concentrations in "batches" for
134 the oil, water and gas tracers respectively. The batches should be
135 initialized with the initial concentration, that must be performed
136 after the concentration values have been assigned. This is done in
137 method prepareTracerBatches() called from eclproblem::finishInit().
138 */
139 void init(bool rst)
140 {
141 this->doInit(rst, simulator_.model().numGridDof(),
142 gasPhaseIdx, oilPhaseIdx, waterPhaseIdx);
143 }
144
146 {
147 for (std::size_t tracerIdx = 0; tracerIdx < this->tracerPhaseIdx_.size(); ++tracerIdx) {
148 if (this->tracerPhaseIdx_[tracerIdx] == FluidSystem::waterPhaseIdx) {
149 if (! FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx)){
150 throw std::runtime_error("Water tracer specified for non-water fluid system: " +
151 this->name(tracerIdx));
152 }
153
154 wat_.addTracer(tracerIdx, this->tracerConcentration_[tracerIdx]);
155 }
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));
160 }
161
162 oil_.addTracer(tracerIdx, this->tracerConcentration_[tracerIdx]);
163 }
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));
168 }
169
170 gas_.addTracer(tracerIdx, this->tracerConcentration_[tracerIdx]);
171 }
172
173 // resize free and solution volume storages
174 vol1_[0][this->tracerPhaseIdx_[tracerIdx]].
175 resize(this->freeTracerConcentration_[tracerIdx].size());
176 vol1_[1][this->tracerPhaseIdx_[tracerIdx]].
177 resize(this->freeTracerConcentration_[tracerIdx].size());
178 dVol_[0][this->tracerPhaseIdx_[tracerIdx]].
179 resize(this->solTracerConcentration_[tracerIdx].size());
180 dVol_[1][this->tracerPhaseIdx_[tracerIdx]].
181 resize(this->solTracerConcentration_[tracerIdx].size());
182 }
183
184 // will be valid after we move out of tracerMatrix_
185 TracerMatrix* base = this->tracerMatrix_.get();
186 for (auto& tr : this->tbatch) {
187 if (tr.numTracer() != 0) {
188 if (this->tracerMatrix_) {
189 tr.mat = std::move(this->tracerMatrix_);
190 }
191 else {
192 tr.mat = std::make_unique<TracerMatrix>(*base);
193 }
194 }
195 }
196 }
197
199 {
200 if (this->numTracers() == 0) {
201 return;
202 }
203
204 OPM_TIMEBLOCK(tracerUpdateCache);
206 }
207
212 {
213 if (this->numTracers() == 0) {
214 return;
215 }
216
217 OPM_TIMEBLOCK(tracerAdvance);
219 }
220
225 template <class Restarter>
226 void serialize(Restarter&)
227 { /* not implemented */ }
228
235 template <class Restarter>
236 void deserialize(Restarter&)
237 { /* not implemented */ }
238
239 template<class Serializer>
240 void serializeOp(Serializer& serializer)
241 {
242 serializer(static_cast<BaseType&>(*this));
243 serializer(tbatch);
244 }
245
246protected:
248 using BaseType::Free;
249 using BaseType::Solution;
250
251 // compute volume associated with free/solution concentration
252 template<TracerTypeIdx Index>
253 Scalar computeVolume_(const int tracerPhaseIdx,
254 const unsigned globalDofIdx,
255 const unsigned timeIdx) const
256 {
257 const auto& intQuants = simulator_.model().intensiveQuantities(globalDofIdx, timeIdx);
258 const auto& fs = intQuants.fluidState();
259 constexpr Scalar min_volume = 1e-10;
260
261 if constexpr (Index == Free) {
262 return std::max(decay<Scalar>(fs.saturation(tracerPhaseIdx)) *
263 decay<Scalar>(fs.invB(tracerPhaseIdx)) *
264 decay<Scalar>(intQuants.porosity()),
265 min_volume);
266 } else {
267 // vaporized oil
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()),
273 min_volume);
274 }
275
276 // dissolved gas
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()),
282 min_volume);
283 }
284
285 return min_volume;
286 }
287 }
288
289 template<TracerTypeIdx Index>
290 std::pair<TracerEvaluation, bool>
291 computeFlux_(const int tracerPhaseIdx,
292 const ElementContext& elemCtx,
293 const unsigned scvfIdx,
294 const unsigned timeIdx) const
295 {
296 const auto& stencil = elemCtx.stencil(timeIdx);
297 const auto& scvf = stencil.interiorFace(scvfIdx);
298
299 const auto& extQuants = elemCtx.extensiveQuantities(scvfIdx, timeIdx);
300 const unsigned inIdx = extQuants.interiorIndex();
301
302 Scalar v;
303 unsigned upIdx;
304
305 if constexpr (Index == Free) {
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));
311 } else {
312 if (tracerPhaseIdx == FluidSystem::oilPhaseIdx && FluidSystem::enableVaporizedOil()) {
313 upIdx = extQuants.upstreamIndex(FluidSystem::gasPhaseIdx);
314
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());
320 }
321 // dissolved gas
322 else if (tracerPhaseIdx == FluidSystem::gasPhaseIdx && FluidSystem::enableDissolvedGas()) {
323 upIdx = extQuants.upstreamIndex(FluidSystem::oilPhaseIdx);
324
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());
330 }
331 else {
332 upIdx = 0;
333 v = 0.0;
334 }
335 }
336
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};
341 }
342
343 template<TracerTypeIdx Index, class TrRe>
344 Scalar storage1_(const TrRe& tr,
345 const unsigned tIdx,
346 const unsigned I,
347 const unsigned I1,
348 const bool cache)
349 {
350 if (cache) {
351 return tr.storageOfTimeIndex1_[tIdx][I][Index];
352 } else {
353 return computeVolume_<Index>(tr.phaseIdx_, I1, 1) *
354 tr.concentration_[tIdx][I1][Index];
355 }
356 }
357
358 template<class TrRe>
360 const ElementContext& elemCtx,
361 const Scalar scvVolume,
362 const Scalar dt,
363 unsigned I,
364 unsigned I1)
365
366 {
367 if (tr.numTracer() == 0) {
368 return;
369 }
370
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) {
376 // Free part
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; // residual + flux
381
382 // Solution part
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; // residual + flux
387 }
388
389 // Derivative matrix
390 (*tr.mat)[I][I][Free][Free] += fVol.derivative(0) * scvVolume/dt;
391 (*tr.mat)[I][I][Solution][Solution] += sVol.derivative(0) * scvVolume/dt;
392 }
393
394 template<class TrRe>
396 const ElementContext& elemCtx,
397 unsigned scvfIdx,
398 unsigned I,
399 unsigned J,
400 const Scalar dt)
401 {
402 if (tr.numTracer() == 0) {
403 return;
404 }
405
406 const auto& [fFlux, isUpF] = computeFlux_<Free>(tr.phaseIdx_, elemCtx, scvfIdx, 0);
407 const auto& [sFlux, isUpS] = computeFlux_<Solution>(tr.phaseIdx_, elemCtx, scvfIdx, 0);
408 dVol_[Solution][tr.phaseIdx_][I] += sFlux.value() * dt;
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) {
413 // Free and solution fluxes
414 tr.residual_[tIdx][I][Free] += fFlux.value()*tr.concentration_[tIdx][fGlobalUpIdx][Free]; // residual + flux
415 tr.residual_[tIdx][I][Solution] += sFlux.value()*tr.concentration_[tIdx][sGlobalUpIdx][Solution]; // residual + flux
416 }
417
418 // Derivative matrix
419 if (isUpF){
420 (*tr.mat)[J][I][Free][Free] = -fFlux.derivative(0);
421 (*tr.mat)[I][I][Free][Free] += fFlux.derivative(0);
422 }
423 if (isUpS) {
424 (*tr.mat)[J][I][Solution][Solution] = -sFlux.derivative(0);
425 (*tr.mat)[I][I][Solution][Solution] += sFlux.derivative(0);
426 }
427 }
428
429 template<class TrRe, class Well>
431 const Well& well)
432 {
433 if (tr.numTracer() == 0) {
434 return;
435 }
436
437 const auto& eclWell = well.wellEcl();
438
439 // Init. well output to zero
440 auto& tracerRate = this->wellTracerRate_[eclWell.seqIndex()];
441 auto& solTracerRate = this->wellSolTracerRate_[eclWell.seqIndex()];
442 auto& freeTracerRate = this->wellFreeTracerRate_[eclWell.seqIndex()];
443 auto* mswTracerRate = eclWell.isMultiSegment()
444 ? &this->mSwTracerRate_[eclWell.seqIndex()]
445 : nullptr;
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);
455 }
456 }
457 }
458
459 std::vector<Scalar> wtracer(tr.numTracer());
460 for (int tIdx = 0; tIdx < tr.numTracer(); ++tIdx) {
461 wtracer[tIdx] = this->currentConcentration_(eclWell, this->name(tr.idx_[tIdx]),
462 simulator_.problem().wellModel().summaryState());
463 }
464
465 const Scalar dt = simulator_.timeStepSize();
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_);
470 Scalar rate_s;
471 if (tr.phaseIdx_ == FluidSystem::oilPhaseIdx && FluidSystem::enableVaporizedOil()) {
472 rate_s = ws.perf_data.phase_mixing_rates[i][ws.vaporized_oil];
473 }
474 else if (tr.phaseIdx_ == FluidSystem::gasPhaseIdx && FluidSystem::enableDissolvedGas()) {
475 rate_s = ws.perf_data.phase_mixing_rates[i][ws.dissolved_gas];
476 }
477 else {
478 rate_s = 0.0;
479 }
480
481 const Scalar rate_f = rate - rate_s;
482 if (rate_f > 0) {
483 for (int tIdx = 0; tIdx < tr.numTracer(); ++tIdx) {
484 const Scalar delta = rate_f * wtracer[tIdx];
485 // Injection of free tracer only
486 tr.residual_[tIdx][I][Free] -= delta;
487
488 // Store _injector_ tracer rate for reporting
489 // (can be done here since WTRACER is constant)
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;
494 }
495 }
496 dVol_[Free][tr.phaseIdx_][I] -= rate_f * dt;
497 }
498 else if (rate_f < 0) {
499 for (int tIdx = 0; tIdx < tr.numTracer(); ++tIdx) {
500 const Scalar delta = rate_f * wtracer[tIdx];
501 // Store _injector_ tracer rate for cross-flowing well connections
502 // (can be done here since WTRACER is constant)
503 tracerRate[tr.idx_[tIdx]].rate += delta;
504 freeTracerRate[tr.idx_[tIdx]].rate += delta;
505
506 // Production of free tracer
507 tr.residual_[tIdx][I][Free] -= rate_f * tr.concentration_[tIdx][I][Free];
508 }
509 dVol_[Free][tr.phaseIdx_][I] -= rate_f * dt;
510
511 // Derivative matrix for free tracer producer
512 (*tr.mat)[I][I][Free][Free] -= rate_f * variable<TracerEvaluation>(1.0, 0).derivative(0);
513 }
514 if (rate_s < 0) {
515 for (int tIdx = 0; tIdx < tr.numTracer(); ++tIdx) {
516 // Production of solution tracer
517 tr.residual_[tIdx][I][Solution] -= rate_s * tr.concentration_[tIdx][I][Solution];
518 }
519 dVol_[Solution][tr.phaseIdx_][I] -= rate_s * dt;
520
521 // Derivative matrix for solution tracer producer
522 (*tr.mat)[I][I][Solution][Solution] -= rate_s * variable<TracerEvaluation>(1.0, 0).derivative(0);
523 }
524 }
525 }
526
527 template<class TrRe>
529 const Scalar dt,
530 unsigned I)
531 {
532 if (tr.numTracer() == 0) {
533 return;
534 }
535
536 // Skip if solution tracers do not exist
537 if (tr.phaseIdx_ == FluidSystem::waterPhaseIdx ||
538 (tr.phaseIdx_ == FluidSystem::gasPhaseIdx && !FluidSystem::enableDissolvedGas()) ||
539 (tr.phaseIdx_ == FluidSystem::oilPhaseIdx && !FluidSystem::enableVaporizedOil()))
540 {
541 return;
542 }
543
544 const Scalar& dsVol = dVol_[Solution][tr.phaseIdx_][I];
545 const Scalar& dfVol = dVol_[Free][tr.phaseIdx_][I];
546
547 // Source term determined by sign of dsVol: if dsVol > 0 then ms -> mf, else mf -> ms
548 for (int tIdx = 0; tIdx < tr.numTracer(); ++tIdx) {
549 if (dsVol >= 0) {
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;
553 }
554 else {
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;
558 }
559 }
560
561 // Derivative matrix
562 if (dsVol >= 0) {
563 const auto delta = (dfVol / dt) * variable<TracerEvaluation>(1.0, 0).derivative(0);
564 (*tr.mat)[I][I][Free][Free] -= delta;
565 (*tr.mat)[I][I][Solution][Free] += delta;
566 }
567 else {
568 const auto delta = (dsVol / dt) * variable<TracerEvaluation>(1.0, 0).derivative(0);
569 (*tr.mat)[I][I][Free][Solution] += delta;
570 (*tr.mat)[I][I][Solution][Solution] -= delta;
571 }
572 }
573
575 {
576 // Note that we formulate the equations in terms of a concentration update
577 // (compared to previous time step) and not absolute concentration.
578 // This implies that current concentration (tr.concentration_[][]) contributes
579 // to the rhs both through storage and flux terms.
580 // Compare also advanceTracerFields(...) below.
581
582 DeferredLogger local_deferredLogger{};
584 {
585 OPM_TIMEBLOCK(tracerAssemble);
586 for (auto& tr : tbatch) {
587 if (tr.numTracer() != 0) {
588 (*tr.mat) = 0.0;
589 for (int tIdx = 0; tIdx < tr.numTracer(); ++tIdx) {
590 tr.residual_[tIdx] = 0.0;
591 }
592 }
593 }
594
595 this->wellTracerRate_.clear();
596 this->wellFreeTracerRate_.clear();
597 this->wellSolTracerRate_.clear();
598
599 // educated guess for new container size
600 const auto num_msw = this->mSwTracerRate_.size();
601 this->mSwTracerRate_.clear();
602
603 // Well terms
604 const auto& wellPtrs = simulator_.problem().wellModel().localNonshutWells();
605 this->wellTracerRate_.reserve(wellPtrs.size());
606 this->wellFreeTracerRate_.reserve(wellPtrs.size());
607 this->wellSolTracerRate_.reserve(wellPtrs.size());
608 this->mSwTracerRate_.reserve(num_msw);
609 for (const auto& wellPtr : wellPtrs) {
610 // Resize vectors of well tracer rates to total number of tracers
611 const auto& eclWell = wellPtr->wellEcl();
612 this->wellTracerRate_[eclWell.seqIndex()].resize(this->numTracers());
613 this->wellFreeTracerRate_[eclWell.seqIndex()].resize(this->numTracers());
614 this->wellSolTracerRate_[eclWell.seqIndex()].resize(this->numTracers());
615 auto* mswTracerRate = eclWell.isMultiSegment()
616 ? &this->mSwTracerRate_[eclWell.seqIndex()]
617 : nullptr;
618 if (mswTracerRate) {
619 mswTracerRate->resize(this->numTracers());
620 }
621 for (auto& tr : tbatch) {
622 this->assembleTracerEquationWell(tr, *wellPtr);
623 }
624 }
625
626 // Parallel loop over element chunks
627 #ifdef _OPENMP
628 #pragma omp parallel for
629 #endif
630 for (const auto& chunk : element_chunks_) {
631 ElementContext elemCtx(simulator_);
632 const Scalar dt = elemCtx.simulator().timeStepSize();
633
634 for (const auto& elem : chunk) {
635 elemCtx.updateStencil(elem);
636 const std::size_t I = elemCtx.globalSpaceIndex(/*dofIdx=*/ 0, /*timeIdx=*/0);
637
638 if (elem.partitionType() != Dune::InteriorEntity) {
639 // Dirichlet boundary conditions for parallel matrix
640 // This is safe as each element has a unique I. So each thread
641 // always writes to different memory locations in the shared arrays.
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.;
646 }
647 }
648 continue;
649 }
650 elemCtx.updateAllIntensiveQuantities();
651 elemCtx.updateAllExtensiveQuantities();
652
653 const Scalar extrusionFactor =
654 elemCtx.intensiveQuantities(/*dofIdx=*/ 0, /*timeIdx=*/0).extrusionFactor();
655 Valgrind::CheckDefined(extrusionFactor);
656 assert(isfinite(extrusionFactor));
657 assert(extrusionFactor > 0.0);
658 const Scalar scvVolume =
659 elemCtx.stencil(/*timeIdx=*/0).subControlVolume(/*dofIdx=*/ 0).volume() * extrusionFactor;
660 const std::size_t I1 = elemCtx.globalSpaceIndex(/*dofIdx=*/ 0, /*timeIdx=*/1);
661
662 // This is safe as each element has a unique I. So each thread
663 // always writes to different memory locations in the shared arrays.
664 for (auto& tr : tbatch) {
665 if (tr.numTracer() == 0) {
666 continue;
667 }
668 this->assembleTracerEquationVolume(tr, elemCtx, scvVolume, dt, I, I1);
669 }
670
671 const std::size_t numInteriorFaces = elemCtx.numInteriorFaces(/*timIdx=*/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(/*dofIdx=*/ j, /*timIdx=*/0);
676 for (auto& tr : tbatch) {
677 if (tr.numTracer() == 0) {
678 continue;
679 }
680 this->assembleTracerEquationFlux(tr, elemCtx, scvfIdx, I, J, dt);
681 }
682 }
683
684 // Source terms (mass transfer between free and solution tracer)
685 for (auto& tr : tbatch) {
686 if (tr.numTracer() == 0) {
687 continue;
688 }
689 this->assembleTracerEquationSource(tr, dt, I);
690 }
691 }
692 }
693 }
694 OPM_END_PARALLEL_TRY_CATCH_LOG(local_deferredLogger,
695 "assembleTracerEquations() failed: ",
696 true, simulator_.gridView().comm())
697
698 // Communicate overlap using grid Communication
699 for (auto& tr : tbatch) {
700 if (tr.numTracer() == 0) {
701 continue;
702 }
704 simulator_.gridView());
705 simulator_.gridView().communicate(handle, Dune::InteriorBorder_All_Interface,
706 Dune::ForwardCommunication);
707 }
708 }
709
710 template<TracerTypeIdx Index, class TrRe>
711 void updateElem(TrRe& tr,
712 const Scalar scvVolume,
713 const unsigned globalDofIdx)
714 {
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];
721 }
722 }
723
725 {
726 for (auto& tr : tbatch) {
727 if (tr.numTracer() != 0) {
728 tr.concentrationInitial_ = tr.concentration_;
729 }
730 }
731
732 // Parallel loop over element chunks
733 #ifdef _OPENMP
734 #pragma omp parallel for
735 #endif
736 for (const auto& chunk : element_chunks_) {
737 ElementContext elemCtx(simulator_);
738
739 for (const auto& elem : chunk) {
740 elemCtx.updatePrimaryStencil(elem);
741 elemCtx.updatePrimaryIntensiveQuantities(/*timeIdx=*/0);
742 const Scalar extrusionFactor = elemCtx.intensiveQuantities(/*dofIdx=*/ 0, /*timeIdx=*/0).extrusionFactor();
743 const Scalar scvVolume = elemCtx.stencil(/*timeIdx=*/0).subControlVolume(/*dofIdx=*/ 0).volume() * extrusionFactor;
744 const unsigned globalDofIdx = elemCtx.globalSpaceIndex(0, /*timeIdx=*/0);
745
746 for (auto& tr : tbatch) {
747 if (tr.numTracer() == 0) {
748 continue;
749 }
750 // This is safe as each element has a unique globalDofIdx. So each thread
751 // always writes to different memory locations in the shared arrays.
752 updateElem<Free>(tr, scvVolume, globalDofIdx);
753 updateElem<Solution>(tr, scvVolume, globalDofIdx);
754 }
755 }
756 }
757 }
758
759 template<TracerTypeIdx Index, class TrRe>
760 void copyForOutput(TrRe& tr,
761 const std::vector<TracerVector>& dx,
762 const Scalar S,
763 const unsigned tIdx,
764 const unsigned globalDofIdx,
765 std::vector<TracerVectorSingle>& sc)
766 {
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;
771 }
772 sc[tr.idx_[tIdx]][globalDofIdx] = tr.concentration_[tIdx][globalDofIdx][Index];
773 }
774
775 template<TracerTypeIdx Index, class TrRe>
776 void assignRates(const TrRe& tr,
777 const Well& eclWell,
778 const std::size_t i,
779 const std::size_t I,
780 const Scalar rate,
781 std::vector<WellTracerRate<Scalar>>& tracerRate,
782 std::vector<MSWellTracerRate<Scalar>>* mswTracerRate,
783 std::vector<WellTracerRate<Scalar>>& splitRate)
784 {
785 if (rate < 0) {
786 for (int tIdx = 0; tIdx < tr.numTracer(); ++tIdx) {
787 // Store _producer_ free tracer rate for reporting
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;
793 }
794 }
795 }
796 }
797
799 {
801
802 for (auto& tr : tbatch) {
803 if (tr.numTracer() == 0) {
804 continue;
805 }
806
807 // Note that we solve for a concentration update (compared to previous time step)
808 // Confer also assembleTracerEquations_(...) above.
809 std::vector<TracerVector> dx(tr.concentration_);
810 for (int tIdx = 0; tIdx < tr.numTracer(); ++tIdx) {
811 dx[tIdx] = 0.0;
812 }
813
814 const bool converged = this->linearSolveBatchwise_(*tr.mat, dx, tr.residual_);
815 if (!converged) {
816 OpmLog::warning("### Tracer model: Linear solver did not converge. ###");
817 }
818
819 OPM_TIMEBLOCK(tracerPost);
820
821 for (int tIdx = 0; tIdx < tr.numTracer(); ++tIdx) {
822 for (std::size_t globalDofIdx = 0; globalDofIdx < tr.concentration_[tIdx].size(); ++globalDofIdx) {
823 // New concetration. Concentrations that are negative or where free/solution phase is not
824 // present are set to zero
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_));
828 Scalar Ss = 0.0;
829
830 if (tr.phaseIdx_ == FluidSystem::gasPhaseIdx && FluidSystem::enableDissolvedGas()) {
831 Ss = decay<Scalar>(fs.saturation(FluidSystem::oilPhaseIdx));
832 }
833 else if (tr.phaseIdx_ == FluidSystem::oilPhaseIdx && FluidSystem::enableVaporizedOil()) {
834 Ss = decay<Scalar>(fs.saturation(FluidSystem::gasPhaseIdx));
835 }
836
837 copyForOutput<Free>(tr, dx, Sf, tIdx, globalDofIdx, this->freeTracerConcentration_);
838 copyForOutput<Solution>(tr, dx, Ss, tIdx, globalDofIdx, this->solTracerConcentration_);
839 }
840 }
841
842 // Store _producer_ tracer rate for reporting
843 const auto& wellPtrs = simulator_.problem().wellModel().localNonshutWells();
844 for (const auto& wellPtr : wellPtrs) {
845 const auto& eclWell = wellPtr->wellEcl();
846
847 // Injection rates already reported during assembly
848 if (!eclWell.isProducer()) {
849 continue;
850 }
851
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);
856 auto& tracerRate = this->wellTracerRate_[eclWell.seqIndex()];
857 auto& freeTracerRate = this->wellFreeTracerRate_[eclWell.seqIndex()];
858 auto& solTracerRate = this->wellSolTracerRate_[eclWell.seqIndex()];
859 auto* mswTracerRate = eclWell.isMultiSegment() ? &this->mSwTracerRate_[eclWell.seqIndex()] : nullptr;
860
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_);
864
865 Scalar rate_s;
866 if (tr.phaseIdx_ == FluidSystem::oilPhaseIdx && FluidSystem::enableVaporizedOil()) {
867 rate_s = ws.perf_data.phase_mixing_rates[i][ws.vaporized_oil];
868 }
869 else if (tr.phaseIdx_ == FluidSystem::gasPhaseIdx && FluidSystem::enableDissolvedGas()) {
870 rate_s = ws.perf_data.phase_mixing_rates[i][ws.dissolved_gas];
871 }
872 else {
873 rate_s = 0.0;
874 }
875
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);
881
882 if (rate < 0) {
883 rateWellNeg += rate;
884 }
885 else {
886 rateWellPos += rate;
887 }
888 }
889
890 // TODO: Some inconsistencies here that perhaps should be clarified.
891 // The "offical" rate as reported below is occasionally significant
892 // different from the sum over connections (as calculated above). Only observed
893 // for small values, neglible for the rate itself, but matters when used to
894 // calculate tracer concentrations.
895 const Scalar official_well_rate_total =
896 simulator_.problem().wellModel().wellState().well(well_index).surface_rates[tr.phaseIdx_];
897
898 const Scalar rateWellTotal = official_well_rate_total;
899
900 if (rateWellTotal > rateWellNeg) { // Cross flow
901 constexpr Scalar bucketPrDay = 10.0 / (1000. * 3600. * 24.); // ... keeps (some) trouble away
902 const Scalar factor = (rateWellTotal < -bucketPrDay) ? rateWellTotal / rateWellNeg : 0.0;
903 for (int tIdx = 0; tIdx < tr.numTracer(); ++tIdx) {
904 tracerRate[tIdx].rate *= factor;
905 }
906 }
907 }
908 }
909 }
910
911 Simulator& simulator_;
912
913 // This struct collects tracers of the same type (i.e, transported in same phase).
914 // The idea being that, under the assumption of linearity, tracers of same type can
915 // be solved in concert, having a common system matrix but separate right-hand-sides.
916
917 // Since oil or gas tracers appears in dual compositions when VAPOIL respectively DISGAS
918 // is active, the template argument is intended to support future extension to these
919 // scenarios by supplying an extended vector type.
920
921 template <typename TV>
923 {
924 std::vector<int> idx_;
925 const int phaseIdx_;
926 std::vector<TV> concentrationInitial_;
927 std::vector<TV> concentration_;
928 std::vector<TV> storageOfTimeIndex1_;
929 std::vector<TV> residual_;
930 std::unique_ptr<TracerMatrix> mat;
931
932 bool operator==(const TracerBatch& rhs) const
933 {
934 return this->concentrationInitial_ == rhs.concentrationInitial_ &&
935 this->concentration_ == rhs.concentration_;
936 }
937
939 {
940 TracerBatch<TV> result(4);
941 result.idx_ = {1,2,3};
942 result.concentrationInitial_ = {5.0, 6.0};
943 result.concentration_ = {7.0, 8.0};
944 result.storageOfTimeIndex1_ = {9.0, 10.0, 11.0};
945 result.residual_ = {12.0, 13.0};
946
947 return result;
948 }
949
950 template<class Serializer>
951 void serializeOp(Serializer& serializer)
952 {
953 serializer(concentrationInitial_);
954 serializer(concentration_);
955 }
956
957 TracerBatch(int phaseIdx = 0) : phaseIdx_(phaseIdx) {}
958
959 int numTracer() const
960 { return idx_.size(); }
961
962 void addTracer(const int idx, const TV& concentration)
963 {
964 const int numGridDof = concentration.size();
965 idx_.emplace_back(idx);
966 concentrationInitial_.emplace_back(concentration);
967 concentration_.emplace_back(concentration);
968 residual_.emplace_back(numGridDof);
969 storageOfTimeIndex1_.emplace_back(numGridDof);
970 }
971 };
972
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_;
979 ElementChunks<GridView, Dune::Partitions::All> element_chunks_;
980};
981
982} // namespace Opm
983
984#endif // OPM_TRACER_MODEL_HPP
#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
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
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