opm-simulators
Loading...
Searching...
No Matches
FlowProblemComp.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 Copyright 2024 SINTEF Digital
5
6 This file is part of the Open Porous Media project (OPM).
7
8 OPM is free software: you can redistribute it and/or modify
9 it under the terms of the GNU General Public License as published by
10 the Free Software Foundation, either version 2 of the License, or
11 (at your option) any later version.
12
13 OPM is distributed in the hope that it will be useful,
14 but WITHOUT ANY WARRANTY; without even the implied warranty of
15 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
16 GNU General Public License for more details.
17
18 You should have received a copy of the GNU General Public License
19 along with OPM. If not, see <http://www.gnu.org/licenses/>.
20
21 Consult the COPYING file in the top-level source directory of this
22 module for the precise wording of the license and the list of
23 copyright holders.
24*/
30#ifndef OPM_FLOW_PROBLEM_COMP_HPP
31#define OPM_FLOW_PROBLEM_COMP_HPP
32
33
37
38#include <opm/material/fluidstates/CompositionalFluidState.hpp>
39
40#include <opm/material/thermal/EclThermalLawManager.hpp>
41
42#include <opm/input/eclipse/EclipseState/Compositional/CompositionalConfig.hpp>
43
44#include <algorithm>
45#include <functional>
46#include <set>
47#include <string>
48#include <vector>
49
50namespace Opm {
51
58template <class TypeTag>
59class FlowProblemComp : public FlowProblem<TypeTag>
60{
61 // TODO: the naming of the Types will be adjusted
62 using FlowProblemType = FlowProblem<TypeTag>;
63
64 using typename FlowProblemType::Scalar;
65 using typename FlowProblemType::Simulator;
66 using typename FlowProblemType::GridView;
67 using typename FlowProblemType::FluidSystem;
68 using typename FlowProblemType::Vanguard;
69
70 // might not be needed
71 using FlowProblemType::dim;
72 using FlowProblemType::dimWorld;
73
74 using FlowProblemType::numPhases;
75 using FlowProblemType::numComponents;
76
77 using FlowProblemType::gasPhaseIdx;
78 using FlowProblemType::oilPhaseIdx;
79 using FlowProblemType::waterPhaseIdx;
80
81 using typename FlowProblemType::Indices;
82 using typename FlowProblemType::PrimaryVariables;
84 using typename FlowProblemType::Evaluation;
85 using typename FlowProblemType::MaterialLaw;
86 using typename FlowProblemType::RateVector;
87
88 using InitialFluidState = CompositionalFluidState<Scalar, FluidSystem>;
90
91public:
94
98 static void registerParameters()
99 {
101
102 EclWriterType::registerParameters();
103
104 // tighter tolerance is needed for compositional modeling here
106 }
107
108 Opm::CompositionalConfig::EOSType getEosType() const
109 {
110 auto& simulator = this->simulator();
111 const auto& eclState = simulator.vanguard().eclState();
112 return eclState.compositionalConfig().eosType(0);
113 }
114
118 explicit FlowProblemComp(Simulator& simulator)
119 : FlowProblemType(simulator)
120 , thresholdPressures_(simulator)
121 {
122 eclWriter_ = std::make_unique<EclWriterType>(simulator);
124 }
125
130 {
131 // TODO: there should be room to remove duplication for this function,
132 // but there is relatively complicated logic in the function calls in this function
133 // some refactoring is needed for this function
134 FlowProblemType::finishInit();
135
136 auto& simulator = this->simulator();
137
138 auto finishTransmissibilities = [updated = false, this]() mutable {
139 if (updated) {
140 return;
141 }
142 this->transmissibilities_.finishInit(
143 [&vg = this->simulator().vanguard()](const unsigned int it) { return vg.gridIdxToEquilGridIdx(it); });
144 updated = true;
145 };
146 // TODO: we might need to do the same with FlowProblemBlackoil for parallel
147
148 finishTransmissibilities();
149
150 if (enableEclOutput_) {
151 eclWriter_->setTransmissibilities(&simulator.problem().eclTransmissibilities());
152 std::function<unsigned int(unsigned int)> equilGridToGrid = [&simulator](unsigned int i) {
153 return simulator.vanguard().gridEquilIdxToGridIdx(i);
154 };
155 eclWriter_->extractOutputTransAndNNC(equilGridToGrid);
156 }
157
158 const auto& eclState = simulator.vanguard().eclState();
159 const auto& schedule = simulator.vanguard().schedule();
160
161 // Set the start time of the simulation
162 simulator.setStartTime(schedule.getStartTime());
163 simulator.setEndTime(schedule.simTime(schedule.size() - 1));
164
165 // We want the episode index to be the same as the report step index to make
166 // things simpler, so we have to set the episode index to -1 because it is
167 // incremented by endEpisode(). The size of the initial time step and
168 // length of the initial episode is set to zero for the same reason.
169 simulator.setEpisodeIndex(-1);
170 simulator.setEpisodeLength(0.0);
171
172 // the "NOGRAV" keyword from Frontsim or setting the EnableGravity to false
173 // disables gravity, else the standard value of the gravity constant at sea level
174 // on earth is used
175 this->gravity_ = 0.0;
177 this->gravity_[dim - 1] = 9.80665;
178 if (!eclState.getInitConfig().hasGravity())
179 this->gravity_[dim - 1] = 0.0;
180
181 if (this->enableTuning_) {
182 // if support for the TUNING keyword is enabled, we get the initial time
183 // steping parameters from it instead of from command line parameters
184 const auto& tuning = schedule[0].tuning();
185 this->initialTimeStepSize_ = tuning.TSINIT.has_value() ? tuning.TSINIT.value() : -1.0;
186 this->maxTimeStepAfterWellEvent_ = tuning.TMAXWC;
187 }
188
189 this->initFluidSystem_();
190
191 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx)
192 && FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx)) {
193 this->maxOilSaturation_.resize(this->model().numGridDof(), 0.0);
194 }
195
196 this->readRockParameters_(simulator.vanguard().cellCenterDepths(), [&simulator](const unsigned idx) {
197 std::array<int, dim> coords;
198 simulator.vanguard().cartesianCoordinate(idx, coords);
199 std::transform(coords.begin(), coords.end(), coords.begin(),
200 [](const auto c) { return c + 1; });
201 return coords;
202 });
203 FlowProblemType::readMaterialParameters_();
204 FlowProblemType::readThermalParameters_();
205
206 // write the static output files (EGRID, INIT)
207 if (enableEclOutput_) {
208 eclWriter_->writeInit();
209 }
210
211 const auto& initconfig = eclState.getInitConfig();
212 if (initconfig.restartRequested())
213 readEclRestartSolution_();
214 else
215 this->readInitialCondition_();
216
217 FlowProblemType::updatePffDofData_();
218
220 const auto& vanguard = this->simulator().vanguard();
221 const auto& gridView = vanguard.gridView();
222 int numElements = gridView.size(/*codim=*/0);
223 this->polymer_.maxAdsorption.resize(numElements, 0.0);
224 }
225
226 /* readBoundaryConditions_();
227
228 // compute and set eq weights based on initial b values
229 computeAndSetEqWeights_();
230
231 if (enableDriftCompensation_) {
232 drift_.resize(this->model().numGridDof());
233 drift_ = 0.0;
234 } */
235
236 // TODO: check wether the following can work with compostional
237 if (this->enableVtkOutput_() && eclState.getIOConfig().initOnly()) {
238 simulator.setTimeStepSize(0.0);
240 }
241
242 // after finishing the initialization and writing the initial solution, we move
243 // to the first "real" episode/report step
244 // for restart the episode index and start is already set
245 if (!initconfig.restartRequested()) {
246 simulator.startNextEpisode(schedule.seconds(1));
247 simulator.setEpisodeIndex(0);
248 simulator.setTimeStepIndex(0);
249 }
250 }
251
255 void endTimeStep() override
256 {
258
259 // after the solution is updated, the values in output module also needs to be updated
260 this->eclWriter_->mutableOutputModule().invalidateLocalData();
261
262 // For CpGrid with LGRs, ecl/vtk output is not supported yet.
263 const auto& grid = this->simulator().vanguard().gridView().grid();
264
265 using GridType = std::remove_cv_t<std::remove_reference_t<decltype(grid)>>;
266 constexpr bool isCpGrid = std::is_same_v<GridType, Dune::CpGrid>;
267 if (!isCpGrid || (grid.maxLevel() == 0)) {
268 this->eclWriter_->evalSummaryState(! this->episodeWillBeOver());
269 }
270 }
271
272 void writeReports(const SimulatorTimer& timer) {
273 if (enableEclOutput_){
274 eclWriter_->writeReports(timer);
275 }
276 }
277
282 void writeOutput(bool verbose) override
283 {
285
286 if (! this->enableEclOutput_) {
287 return;
288 }
289
290 const auto isSubStep = !this->episodeWillBeOver();
291
293 auto localCellData = data::Solution {};
294
295 this->eclWriter_->writeOutput(std::move(localCellData), isSubStep);
296 }
297 }
298
304 template <class Context>
305 void boundary(BoundaryRateVector& values,
306 const Context& context,
307 unsigned spaceIdx,
308 unsigned /* timeIdx */) const
309 {
310 OPM_TIMEBLOCK_LOCAL(eclProblemBoundary, Subsystem::Assembly);
311 if (!context.intersection(spaceIdx).boundary())
312 return;
313
314 values.setNoFlow();
315
316 if (this->nonTrivialBoundaryConditions()) {
317 throw std::logic_error("boundary condition is not supported by compostional modeling yet");
318 }
319 }
320
327 template <class Context>
328 void initial(PrimaryVariables& values, const Context& context, unsigned spaceIdx, unsigned timeIdx) const
329 {
330 const unsigned globalDofIdx = context.globalSpaceIndex(spaceIdx, timeIdx);
331 const auto& initial_fs = initialFluidStates_[globalDofIdx];
332 Opm::CompositionalFluidState<Scalar, FluidSystem> fs;
333 for (unsigned p = 0; p < numPhases; ++p) { // TODO: assuming the phaseidx continuous
334 // pressure
335 fs.setPressure(p, initial_fs.pressure(p));
336
337 // saturation
338 fs.setSaturation(p, initial_fs.saturation(p));
339
340 // temperature
341 fs.setTemperature(initial_fs.temperature(p));
342 }
343
344
345 if (!zmf_initialization_) {
346 for (unsigned p = 0; p < numPhases; ++p) {
347 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
348 fs.setMoleFraction(p, compIdx, initial_fs.moleFraction(p, compIdx));
349 }
350 }
351
352 {
353 const auto& eos_type = getEosType();
354 typename FluidSystem::template ParameterCache<Scalar> paramCache(eos_type);
355 paramCache.updatePhase(fs, FluidSystem::oilPhaseIdx);
356 paramCache.updatePhase(fs, FluidSystem::gasPhaseIdx);
357 fs.setDensity(FluidSystem::oilPhaseIdx, FluidSystem::density(fs, paramCache, FluidSystem::oilPhaseIdx));
358 fs.setDensity(FluidSystem::gasPhaseIdx, FluidSystem::density(fs, paramCache, FluidSystem::gasPhaseIdx));
359 }
360 // determine the component fractions
361 Dune::FieldVector<Scalar, numComponents> z(0.0);
362 Scalar sumMoles = 0.0;
363 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
364 if (Indices::waterEnabled && phaseIdx == static_cast<unsigned int>(waterPhaseIdx)){
365 continue;
366 }
367 const auto saturation = fs.saturation(phaseIdx);
368 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
369 Scalar tmp = fs.molarity(phaseIdx, compIdx) * saturation;
370 tmp = max(tmp, 1e-8);
371 z[compIdx] += tmp;
372 sumMoles += tmp;
373 }
374 }
375 z /= sumMoles;
376 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
377 fs.setMoleFraction(compIdx, z[compIdx]);
378 }
379 } else {
380 // TODO: should we normalize the input?
381 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
382 fs.setMoleFraction(compIdx, initial_fs.moleFraction(compIdx));
383 }
384 }
385
386 // Set initial K and L
387 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
388 const auto& Ktmp = fs.wilsonK_(compIdx);
389 fs.setKvalue(compIdx, Ktmp);
390 }
391
392 const Scalar& Ltmp = -1.0;
393 fs.setLvalue(Ltmp);
394
395 values.assignNaive(fs);
396 }
397
398 void addToSourceDense(RateVector&, unsigned, unsigned) const override
399 {
400 // we do nothing for now
401 }
402
403 const InitialFluidState& initialFluidState(unsigned globalDofIdx) const
404 { return initialFluidStates_[globalDofIdx]; }
405
406 std::vector<InitialFluidState>& initialFluidStates()
407 { return initialFluidStates_; }
408
409 const std::vector<InitialFluidState>& initialFluidStates() const
410 { return initialFluidStates_; }
411
412 const FlowThresholdPressure<TypeTag>& thresholdPressure() const
413 {
414 assert( !thresholdPressures_.enableThresholdPressure() &&
415 " Threshold Pressures are not supported by compostional simulation ");
416 return thresholdPressures_;
417 }
418
419 const EclWriterType& eclWriter() const
420 { return *eclWriter_; }
421
422 EclWriterType& eclWriter()
423 { return *eclWriter_; }
424
425 // TODO: do we need this one?
426 template<class Serializer>
427 void serializeOp(Serializer& serializer)
428 {
429 serializer(static_cast<FlowProblemType&>(*this));
430 serializer(*eclWriter_);
431 }
432protected:
433
434 void updateExplicitQuantities_(int /* episodeIdx*/, int /* timeStepSize */, bool /* first_step_after_restart */) override
435 {
436 // we do nothing here for now
437 }
438
439 void readEquilInitialCondition_() override
440 {
441 throw std::logic_error("Equilibration is not supported by compositional modeling yet");
442 }
443
444 void readEclRestartSolution_()
445 {
446 throw std::logic_error("Restarting is not supported by compositional modeling yet");
447 }
448
449 void readExplicitInitialCondition_() override
450 {
451 readExplicitInitialConditionCompositional_();
452 }
453
454 void readExplicitInitialConditionCompositional_()
455 {
456 const auto& simulator = this->simulator();
457 const auto& vanguard = simulator.vanguard();
458 const auto& eclState = vanguard.eclState();
459 const auto& fp = eclState.fieldProps();
460 const bool has_pressure = fp.has_double("PRESSURE");
461 if (!has_pressure)
462 throw std::runtime_error("The ECL input file requires the presence of the PRESSURE "
463 "keyword if the model is initialized explicitly");
464
465 const bool has_xmf = fp.has_double("XMF");
466 const bool has_ymf = fp.has_double("YMF");
467 const bool has_zmf = fp.has_double("ZMF");
468 if ( !has_zmf && !(has_xmf && has_ymf) ) {
469 throw std::runtime_error("The ECL input file requires the presence of ZMF or XMF and YMF "
470 "keyword if the model is initialized explicitly");
471 }
472
473 if (has_zmf && (has_xmf || has_ymf)) {
474 throw std::runtime_error("The ECL input file can not handle explicit initialization "
475 "with both ZMF and XMF or YMF");
476 }
477
478 if (has_xmf != has_ymf) {
479 throw std::runtime_error("The ECL input file needs XMF and YMF combined to do the explicit "
480 "initializtion when using XMF or YMF");
481 }
482
483 const bool has_temp = fp.has_double("TEMPI");
484
485 // const bool has_gas = fp.has_double("SGAS");
486 assert(fp.has_double("SGAS"));
487
488 std::size_t numDof = this->model().numGridDof();
489
490 initialFluidStates_.resize(numDof);
491
492 std::vector<double> waterSaturationData;
493 std::vector<double> gasSaturationData;
494 std::vector<double> soilData;
495 std::vector<double> pressureData;
496 std::vector<double> tempiData;
497
498 const bool water_active = FluidSystem::phaseIsActive(waterPhaseIdx);
499 const bool gas_active = FluidSystem::phaseIsActive(gasPhaseIdx);
500 const bool oil_active = FluidSystem::phaseIsActive(oilPhaseIdx);
501
502 if (water_active && Indices::numPhases > 2)
503 waterSaturationData = fp.get_double("SWAT");
504 else
505 waterSaturationData.resize(numDof);
506
507 pressureData = fp.get_double("PRESSURE");
508
509 if (has_temp) {
510 tempiData = fp.get_double("TEMPI");
511 } else {
512 ; // TODO: throw?
513 }
514
515 if (gas_active) // && FluidSystem::phaseIsActive(oilPhaseIdx))
516 gasSaturationData = fp.get_double("SGAS");
517 else
518 gasSaturationData.resize(numDof);
519
520 for (std::size_t dofIdx = 0; dofIdx < numDof; ++dofIdx) {
521 auto& dofFluidState = initialFluidStates_[dofIdx];
522 // dofFluidState.setPvtRegionIndex(pvtRegionIndex(dofIdx));
523
524 Scalar temperatureLoc = tempiData[dofIdx];
525 assert(std::isfinite(temperatureLoc) && temperatureLoc > 0);
526 dofFluidState.setTemperature(temperatureLoc);
527
528 if (gas_active) {
529 dofFluidState.setSaturation(FluidSystem::gasPhaseIdx,
530 gasSaturationData[dofIdx]);
531 }
532 if (oil_active) {
533 dofFluidState.setSaturation(FluidSystem::oilPhaseIdx,
534 1.0
535 - waterSaturationData[dofIdx]
536 - gasSaturationData[dofIdx]);
537 }
538 if (water_active) {
539 dofFluidState.setSaturation(FluidSystem::waterPhaseIdx,
540 waterSaturationData[dofIdx]);
541 }
542
544 // set phase pressures
546 const Scalar pressure = pressureData[dofIdx]; // oil pressure (or gas pressure for water-gas system or water pressure for single phase)
547
548 // TODO: zero capillary pressure for now
549 const std::array<Scalar, numPhases> pc = {0};
550 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
551 if (!FluidSystem::phaseIsActive(phaseIdx))
552 continue;
553
554 if (Indices::oilEnabled)
555 dofFluidState.setPressure(phaseIdx, pressure + (pc[phaseIdx] - pc[oilPhaseIdx]));
556 else if (Indices::gasEnabled)
557 dofFluidState.setPressure(phaseIdx, pressure + (pc[phaseIdx] - pc[gasPhaseIdx]));
558 else if (Indices::waterEnabled)
559 // single (water) phase
560 dofFluidState.setPressure(phaseIdx, pressure);
561 }
562
563 if (has_xmf && has_ymf) {
564 const auto& xmfData = fp.get_double("XMF");
565 const auto& ymfData = fp.get_double("YMF");
566 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
567 const std::size_t data_idx = compIdx * numDof + dofIdx;
568 const Scalar xmf = xmfData[data_idx];
569 const Scalar ymf = ymfData[data_idx];
570
571 dofFluidState.setMoleFraction(FluidSystem::oilPhaseIdx, compIdx, xmf);
572 dofFluidState.setMoleFraction(FluidSystem::gasPhaseIdx, compIdx, ymf);
573 }
574 }
575
576 if (has_zmf) {
577 zmf_initialization_ = true;
578 const auto& zmfData = fp.get_double("ZMF");
579 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
580 const std::size_t data_idx = compIdx * numDof + dofIdx;
581 const Scalar zmf = zmfData[data_idx];
582 dofFluidState.setMoleFraction(compIdx, zmf);
583
584 if (gas_active) {
585 const auto ymf = (dofFluidState.saturation(FluidSystem::gasPhaseIdx) > 0.) ? zmf : Scalar{0};
586 dofFluidState.setMoleFraction(FluidSystem::gasPhaseIdx, compIdx, ymf);
587 }
588 if (oil_active) {
589 const auto xmf = (dofFluidState.saturation(FluidSystem::oilPhaseIdx) > 0.) ? zmf : Scalar{0};
590 dofFluidState.setMoleFraction(FluidSystem::oilPhaseIdx, compIdx, xmf);
591 }
592 }
593 }
594 }
595 }
596
597private:
598
599 void handleSolventBC(const BCProp::BCFace& /* bc */, RateVector& /* rate */) const override
600 {
601 throw std::logic_error("solvent is disabled for compositional modeling and you're trying to add solvent to BC");
602 }
603
604 void handlePolymerBC(const BCProp::BCFace& /* bc */, RateVector& /* rate */) const override
605 {
606 throw std::logic_error("polymer is disabled for compositional modeling and you're trying to add polymer to BC");
607 }
608
609 void handleMicrBC(const BCProp::BCFace& /* bc */, RateVector& /* rate */) const override
610 {
611 throw std::logic_error("MICP is disabled for compositional modeling and you're trying to add microbes to BC");
612 }
613
614 void handleOxygBC(const BCProp::BCFace& /* bc */, RateVector& /* rate */) const override
615 {
616 throw std::logic_error("MICP is disabled for compositional modeling and you're trying to add oxygen to BC");
617 }
618
619 void handleUreaBC(const BCProp::BCFace& /* bc */, RateVector& /* rate */) const override
620 {
621 throw std::logic_error("MICP is disabled for compositional modeling and you're trying to add urea to BC");
622 }
623
624 FlowThresholdPressure<TypeTag> thresholdPressures_;
625
626 std::vector<InitialFluidState> initialFluidStates_;
627
628 bool zmf_initialization_ {false};
629
630 bool enableEclOutput_{false};
631 std::unique_ptr<EclWriterType> eclWriter_;
632};
633
634} // namespace Opm
635
636#endif // OPM_FLOW_PROBLEM_COMP_HPP
This problem simulates an input file given in the data format used by the commercial ECLiPSE simulato...
This class calculates the threshold pressure for grid faces according to the Eclipse Reference Manual...
Output module for the results black oil model writing in ECL binary format.
Collects necessary output values and pass it to opm-common's ECL output.
Definition EclWriter.hpp:114
void writeOutput(bool verbose) override
Write the requested quantities of the current solution into the output files.
Definition FlowProblemComp.hpp:282
void finishInit()
Called by the Opm::Simulator in order to initialize the problem.
Definition FlowProblemComp.hpp:129
FlowProblemComp(Simulator &simulator)
Definition FlowProblemComp.hpp:118
void endTimeStep() override
Called by the simulator after each time integration.
Definition FlowProblemComp.hpp:255
void boundary(BoundaryRateVector &values, const Context &context, unsigned spaceIdx, unsigned) const
Evaluate the boundary conditions for a boundary segment.
Definition FlowProblemComp.hpp:305
void initial(PrimaryVariables &values, const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Evaluate the initial value for a control volume.
Definition FlowProblemComp.hpp:328
static void registerParameters()
Registers all available parameters for the problem and the model.
Definition FlowProblemComp.hpp:98
virtual bool episodeWillBeOver() const
Whether or not the current episode will end at the end of the current time step.
Definition FlowProblem.hpp:1733
virtual void writeOutput(bool verbose)
Write the requested quantities of the current solution into the output files.
Definition FlowProblem.hpp:479
unsigned pvtRegionIndex(const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Returns the index of the relevant region for thermodynmic properties.
Definition FlowProblem.hpp:826
Scalar porosity(const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Definition FlowProblem.hpp:659
FlowProblem(Simulator &simulator)
Definition FlowProblem.hpp:216
static void registerParameters()
Registers all available parameters for the problem and the model.
Definition FlowProblem.hpp:181
virtual void endTimeStep()
Called by the simulator after each time integration.
Definition FlowProblem.hpp:404
This class calculates the threshold pressure for grid faces according to the Eclipse Reference Manual...
Definition FlowThresholdPressure.hpp:59
Definition SimulatorTimer.hpp:39
This file contains a set of helper functions used by VFPProd / VFPInj.
Definition blackoilbioeffectsmodules.hh:43
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property (equivalent to old macro GET_PROP_TYPE(....
Definition propertysystem.hh:233
constexpr auto getPropValue()
get the value data member of a property
Definition propertysystem.hh:240
void SetDefault(decltype(Param::value) new_value)
Set a runtime parameter.
Definition parametersystem.hpp:212
auto Get(bool errorIfNotRegistered=true)
Retrieve a runtime parameter.
Definition parametersystem.hpp:187