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
33#include <opm/models/utils/propertysystem.hh>
34
37
38#include <array>
39#include <cstddef>
40#include <memory>
41#include <stdexcept>
42#include <string>
43#include <vector>
44
45namespace Opm::Properties {
46
47template<class TypeTag, class MyTypeTag>
49 using type = UndefinedProperty;
50};
51
52} // namespace Opm::Properties
53
54namespace Opm {
55
61template <class TypeTag>
62class TracerModel : public GenericTracerModel<GetPropType<TypeTag, Properties::Grid>,
63 GetPropType<TypeTag, Properties::GridView>,
64 GetPropType<TypeTag, Properties::DofMapper>,
65 GetPropType<TypeTag, Properties::Stencil>,
66 GetPropType<TypeTag, Properties::Scalar>>
67{
69 GetPropType<TypeTag, Properties::GridView>,
70 GetPropType<TypeTag, Properties::DofMapper>,
71 GetPropType<TypeTag, Properties::Stencil>,
72 GetPropType<TypeTag, Properties::Scalar>>;
73 using Simulator = GetPropType<TypeTag, Properties::Simulator>;
74 using GridView = GetPropType<TypeTag, Properties::GridView>;
75 using Grid = GetPropType<TypeTag, Properties::Grid>;
76 using Scalar = GetPropType<TypeTag, Properties::Scalar>;
77 using Stencil = GetPropType<TypeTag, Properties::Stencil>;
78 using FluidSystem = GetPropType<TypeTag, Properties::FluidSystem>;
79 using ElementContext = GetPropType<TypeTag, Properties::ElementContext>;
80 using RateVector = GetPropType<TypeTag, Properties::RateVector>;
81 using Indices = GetPropType<TypeTag, Properties::Indices>;
82
83 using TracerEvaluation = DenseAd::Evaluation<Scalar,1>;
84
85 using TracerMatrix = typename BaseType::TracerMatrix;
86 using TracerVector = typename BaseType::TracerVector;
87
88 enum { numEq = getPropValue<TypeTag, Properties::NumEq>() };
89 enum { numPhases = FluidSystem::numPhases };
90 enum { waterPhaseIdx = FluidSystem::waterPhaseIdx };
91 enum { oilPhaseIdx = FluidSystem::oilPhaseIdx };
92 enum { gasPhaseIdx = FluidSystem::gasPhaseIdx };
93
94public:
95 TracerModel(Simulator& simulator)
96 : BaseType(simulator.vanguard().gridView(),
97 simulator.vanguard().eclState(),
98 simulator.vanguard().cartesianIndexMapper(),
99 simulator.model().dofMapper(),
100 simulator.vanguard().cellCentroids())
101 , simulator_(simulator)
102 , tbatch({waterPhaseIdx, oilPhaseIdx, gasPhaseIdx})
103 , wat_(tbatch[0])
104 , oil_(tbatch[1])
105 , gas_(tbatch[2])
106 { }
107
108
109 /*
110 The initialization of the tracer model is a three step process:
111
112 1. The init() method is called. This will allocate buffers and initialize
113 some phase index stuff. If this is a normal run the initial tracer
114 concentrations will be assigned from the TBLK or TVDPF keywords.
115
116 2. [Restart only:] The tracer concenntration are read from the restart
117 file and the concentrations are applied with repeated calls to the
118 setTracerConcentration() method. This is currently done in the
119 eclwriter::beginRestart() method.
120
121 3. Internally the tracer model manages the concentrations in "batches" for
122 the oil, water and gas tracers respectively. The batches should be
123 initialized with the initial concentration, that must be performed
124 after the concentration values have been assigned. This is done in
125 method prepareTracerBatches() called from eclproblem::finishInit().
126 */
127 void init(bool rst)
128 {
129 this->doInit(rst, simulator_.model().numGridDof(),
130 gasPhaseIdx, oilPhaseIdx, waterPhaseIdx);
131 }
132
134 {
135 for (std::size_t tracerIdx = 0; tracerIdx < this->tracerPhaseIdx_.size(); ++tracerIdx) {
136 if (this->tracerPhaseIdx_[tracerIdx] == FluidSystem::waterPhaseIdx) {
137 if (! FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx)){
138 throw std::runtime_error("Water tracer specified for non-water fluid system:" + this->name(tracerIdx));
139 }
140
141 wat_.addTracer(tracerIdx, this->tracerConcentration_[tracerIdx]);
142 }
143 else if (this->tracerPhaseIdx_[tracerIdx] == FluidSystem::oilPhaseIdx) {
144 if (! FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx)){
145 throw std::runtime_error("Oil tracer specified for non-oil fluid system:" + this->name(tracerIdx));
146 }
147 if (FluidSystem::enableVaporizedOil()) {
148 throw std::runtime_error("Oil tracer in combination with kw VAPOIL is not supported: " + this->name(tracerIdx));
149 }
150
151 oil_.addTracer(tracerIdx, this->tracerConcentration_[tracerIdx]);
152 }
153 else if (this->tracerPhaseIdx_[tracerIdx] == FluidSystem::gasPhaseIdx) {
154 if (! FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx)){
155 throw std::runtime_error("Gas tracer specified for non-gas fluid system:" + this->name(tracerIdx));
156 }
157 if (FluidSystem::enableDissolvedGas()) {
158 throw std::runtime_error("Gas tracer in combination with kw DISGAS is not supported: " + this->name(tracerIdx));
159 }
160
161 gas_.addTracer(tracerIdx, this->tracerConcentration_[tracerIdx]);
162 }
163 }
164
165 // will be valid after we move out of tracerMatrix_
166 TracerMatrix* base = this->tracerMatrix_.get();
167 for (auto& tr : this->tbatch) {
168 if (tr.numTracer() != 0) {
169 if (this->tracerMatrix_)
170 tr.mat = std::move(this->tracerMatrix_);
171 else
172 tr.mat = std::make_unique<TracerMatrix>(*base);
173 }
174 }
175 }
176
178 {
179 if (this->numTracers() == 0)
180 return;
181
183 }
184
189 {
190 if (this->numTracers() == 0)
191 return;
192
194 }
195
200 template <class Restarter>
201 void serialize(Restarter&)
202 { /* not implemented */ }
203
210 template <class Restarter>
211 void deserialize(Restarter&)
212 { /* not implemented */ }
213
214 template<class Serializer>
215 void serializeOp(Serializer& serializer)
216 {
217 serializer(static_cast<BaseType&>(*this));
218 serializer(tbatch);
219 }
220
221protected:
222
223 // evaluate water storage volume(s) in a single cell
224 template <class LhsEval>
225 void computeVolume_(LhsEval& freeVolume,
226 const int tracerPhaseIdx,
227 const ElementContext& elemCtx,
228 unsigned scvIdx,
229 unsigned timeIdx)
230 {
231 const auto& intQuants = elemCtx.intensiveQuantities(scvIdx, timeIdx);
232 const auto& fs = intQuants.fluidState();
233 Scalar phaseVolume =
234 decay<Scalar>(fs.saturation(tracerPhaseIdx))
235 *decay<Scalar>(fs.invB(tracerPhaseIdx))
236 *decay<Scalar>(intQuants.porosity());
237
238 // avoid singular matrix if no water is present.
239 phaseVolume = max(phaseVolume, 1e-10);
240
241 if (std::is_same<LhsEval, Scalar>::value)
242 freeVolume = phaseVolume;
243 else
244 freeVolume = phaseVolume * variable<LhsEval>(1.0, 0);
245
246 }
247
248 // evaluate the flux(es) over one face
249 void computeFlux_(TracerEvaluation & freeFlux,
250 bool & isUpFree,
251 const int tracerPhaseIdx,
252 const ElementContext& elemCtx,
253 unsigned scvfIdx,
254 unsigned timeIdx)
255
256 {
257 const auto& stencil = elemCtx.stencil(timeIdx);
258 const auto& scvf = stencil.interiorFace(scvfIdx);
259
260 const auto& extQuants = elemCtx.extensiveQuantities(scvfIdx, timeIdx);
261 unsigned inIdx = extQuants.interiorIndex();
262
263 unsigned upIdx = extQuants.upstreamIndex(tracerPhaseIdx);
264
265 const auto& intQuants = elemCtx.intensiveQuantities(upIdx, timeIdx);
266 const auto& fs = intQuants.fluidState();
267
268 Scalar A = scvf.area();
269 Scalar v = decay<Scalar>(extQuants.volumeFlux(tracerPhaseIdx));
270 Scalar b = decay<Scalar>(fs.invB(tracerPhaseIdx));
271
272 if (inIdx == upIdx) {
273 freeFlux = A*v*b*variable<TracerEvaluation>(1.0, 0);
274 isUpFree = true;
275 }
276 else {
277 freeFlux = A*v*b*1.0;
278 isUpFree = false;
279 }
280 }
281
282 template<class TrRe>
284 const ElementContext& elemCtx,
285 const Scalar scvVolume,
286 const Scalar dt,
287 unsigned I,
288 unsigned I1)
289
290 {
291 if (tr.numTracer() == 0)
292 return;
293
294 std::vector<Scalar> storageOfTimeIndex1(tr.numTracer());
295 if (elemCtx.enableStorageCache()) {
296 for (int tIdx = 0; tIdx < tr.numTracer(); ++tIdx) {
297 storageOfTimeIndex1[tIdx] = tr.storageOfTimeIndex1_[tIdx][I];
298 }
299 }
300 else {
301 Scalar fVolume1;
302 computeVolume_(fVolume1, tr.phaseIdx_, elemCtx, 0, /*timeIdx=*/1);
303 for (int tIdx = 0; tIdx < tr.numTracer(); ++tIdx) {
304 storageOfTimeIndex1[tIdx] = fVolume1*tr.concentrationInitial_[tIdx][I1];
305 }
306 }
307
308 TracerEvaluation fVolume;
309 computeVolume_(fVolume, tr.phaseIdx_, elemCtx, 0, /*timeIdx=*/0);
310 for (int tIdx = 0; tIdx < tr.numTracer(); ++tIdx) {
311 Scalar storageOfTimeIndex0 = fVolume.value()*tr.concentration_[tIdx][I];
312 Scalar localStorage = (storageOfTimeIndex0 - storageOfTimeIndex1[tIdx]) * scvVolume/dt;
313 tr.residual_[tIdx][I][0] += localStorage; //residual + flux
314 }
315 (*tr.mat)[I][I][0][0] += fVolume.derivative(0) * scvVolume/dt;
316 }
317
318 template<class TrRe>
320 const ElementContext& elemCtx,
321 unsigned scvfIdx,
322 unsigned I,
323 unsigned J)
324 {
325 if (tr.numTracer() == 0)
326 return;
327
328 TracerEvaluation flux;
329 bool isUpF;
330 computeFlux_(flux, isUpF, tr.phaseIdx_, elemCtx, scvfIdx, 0);
331 int globalUpIdx = isUpF ? I : J;
332 for (int tIdx =0; tIdx < tr.numTracer(); ++tIdx) {
333 tr.residual_[tIdx][I][0] += flux.value()*tr.concentration_[tIdx][globalUpIdx]; //residual + flux
334 }
335 if (isUpF) {
336 (*tr.mat)[J][I][0][0] = -flux.derivative(0);
337 (*tr.mat)[I][I][0][0] += flux.derivative(0);
338 }
339 }
340
341 template<class TrRe, class Well>
343 const Well& well)
344 {
345 if (tr.numTracer() == 0)
346 return;
347
348 const auto& eclWell = well.wellEcl();
349 for (int tIdx = 0; tIdx < tr.numTracer(); ++tIdx) {
350 this->wellTracerRate_[std::make_pair(eclWell.name(), this->name(tr.idx_[tIdx]))] = 0.0;
351 }
352
353 std::vector<double> wtracer(tr.numTracer());
354 for (int tIdx = 0; tIdx < tr.numTracer(); ++tIdx) {
355 wtracer[tIdx] = this->currentConcentration_(eclWell, this->name(tr.idx_[tIdx]));
356 }
357
358 for (auto& perfData : well.perforationData()) {
359 const auto I = perfData.cell_index;
360 Scalar rate = well.volumetricSurfaceRateForConnection(I, tr.phaseIdx_);
361 if (rate > 0) {
362 for (int tIdx = 0; tIdx < tr.numTracer(); ++tIdx) {
363 tr.residual_[tIdx][I][0] -= rate*wtracer[tIdx];
364 // Store _injector_ tracer rate for reporting
365 this->wellTracerRate_.at(std::make_pair(eclWell.name(),this->name(tr.idx_[tIdx]))) += rate*wtracer[tIdx];
366 }
367 }
368 else if (rate < 0) {
369 for (int tIdx = 0; tIdx < tr.numTracer(); ++tIdx) {
370 tr.residual_[tIdx][I][0] -= rate*tr.concentration_[tIdx][I];
371 }
372 (*tr.mat)[I][I][0][0] -= rate*variable<TracerEvaluation>(1.0, 0).derivative(0);
373 }
374 }
375 }
376
378 {
379 // Note that we formulate the equations in terms of a concentration update
380 // (compared to previous time step) and not absolute concentration.
381 // This implies that current concentration (tr.concentration_[][]) contributes
382 // to the rhs both through storrage and flux terms.
383 // Compare also advanceTracerFields(...) below.
384
385 for (auto& tr : tbatch) {
386 if (tr.numTracer() != 0) {
387 (*tr.mat) = 0.0;
388 for (int tIdx = 0; tIdx < tr.numTracer(); ++tIdx)
389 tr.residual_[tIdx] = 0.0;
390 }
391 }
392
393 ElementContext elemCtx(simulator_);
394 for (const auto& elem : elements(simulator_.gridView())) {
395 elemCtx.updateStencil(elem);
396
397 std::size_t I = elemCtx.globalSpaceIndex(/*dofIdx=*/ 0, /*timeIdx=*/0);
398
399 if (elem.partitionType() != Dune::InteriorEntity)
400 {
401 // Dirichlet boundary conditions needed for the parallel matrix
402 for (auto& tr : tbatch) {
403 if (tr.numTracer() != 0)
404 (*tr.mat)[I][I][0][0] = 1.;
405 }
406 continue;
407 }
408 elemCtx.updateAllIntensiveQuantities();
409 elemCtx.updateAllExtensiveQuantities();
410
411 Scalar extrusionFactor =
412 elemCtx.intensiveQuantities(/*dofIdx=*/ 0, /*timeIdx=*/0).extrusionFactor();
413 Valgrind::CheckDefined(extrusionFactor);
414 assert(isfinite(extrusionFactor));
415 assert(extrusionFactor > 0.0);
416 Scalar scvVolume =
417 elemCtx.stencil(/*timeIdx=*/0).subControlVolume(/*dofIdx=*/ 0).volume()
418 * extrusionFactor;
419 Scalar dt = elemCtx.simulator().timeStepSize();
420
421 std::size_t I1 = elemCtx.globalSpaceIndex(/*dofIdx=*/ 0, /*timeIdx=*/1);
422
423 for (auto& tr : tbatch) {
424 this->assembleTracerEquationVolume(tr, elemCtx, scvVolume, dt, I, I1);
425 }
426
427 std::size_t numInteriorFaces = elemCtx.numInteriorFaces(/*timIdx=*/0);
428 for (unsigned scvfIdx = 0; scvfIdx < numInteriorFaces; scvfIdx++) {
429 const auto& face = elemCtx.stencil(0).interiorFace(scvfIdx);
430 unsigned j = face.exteriorIndex();
431 unsigned J = elemCtx.globalSpaceIndex(/*dofIdx=*/ j, /*timIdx=*/0);
432 for (auto& tr : tbatch) {
433 this->assembleTracerEquationFlux(tr, elemCtx, scvfIdx, I, J);
434 }
435 }
436 }
437
438 // Well terms
439 const auto& wellPtrs = simulator_.problem().wellModel().localNonshutWells();
440 for (const auto& wellPtr : wellPtrs) {
441 for (auto& tr : tbatch) {
442 this->assembleTracerEquationWell(tr, *wellPtr);
443 }
444 }
445
446 // Communicate overlap using grid Communication
447 for (auto& tr : tbatch) {
448 if (tr.numTracer() == 0)
449 continue;
451 simulator_.gridView());
452 simulator_.gridView().communicate(handle, Dune::InteriorBorder_All_Interface,
453 Dune::ForwardCommunication);
454 }
455 }
456
458 {
459 for (auto& tr : tbatch) {
460 if (tr.numTracer() != 0)
461 tr.concentrationInitial_ = tr.concentration_;
462 }
463
464 ElementContext elemCtx(simulator_);
465 for (const auto& elem : elements(simulator_.gridView())) {
466 elemCtx.updatePrimaryStencil(elem);
467 elemCtx.updatePrimaryIntensiveQuantities(/*timeIdx=*/0);
468 int globalDofIdx = elemCtx.globalSpaceIndex(0, /*timeIdx=*/0);
469 for (auto& tr : tbatch) {
470 if (tr.numTracer() == 0)
471 continue;
472 Scalar fVolume;
473 computeVolume_(fVolume, tr.phaseIdx_, elemCtx, 0, /*timeIdx=*/0);
474 for (int tIdx = 0; tIdx < tr.numTracer(); ++tIdx) {
475 tr.storageOfTimeIndex1_[tIdx][globalDofIdx] = fVolume*tr.concentrationInitial_[tIdx][globalDofIdx];
476 }
477 }
478 }
479 }
480
482 {
484
485 for (auto& tr : tbatch) {
486 if (tr.numTracer() == 0)
487 continue;
488
489 // Note that we solve for a concentration update (compared to previous time step)
490 // Confer also assembleTracerEquations_(...) above.
491 std::vector<TracerVector> dx(tr.concentration_);
492 for (int tIdx = 0; tIdx < tr.numTracer(); ++tIdx)
493 dx[tIdx] = 0.0;
494
495 bool converged = this->linearSolveBatchwise_(*tr.mat, dx, tr.residual_);
496 if (!converged) {
497 OpmLog::warning("### Tracer model: Linear solver did not converge. ###");
498 }
499
500 for (int tIdx = 0; tIdx < tr.numTracer(); ++tIdx) {
501 tr.concentration_[tIdx] -= dx[tIdx];
502 // Tracer concentrations for restart report
503 this->tracerConcentration_[tr.idx_[tIdx]] = tr.concentration_[tIdx];
504 }
505
506 // Store _producer_ tracer rate for reporting
507 const auto& wellPtrs = simulator_.problem().wellModel().localNonshutWells();
508 for (const auto& wellPtr : wellPtrs) {
509 const auto& well = wellPtr->wellEcl();
510
511 if (!well.isProducer()) //Injection rates already reported during assembly
512 continue;
513
514 Scalar rateWellPos = 0.0;
515 Scalar rateWellNeg = 0.0;
516 for (auto& perfData : wellPtr->perforationData()) {
517 const int I = perfData.cell_index;
518 Scalar rate = wellPtr->volumetricSurfaceRateForConnection(I, tr.phaseIdx_);
519 if (rate < 0) {
520 rateWellNeg += rate;
521 for (int tIdx = 0; tIdx < tr.numTracer(); ++tIdx) {
522 this->wellTracerRate_.at(std::make_pair(well.name(),this->name(tr.idx_[tIdx]))) += rate*tr.concentration_[tIdx][I];
523 }
524 }
525 else {
526 rateWellPos += rate;
527 }
528 }
529
530 Scalar rateWellTotal = rateWellNeg + rateWellPos;
531
532 // TODO: Some inconsistencies here that perhaps should be clarified. The "offical" rate as reported below is
533 // occasionally significant different from the sum over connections (as calculated above). Only observed
534 // for small values, neglible for the rate itself, but matters when used to calculate tracer concentrations.
535 std::size_t well_index = simulator_.problem().wellModel().wellState().index(well.name()).value();
536 Scalar official_well_rate_total = simulator_.problem().wellModel().wellState().well(well_index).surface_rates[tr.phaseIdx_];
537
538 rateWellTotal = official_well_rate_total;
539
540 if (rateWellTotal > rateWellNeg) { // Cross flow
541 const Scalar bucketPrDay = 10.0/(1000.*3600.*24.); // ... keeps (some) trouble away
542 const Scalar factor = (rateWellTotal < -bucketPrDay) ? rateWellTotal/rateWellNeg : 0.0;
543 for (int tIdx = 0; tIdx < tr.numTracer(); ++tIdx) {
544 this->wellTracerRate_.at(std::make_pair(well.name(),this->name(tr.idx_[tIdx]))) *= factor;
545 }
546 }
547 }
548 }
549 }
550
551
552 Simulator& simulator_;
553
554 // This struct collects tracers of the same type (i.e, transported in same phase).
555 // The idea being that, under the assumption of linearity, tracers of same type can
556 // be solved in concert, having a common system matrix but separate right-hand-sides.
557
558 // Since oil or gas tracers appears in dual compositions when VAPOIL respectively DISGAS
559 // is active, the template argument is intended to support future extension to these
560 // scenarios by supplying an extended vector type.
561
562 template <typename TV>
563 struct TracerBatch {
564 std::vector<int> idx_;
565 const int phaseIdx_;
566 std::vector<TV> concentrationInitial_;
567 std::vector<TV> concentration_;
568 std::vector<TV> storageOfTimeIndex1_;
569 std::vector<TV> residual_;
570 std::unique_ptr<TracerMatrix> mat;
571
572 bool operator==(const TracerBatch& rhs) const
573 {
574 return this->concentrationInitial_ == rhs.concentrationInitial_ &&
575 this->concentration_ == rhs.concentration_;
576 }
577
579 {
580 TracerBatch<TV> result(4);
581 result.idx_ = {1,2,3};
582 result.concentrationInitial_ = {5.0, 6.0};
583 result.concentration_ = {7.0, 8.0};
584 result.storageOfTimeIndex1_ = {9.0, 10.0, 11.0};
585 result.residual_ = {12.0, 13.0};
586
587 return result;
588 }
589
590 template<class Serializer>
591 void serializeOp(Serializer& serializer)
592 {
593 serializer(concentrationInitial_);
594 serializer(concentration_);
595 }
596
597 TracerBatch(int phaseIdx = 0) : phaseIdx_(phaseIdx) {}
598
599 int numTracer() const {return idx_.size(); }
600
601 void addTracer(const int idx, const TV & concentration)
602 {
603 int numGridDof = concentration.size();
604 idx_.emplace_back(idx);
605 concentrationInitial_.emplace_back(concentration);
606 concentration_.emplace_back(concentration);
607 storageOfTimeIndex1_.emplace_back(numGridDof);
608 residual_.emplace_back(numGridDof);
609 }
610 };
611
612 std::array<TracerBatch<TracerVector>,3> tbatch;
616};
617
618} // namespace Opm
619
620#endif // OPM_TRACER_MODEL_HPP
A datahandle sending data located in multiple vectors.
Definition: GenericTracerModel.hpp:53
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:162
A class which handles tracers as specified in by ECL.
Definition: TracerModel.hpp:67
std::array< TracerBatch< TracerVector >, 3 > tbatch
Definition: TracerModel.hpp:612
void updateStorageCache()
Definition: TracerModel.hpp:457
TracerBatch< TracerVector > & oil_
Definition: TracerModel.hpp:614
void advanceTracerFields()
Definition: TracerModel.hpp:481
void init(bool rst)
Definition: TracerModel.hpp:127
TracerModel(Simulator &simulator)
Definition: TracerModel.hpp:95
void computeFlux_(TracerEvaluation &freeFlux, bool &isUpFree, const int tracerPhaseIdx, const ElementContext &elemCtx, unsigned scvfIdx, unsigned timeIdx)
Definition: TracerModel.hpp:249
void assembleTracerEquationWell(TrRe &tr, const Well &well)
Definition: TracerModel.hpp:342
void assembleTracerEquationVolume(TrRe &tr, const ElementContext &elemCtx, const Scalar scvVolume, const Scalar dt, unsigned I, unsigned I1)
Definition: TracerModel.hpp:283
void serialize(Restarter &)
This method writes the complete state of all tracer to the hard disk.
Definition: TracerModel.hpp:201
TracerBatch< TracerVector > & gas_
Definition: TracerModel.hpp:615
void assembleTracerEquations_()
Definition: TracerModel.hpp:377
void computeVolume_(LhsEval &freeVolume, const int tracerPhaseIdx, const ElementContext &elemCtx, unsigned scvIdx, unsigned timeIdx)
Definition: TracerModel.hpp:225
void beginTimeStep()
Definition: TracerModel.hpp:177
void serializeOp(Serializer &serializer)
Definition: TracerModel.hpp:215
void prepareTracerBatches()
Definition: TracerModel.hpp:133
void assembleTracerEquationFlux(TrRe &tr, const ElementContext &elemCtx, unsigned scvfIdx, unsigned I, unsigned J)
Definition: TracerModel.hpp:319
TracerBatch< TracerVector > & wat_
Definition: TracerModel.hpp:613
void endTimeStep()
Informs the tracer model that a time step has just been finished.
Definition: TracerModel.hpp:188
Simulator & simulator_
Definition: TracerModel.hpp:552
void deserialize(Restarter &)
This method restores the complete state of the tracer from disk.
Definition: TracerModel.hpp:211
A data handle sending multiple data store in vectors attached to cells.
Definition: VectorVectorDataHandle.hpp:49
Definition: AluGridVanguard.hpp:57
Definition: BlackoilPhases.hpp:27
Definition: TracerModel.hpp:48
UndefinedProperty type
Definition: TracerModel.hpp:49
Definition: TracerModel.hpp:563
const int phaseIdx_
Definition: TracerModel.hpp:565
void addTracer(const int idx, const TV &concentration)
Definition: TracerModel.hpp:601
bool operator==(const TracerBatch &rhs) const
Definition: TracerModel.hpp:572
static TracerBatch serializationTestObject()
Definition: TracerModel.hpp:578
std::vector< TV > concentrationInitial_
Definition: TracerModel.hpp:566
void serializeOp(Serializer &serializer)
Definition: TracerModel.hpp:591
TracerBatch(int phaseIdx=0)
Definition: TracerModel.hpp:597
int numTracer() const
Definition: TracerModel.hpp:599
std::vector< int > idx_
Definition: TracerModel.hpp:564
std::unique_ptr< TracerMatrix > mat
Definition: TracerModel.hpp:570
std::vector< TV > storageOfTimeIndex1_
Definition: TracerModel.hpp:568
std::vector< TV > residual_
Definition: TracerModel.hpp:569
std::vector< TV > concentration_
Definition: TracerModel.hpp:567