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::ranges::transform(coords, 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 this->simulator().vanguard().schedule()
297 .exitStatus().has_value());
298 }
299 }
300
306 template <class Context>
307 void boundary(BoundaryRateVector& values,
308 const Context& context,
309 unsigned spaceIdx,
310 unsigned /* timeIdx */) const
311 {
312 OPM_TIMEBLOCK_LOCAL(eclProblemBoundary, Subsystem::Assembly);
313 if (!context.intersection(spaceIdx).boundary())
314 return;
315
316 values.setNoFlow();
317
318 if (this->nonTrivialBoundaryConditions()) {
319 throw std::logic_error("boundary condition is not supported by compostional modeling yet");
320 }
321 }
322
329 template <class Context>
330 void initial(PrimaryVariables& values, const Context& context, unsigned spaceIdx, unsigned timeIdx) const
331 {
332 const unsigned globalDofIdx = context.globalSpaceIndex(spaceIdx, timeIdx);
333 const auto& initial_fs = initialFluidStates_[globalDofIdx];
334 Opm::CompositionalFluidState<Scalar, FluidSystem> fs;
335 for (unsigned p = 0; p < numPhases; ++p) { // TODO: assuming the phaseidx continuous
336 // pressure
337 fs.setPressure(p, initial_fs.pressure(p));
338
339 // saturation
340 fs.setSaturation(p, initial_fs.saturation(p));
341
342 // temperature
343 fs.setTemperature(initial_fs.temperature(p));
344 }
345
346
347 if (!zmf_initialization_) {
348 for (unsigned p = 0; p < numPhases; ++p) {
349 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
350 fs.setMoleFraction(p, compIdx, initial_fs.moleFraction(p, compIdx));
351 }
352 }
353
354 {
355 const auto& eos_type = getEosType();
356 typename FluidSystem::template ParameterCache<Scalar> paramCache(eos_type);
357 paramCache.updatePhase(fs, FluidSystem::oilPhaseIdx);
358 paramCache.updatePhase(fs, FluidSystem::gasPhaseIdx);
359 fs.setDensity(FluidSystem::oilPhaseIdx, FluidSystem::density(fs, paramCache, FluidSystem::oilPhaseIdx));
360 fs.setDensity(FluidSystem::gasPhaseIdx, FluidSystem::density(fs, paramCache, FluidSystem::gasPhaseIdx));
361 }
362 // determine the component fractions
363 Dune::FieldVector<Scalar, numComponents> z(0.0);
364 Scalar sumMoles = 0.0;
365 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
366 if (Indices::waterEnabled && phaseIdx == static_cast<unsigned int>(waterPhaseIdx)){
367 continue;
368 }
369 const auto saturation = fs.saturation(phaseIdx);
370 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
371 Scalar tmp = fs.molarity(phaseIdx, compIdx) * saturation;
372 tmp = max(tmp, 1e-8);
373 z[compIdx] += tmp;
374 sumMoles += tmp;
375 }
376 }
377 z /= sumMoles;
378 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
379 fs.setMoleFraction(compIdx, z[compIdx]);
380 }
381 } else {
382 // TODO: should we normalize the input?
383 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
384 fs.setMoleFraction(compIdx, initial_fs.moleFraction(compIdx));
385 }
386 }
387
388 // Set initial K and L
389 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
390 const auto& Ktmp = fs.wilsonK_(compIdx);
391 fs.setKvalue(compIdx, Ktmp);
392 }
393
394 const Scalar& Ltmp = -1.0;
395 fs.setLvalue(Ltmp);
396
397 values.assignNaive(fs);
398 }
399
400 void addToSourceDense(RateVector&, unsigned, unsigned) const override
401 {
402 // we do nothing for now
403 }
404
405 const InitialFluidState& initialFluidState(unsigned globalDofIdx) const
406 { return initialFluidStates_[globalDofIdx]; }
407
408 std::vector<InitialFluidState>& initialFluidStates()
409 { return initialFluidStates_; }
410
411 const std::vector<InitialFluidState>& initialFluidStates() const
412 { return initialFluidStates_; }
413
414 const FlowThresholdPressure<TypeTag>& thresholdPressure() const
415 {
416 assert( !thresholdPressures_.enableThresholdPressure() &&
417 " Threshold Pressures are not supported by compostional simulation ");
418 return thresholdPressures_;
419 }
420
421 const EclWriterType& eclWriter() const
422 { return *eclWriter_; }
423
424 EclWriterType& eclWriter()
425 { return *eclWriter_; }
426
427 // TODO: do we need this one?
428 template<class Serializer>
429 void serializeOp(Serializer& serializer)
430 {
431 serializer(static_cast<FlowProblemType&>(*this));
432 serializer(*eclWriter_);
433 }
434protected:
435
436 void updateExplicitQuantities_(int /* episodeIdx*/, int /* timeStepSize */, bool /* first_step_after_restart */) override
437 {
438 // we do nothing here for now
439 }
440
441 void readEquilInitialCondition_() override
442 {
443 throw std::logic_error("Equilibration is not supported by compositional modeling yet");
444 }
445
446 void readEclRestartSolution_()
447 {
448 throw std::logic_error("Restarting is not supported by compositional modeling yet");
449 }
450
451 void readExplicitInitialCondition_() override
452 {
453 readExplicitInitialConditionCompositional_();
454 }
455
456 void readExplicitInitialConditionCompositional_()
457 {
458 const auto& simulator = this->simulator();
459 const auto& vanguard = simulator.vanguard();
460 const auto& eclState = vanguard.eclState();
461 const auto& fp = eclState.fieldProps();
462 const bool has_pressure = fp.has_double("PRESSURE");
463 if (!has_pressure)
464 throw std::runtime_error("The ECL input file requires the presence of the PRESSURE "
465 "keyword if the model is initialized explicitly");
466
467 const bool has_xmf = fp.has_double("XMF");
468 const bool has_ymf = fp.has_double("YMF");
469 const bool has_zmf = fp.has_double("ZMF");
470 if ( !has_zmf && !(has_xmf && has_ymf) ) {
471 throw std::runtime_error("The ECL input file requires the presence of ZMF or XMF and YMF "
472 "keyword if the model is initialized explicitly");
473 }
474
475 if (has_zmf && (has_xmf || has_ymf)) {
476 throw std::runtime_error("The ECL input file can not handle explicit initialization "
477 "with both ZMF and XMF or YMF");
478 }
479
480 if (has_xmf != has_ymf) {
481 throw std::runtime_error("The ECL input file needs XMF and YMF combined to do the explicit "
482 "initializtion when using XMF or YMF");
483 }
484
485 const bool has_temp = fp.has_double("TEMPI");
486
487 // const bool has_gas = fp.has_double("SGAS");
488 assert(fp.has_double("SGAS"));
489
490 std::size_t numDof = this->model().numGridDof();
491
492 initialFluidStates_.resize(numDof);
493
494 std::vector<double> waterSaturationData;
495 std::vector<double> gasSaturationData;
496 std::vector<double> soilData;
497 std::vector<double> pressureData;
498 std::vector<double> tempiData;
499
500 const bool water_active = FluidSystem::phaseIsActive(waterPhaseIdx);
501 const bool gas_active = FluidSystem::phaseIsActive(gasPhaseIdx);
502 const bool oil_active = FluidSystem::phaseIsActive(oilPhaseIdx);
503
504 if (water_active && Indices::numPhases > 2)
505 waterSaturationData = fp.get_double("SWAT");
506 else
507 waterSaturationData.resize(numDof);
508
509 pressureData = fp.get_double("PRESSURE");
510
511 if (has_temp) {
512 tempiData = fp.get_double("TEMPI");
513 } else {
514 ; // TODO: throw?
515 }
516
517 if (gas_active) // && FluidSystem::phaseIsActive(oilPhaseIdx))
518 gasSaturationData = fp.get_double("SGAS");
519 else
520 gasSaturationData.resize(numDof);
521
522 for (std::size_t dofIdx = 0; dofIdx < numDof; ++dofIdx) {
523 auto& dofFluidState = initialFluidStates_[dofIdx];
524 // dofFluidState.setPvtRegionIndex(pvtRegionIndex(dofIdx));
525
526 Scalar temperatureLoc = tempiData[dofIdx];
527 assert(std::isfinite(temperatureLoc) && temperatureLoc > 0);
528 dofFluidState.setTemperature(temperatureLoc);
529
530 if (gas_active) {
531 dofFluidState.setSaturation(FluidSystem::gasPhaseIdx,
532 gasSaturationData[dofIdx]);
533 }
534 if (oil_active) {
535 dofFluidState.setSaturation(FluidSystem::oilPhaseIdx,
536 1.0
537 - waterSaturationData[dofIdx]
538 - gasSaturationData[dofIdx]);
539 }
540 if (water_active) {
541 dofFluidState.setSaturation(FluidSystem::waterPhaseIdx,
542 waterSaturationData[dofIdx]);
543 }
544
546 // set phase pressures
548 const Scalar pressure = pressureData[dofIdx]; // oil pressure (or gas pressure for water-gas system or water pressure for single phase)
549
550 // TODO: zero capillary pressure for now
551 const std::array<Scalar, numPhases> pc = {0};
552 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
553 if (!FluidSystem::phaseIsActive(phaseIdx))
554 continue;
555
556 if (Indices::oilEnabled)
557 dofFluidState.setPressure(phaseIdx, pressure + (pc[phaseIdx] - pc[oilPhaseIdx]));
558 else if (Indices::gasEnabled)
559 dofFluidState.setPressure(phaseIdx, pressure + (pc[phaseIdx] - pc[gasPhaseIdx]));
560 else if (Indices::waterEnabled)
561 // single (water) phase
562 dofFluidState.setPressure(phaseIdx, pressure);
563 }
564
565 if (has_xmf && has_ymf) {
566 const auto& xmfData = fp.get_double("XMF");
567 const auto& ymfData = fp.get_double("YMF");
568 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
569 const std::size_t data_idx = compIdx * numDof + dofIdx;
570 const Scalar xmf = xmfData[data_idx];
571 const Scalar ymf = ymfData[data_idx];
572
573 dofFluidState.setMoleFraction(FluidSystem::oilPhaseIdx, compIdx, xmf);
574 dofFluidState.setMoleFraction(FluidSystem::gasPhaseIdx, compIdx, ymf);
575 }
576 }
577
578 if (has_zmf) {
579 zmf_initialization_ = true;
580 const auto& zmfData = fp.get_double("ZMF");
581 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
582 const std::size_t data_idx = compIdx * numDof + dofIdx;
583 const Scalar zmf = zmfData[data_idx];
584 dofFluidState.setMoleFraction(compIdx, zmf);
585
586 if (gas_active) {
587 const auto ymf = (dofFluidState.saturation(FluidSystem::gasPhaseIdx) > 0.) ? zmf : Scalar{0};
588 dofFluidState.setMoleFraction(FluidSystem::gasPhaseIdx, compIdx, ymf);
589 }
590 if (oil_active) {
591 const auto xmf = (dofFluidState.saturation(FluidSystem::oilPhaseIdx) > 0.) ? zmf : Scalar{0};
592 dofFluidState.setMoleFraction(FluidSystem::oilPhaseIdx, compIdx, xmf);
593 }
594 }
595 }
596 }
597 }
598
599private:
600
601 void handleSolventBC(const BCProp::BCFace& /* bc */, RateVector& /* rate */) const override
602 {
603 throw std::logic_error("solvent is disabled for compositional modeling and you're trying to add solvent to BC");
604 }
605
606 void handlePolymerBC(const BCProp::BCFace& /* bc */, RateVector& /* rate */) const override
607 {
608 throw std::logic_error("polymer is disabled for compositional modeling and you're trying to add polymer to BC");
609 }
610
611 void handleMicrBC(const BCProp::BCFace& /* bc */, RateVector& /* rate */) const override
612 {
613 throw std::logic_error("MICP is disabled for compositional modeling and you're trying to add microbes to BC");
614 }
615
616 void handleOxygBC(const BCProp::BCFace& /* bc */, RateVector& /* rate */) const override
617 {
618 throw std::logic_error("MICP is disabled for compositional modeling and you're trying to add oxygen to BC");
619 }
620
621 void handleUreaBC(const BCProp::BCFace& /* bc */, RateVector& /* rate */) const override
622 {
623 throw std::logic_error("MICP is disabled for compositional modeling and you're trying to add urea to BC");
624 }
625
626 FlowThresholdPressure<TypeTag> thresholdPressures_;
627
628 std::vector<InitialFluidState> initialFluidStates_;
629
630 bool zmf_initialization_ {false};
631
632 bool enableEclOutput_{false};
633 std::unique_ptr<EclWriterType> eclWriter_;
634};
635
636} // namespace Opm
637
638#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:115
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:307
void initial(PrimaryVariables &values, const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Evaluate the initial value for a control volume.
Definition FlowProblemComp.hpp:330
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:1863
virtual void writeOutput(bool verbose)
Write the requested quantities of the current solution into the output files.
Definition FlowProblem.hpp:498
unsigned pvtRegionIndex(const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Returns the index of the relevant region for thermodynmic properties.
Definition FlowProblem.hpp:878
Scalar porosity(const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Definition FlowProblem.hpp:681
FlowProblem(Simulator &simulator)
Definition FlowProblem.hpp:220
static void registerParameters()
Registers all available parameters for the problem and the model.
Definition FlowProblem.hpp:185
virtual void endTimeStep()
Called by the simulator after each time integration.
Definition FlowProblem.hpp:418
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: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
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