BlackoilModel.hpp
Go to the documentation of this file.
1/*
2 Copyright 2013, 2015 SINTEF ICT, Applied Mathematics.
3 Copyright 2014, 2015 Dr. Blatt - HPC-Simulation-Software & Services
4 Copyright 2014, 2015 Statoil ASA.
5 Copyright 2015 NTNU
6 Copyright 2015, 2016, 2017 IRIS AS
7
8 This file is part of the Open Porous Media project (OPM).
9
10 OPM is free software: you can redistribute it and/or modify
11 it under the terms of the GNU General Public License as published by
12 the Free Software Foundation, either version 3 of the License, or
13 (at your option) any later version.
14
15 OPM is distributed in the hope that it will be useful,
16 but WITHOUT ANY WARRANTY; without even the implied warranty of
17 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
18 GNU General Public License for more details.
19
20 You should have received a copy of the GNU General Public License
21 along with OPM. If not, see <http://www.gnu.org/licenses/>.
22*/
23
24#ifndef OPM_BLACKOILMODEL_HEADER_INCLUDED
25#define OPM_BLACKOILMODEL_HEADER_INCLUDED
26
28
34
36
40
42
44
45#include <memory>
46#include <tuple>
47#include <vector>
48
49#include <fmt/format.h>
50
51namespace Opm {
52
59template <class TypeTag>
61{
62public:
63 // --------- Types and enums ---------
77
78 static constexpr bool enableSaltPrecipitation = getPropValue<TypeTag, Properties::EnableSaltPrecipitation>();
79
80 static constexpr int numEq = Indices::numEq;
81 static constexpr int contiSolventEqIdx = Indices::contiSolventEqIdx;
82 static constexpr int contiZfracEqIdx = Indices::contiZfracEqIdx;
83 static constexpr int contiPolymerEqIdx = Indices::contiPolymerEqIdx;
84 static constexpr int contiEnergyEqIdx = Indices::contiEnergyEqIdx;
85 static constexpr int contiPolymerMWEqIdx = Indices::contiPolymerMWEqIdx;
86 static constexpr int contiFoamEqIdx = Indices::contiFoamEqIdx;
87 static constexpr int contiBrineEqIdx = Indices::contiBrineEqIdx;
88 static constexpr int contiMicrobialEqIdx = Indices::contiMicrobialEqIdx;
89 static constexpr int contiOxygenEqIdx = Indices::contiOxygenEqIdx;
90 static constexpr int contiUreaEqIdx = Indices::contiUreaEqIdx;
91 static constexpr int contiBiofilmEqIdx = Indices::contiBiofilmEqIdx;
92 static constexpr int contiCalciteEqIdx = Indices::contiCalciteEqIdx;
93 static constexpr int solventSaturationIdx = Indices::solventSaturationIdx;
94 static constexpr int zFractionIdx = Indices::zFractionIdx;
95 static constexpr int polymerConcentrationIdx = Indices::polymerConcentrationIdx;
96 static constexpr int polymerMoleWeightIdx = Indices::polymerMoleWeightIdx;
97 static constexpr int temperatureIdx = Indices::temperatureIdx;
98 static constexpr int foamConcentrationIdx = Indices::foamConcentrationIdx;
99 static constexpr int saltConcentrationIdx = Indices::saltConcentrationIdx;
100 static constexpr int microbialConcentrationIdx = Indices::microbialConcentrationIdx;
101 static constexpr int oxygenConcentrationIdx = Indices::oxygenConcentrationIdx;
102 static constexpr int ureaConcentrationIdx = Indices::ureaConcentrationIdx;
103 static constexpr int biofilmVolumeFractionIdx = Indices::biofilmVolumeFractionIdx;
104 static constexpr int calciteVolumeFractionIdx = Indices::calciteVolumeFractionIdx;
105
106 using VectorBlockType = Dune::FieldVector<Scalar, numEq>;
107 using MatrixBlockType = typename SparseMatrixAdapter::MatrixBlock;
108 using Mat = typename SparseMatrixAdapter::IstlMatrix;
109 using BVector = Dune::BlockVector<VectorBlockType>;
110
112
113 // Helper structs
115 std::pair<std::vector<double>, std::vector<int>> cnvPvSplit;
116 std::vector<unsigned> ixCells;
117 };
118
124 };
125
126 // Output debug flags for which tolerances used
127 enum class DebugFlags {
128 STRICT = 0,
129 RELAXED = 1,
130 TUNINGDP = 2
131 };
132
133 // --------- Public methods ---------
134
143 const ModelParameters& param,
144 BlackoilWellModel<TypeTag>& well_model,
145 const bool terminal_output);
146
147 bool isParallel() const
148 { return grid_.comm().size() > 1; }
149
150 const EclipseState& eclState() const
151 { return simulator_.vanguard().eclState(); }
152
156
158 const int iteration,
159 const int minIter,
160 const int maxIter,
161 const SimulatorTimerInterface& timer);
162
170 template <class NonlinearSolverType>
171 SimulatorReportSingle nonlinearIteration(const int iteration,
172 const SimulatorTimerInterface& timer,
173 NonlinearSolverType& nonlinear_solver);
174
175 template <class NonlinearSolverType>
177 const SimulatorTimerInterface& timer,
178 NonlinearSolverType& nonlinear_solver);
179
182 const int iterationIdx);
183
184 // compute the "relative" change of the solution between time steps
185 Scalar relativeChange() const;
186
189 { return simulator_.model().newtonMethod().linearSolver().iterations (); }
190
191 // Obtain reference to linear solver setup time
193 { return linear_solve_setup_time_; }
194
198
200 void updateSolution(const BVector& dx);
201
204 void storeSolutionUpdate(const BVector& dx);
205 MaxSolutionUpdateData getMaxSolutionUpdate(const std::vector<unsigned>& ixCells);
206
209 { return terminal_output_; }
210
211 std::tuple<Scalar,Scalar>
213 const Scalar pvSumLocal,
214 const Scalar numAquiferPvSumLocal,
215 std::vector<Scalar>& R_sum,
216 std::vector<Scalar>& maxCoeff,
217 std::vector<Scalar>& B_avg);
218
222 std::pair<Scalar,Scalar>
223 localConvergenceData(std::vector<Scalar>& R_sum,
224 std::vector<Scalar>& maxCoeff,
225 std::vector<Scalar>& B_avg,
226 std::vector<int>& maxCoeffCell);
227
232 CnvPvSplitData characteriseCnvPvSplit(const std::vector<Scalar>& B_avg, const double dt);
233
236 void convergencePerCell(const std::vector<Scalar>& B_avg,
237 const double dt,
238 const double tol_cnv,
239 const double tol_cnv_energy,
240 const int iteration);
241
242 void updateTUNING(const Tuning& tuning);
243 void updateTUNINGDP(const TuningDp& tuning_dp);
244
246 getReservoirConvergence(const double reportTime,
247 const double dt,
248 const int iteration,
249 const int maxIter,
250 std::vector<Scalar>& B_avg,
251 std::vector<Scalar>& residual_norms);
252
261 const int iteration,
262 const int maxIter,
263 std::vector<Scalar>& residual_norms);
264
266 int numPhases() const
267 { return Indices::numPhases; }
268
270 template<class T>
271 std::vector<std::vector<Scalar> >
272 computeFluidInPlace(const T&, const std::vector<int>& fipnum) const
273 { return this->computeFluidInPlace(fipnum); }
274
276 std::vector<std::vector<Scalar> >
277 computeFluidInPlace(const std::vector<int>& /*fipnum*/) const;
278
279 const Simulator& simulator() const
280 { return simulator_; }
281
283 { return simulator_; }
284
287 { return failureReport_; }
288
291
293 const std::vector<SimulatorReport>& domainAccumulatedReports() const;
294
296 void writeNonlinearIterationsPerCell(const std::filesystem::path& odir) const;
297
298 const std::vector<StepReport>& stepReports() const
299 { return convergence_reports_; }
300
301 void writePartitions(const std::filesystem::path& odir) const;
302
306 { return well_model_; }
307
309 wellModel() const
310 { return well_model_; }
311
313 { simulator_.problem().beginEpisode(); }
314
316 { simulator_.problem().endEpisode(); }
317
318 template<class FluidState, class Residual>
319 void getMaxCoeff(const unsigned cell_idx,
320 const IntensiveQuantities& intQuants,
321 const FluidState& fs,
322 const Residual& modelResid,
323 const Scalar pvValue,
324 std::vector<Scalar>& B_avg,
325 std::vector<Scalar>& R_sum,
326 std::vector<Scalar>& maxCoeff,
327 std::vector<int>& maxCoeffCell);
328
330 const ModelParameters& param() const
331 { return param_; }
332
335 { return compNames_; }
336
338 bool hasNlddSolver() const
339 { return nlddSolver_ != nullptr; }
340
341protected:
342 // --------- Data members ---------
344 const Grid& grid_;
345 static constexpr bool has_solvent_ = getPropValue<TypeTag, Properties::EnableSolvent>();
346 static constexpr bool has_extbo_ = getPropValue<TypeTag, Properties::EnableExtbo>();
347 static constexpr bool has_polymer_ = getPropValue<TypeTag, Properties::EnablePolymer>();
348 static constexpr bool has_polymermw_ = getPropValue<TypeTag, Properties::EnablePolymerMW>();
349 static constexpr bool has_energy_ = getPropValue<TypeTag, Properties::EnergyModuleType>() == EnergyModules::FullyImplicitThermal;
350 static constexpr bool has_foam_ = getPropValue<TypeTag, Properties::EnableFoam>();
351 static constexpr bool has_brine_ = getPropValue<TypeTag, Properties::EnableBrine>();
352 static constexpr bool has_bioeffects_ = getPropValue<TypeTag, Properties::EnableBioeffects>();
353 static constexpr bool has_micp_ = Indices::enableMICP;
354
357
358 // Well Model
360
364 long int global_nc_;
365
366 std::vector<std::vector<Scalar>> residual_norms_history_;
369
371
372 std::vector<StepReport> convergence_reports_;
374
375 std::unique_ptr<BlackoilModelNldd<TypeTag>> nlddSolver_;
377
378private:
379 Scalar dpMaxRel() const { return param_.dp_max_rel_; }
380 Scalar dsMax() const { return param_.ds_max_; }
381 Scalar drMaxRel() const { return param_.dr_max_rel_; }
382 Scalar maxResidualAllowed() const { return param_.max_residual_allowed_; }
383 double linear_solve_setup_time_;
384 std::vector<bool> wasSwitched_;
385};
386
387} // namespace Opm
388
390
391#endif // OPM_BLACKOILMODEL_HEADER_INCLUDED
Implementation of penalty cards for three-phase black oil.
Definition: BlackoilModelConvergenceMonitor.hpp:35
Definition: BlackoilModel.hpp:61
static constexpr int contiEnergyEqIdx
Definition: BlackoilModel.hpp:84
SimulatorReportSingle failureReport_
Definition: BlackoilModel.hpp:356
BlackoilModel(Simulator &simulator, const ModelParameters &param, BlackoilWellModel< TypeTag > &well_model, const bool terminal_output)
Definition: BlackoilModel_impl.hpp:51
static constexpr int biofilmVolumeFractionIdx
Definition: BlackoilModel.hpp:103
std::vector< std::vector< Scalar > > computeFluidInPlace(const T &, const std::vector< int > &fipnum) const
Wrapper required due to not following generic API.
Definition: BlackoilModel.hpp:272
static constexpr bool has_energy_
Definition: BlackoilModel.hpp:349
BlackoilModelParameters< Scalar > ModelParameters
Definition: BlackoilModel.hpp:76
ModelParameters param_
Definition: BlackoilModel.hpp:355
int linearIterationsLastSolve() const
Number of linear iterations used in last call to solveJacobianSystem().
Definition: BlackoilModel.hpp:188
SimulatorReportSingle nonlinearIteration(const int iteration, const SimulatorTimerInterface &timer, NonlinearSolverType &nonlinear_solver)
Definition: BlackoilModel_impl.hpp:214
void writeNonlinearIterationsPerCell(const std::filesystem::path &odir) const
Write the number of nonlinear iterations per cell to a file in ResInsight compatible format.
Definition: BlackoilModel_impl.hpp:1204
static constexpr bool has_extbo_
Definition: BlackoilModel.hpp:346
SimulatorReportSingle assembleReservoir(const SimulatorTimerInterface &, const int iterationIdx)
Assemble the residual and Jacobian of the nonlinear system.
Definition: BlackoilModel_impl.hpp:344
GetPropType< TypeTag, Properties::SparseMatrixAdapter > SparseMatrixAdapter
Definition: BlackoilModel.hpp:68
ConvergenceReport getConvergence(const SimulatorTimerInterface &timer, const int iteration, const int maxIter, std::vector< Scalar > &residual_norms)
Definition: BlackoilModel_impl.hpp:1146
ConvergenceReport getReservoirConvergence(const double reportTime, const double dt, const int iteration, const int maxIter, std::vector< Scalar > &B_avg, std::vector< Scalar > &residual_norms)
Definition: BlackoilModel_impl.hpp:851
GetPropType< TypeTag, Properties::IntensiveQuantities > IntensiveQuantities
Definition: BlackoilModel.hpp:67
static constexpr bool has_foam_
Definition: BlackoilModel.hpp:350
static constexpr int foamConcentrationIdx
Definition: BlackoilModel.hpp:98
Scalar current_relaxation_
Definition: BlackoilModel.hpp:367
GetPropType< TypeTag, Properties::FluidSystem > FluidSystem
Definition: BlackoilModel.hpp:71
static constexpr int microbialConcentrationIdx
Definition: BlackoilModel.hpp:100
static constexpr int numEq
Definition: BlackoilModel.hpp:80
const std::vector< SimulatorReport > & domainAccumulatedReports() const
return the statistics of local solves accumulated for each domain on this rank
Definition: BlackoilModel_impl.hpp:1194
static constexpr bool has_polymermw_
Definition: BlackoilModel.hpp:348
static constexpr int contiFoamEqIdx
Definition: BlackoilModel.hpp:86
static constexpr int polymerMoleWeightIdx
Definition: BlackoilModel.hpp:96
void endReportStep()
Definition: BlackoilModel.hpp:315
BlackoilModelConvergenceMonitor< Scalar > conv_monitor_
Definition: BlackoilModel.hpp:376
static constexpr int contiCalciteEqIdx
Definition: BlackoilModel.hpp:92
std::pair< Scalar, Scalar > localConvergenceData(std::vector< Scalar > &R_sum, std::vector< Scalar > &maxCoeff, std::vector< Scalar > &B_avg, std::vector< int > &maxCoeffCell)
Get reservoir quantities on this process needed for convergence calculations.
Definition: BlackoilModel_impl.hpp:694
GetPropType< TypeTag, Properties::ElementContext > ElementContext
Definition: BlackoilModel.hpp:66
bool isParallel() const
Definition: BlackoilModel.hpp:147
const std::vector< StepReport > & stepReports() const
Definition: BlackoilModel.hpp:298
std::vector< StepReport > convergence_reports_
Definition: BlackoilModel.hpp:372
static constexpr int contiMicrobialEqIdx
Definition: BlackoilModel.hpp:88
bool terminalOutputEnabled() const
Return true if output to cout is wanted.
Definition: BlackoilModel.hpp:208
double & linearSolveSetupTime()
Definition: BlackoilModel.hpp:192
int numPhases() const
The number of active fluid phases in the model.
Definition: BlackoilModel.hpp:266
const SimulatorReport & localAccumulatedReports() const
return the statistics of local solves accumulated for this rank
Definition: BlackoilModel_impl.hpp:1183
std::tuple< Scalar, Scalar > convergenceReduction(Parallel::Communication comm, const Scalar pvSumLocal, const Scalar numAquiferPvSumLocal, std::vector< Scalar > &R_sum, std::vector< Scalar > &maxCoeff, std::vector< Scalar > &B_avg)
Definition: BlackoilModel_impl.hpp:634
static constexpr int contiBiofilmEqIdx
Definition: BlackoilModel.hpp:91
const Grid & grid_
Definition: BlackoilModel.hpp:344
static constexpr int contiUreaEqIdx
Definition: BlackoilModel.hpp:90
static constexpr int temperatureIdx
Definition: BlackoilModel.hpp:97
GetPropType< TypeTag, Properties::SolutionVector > SolutionVector
Definition: BlackoilModel.hpp:69
const ComponentName & compNames() const
Returns const reference to component names.
Definition: BlackoilModel.hpp:334
GetPropType< TypeTag, Properties::PrimaryVariables > PrimaryVariables
Definition: BlackoilModel.hpp:70
void initialLinearization(SimulatorReportSingle &report, const int iteration, const int minIter, const int maxIter, const SimulatorTimerInterface &timer)
Definition: BlackoilModel_impl.hpp:155
void updateSolution(const BVector &dx)
Apply an update to the primary variables.
Definition: BlackoilModel_impl.hpp:507
long int global_nc_
The number of cells of the global grid.
Definition: BlackoilModel.hpp:364
DebugFlags
Definition: BlackoilModel.hpp:127
GetPropType< TypeTag, Properties::Indices > Indices
Definition: BlackoilModel.hpp:72
void updateTUNING(const Tuning &tuning)
Definition: BlackoilModel_impl.hpp:826
typename SparseMatrixAdapter::MatrixBlock MatrixBlockType
Definition: BlackoilModel.hpp:107
const ModelParameters & param() const
Returns const reference to model parameters.
Definition: BlackoilModel.hpp:330
static constexpr int contiOxygenEqIdx
Definition: BlackoilModel.hpp:89
void getMaxCoeff(const unsigned cell_idx, const IntensiveQuantities &intQuants, const FluidState &fs, const Residual &modelResid, const Scalar pvValue, std::vector< Scalar > &B_avg, std::vector< Scalar > &R_sum, std::vector< Scalar > &maxCoeff, std::vector< int > &maxCoeffCell)
Definition: BlackoilModel_impl.hpp:1241
std::unique_ptr< BlackoilModelNldd< TypeTag > > nlddSolver_
Non-linear DD solver.
Definition: BlackoilModel.hpp:375
SimulatorReportSingle nonlinearIterationNewton(const int iteration, const SimulatorTimerInterface &timer, NonlinearSolverType &nonlinear_solver)
Definition: BlackoilModel_impl.hpp:252
void writePartitions(const std::filesystem::path &odir) const
Definition: BlackoilModel_impl.hpp:1214
static constexpr int solventSaturationIdx
Definition: BlackoilModel.hpp:93
bool hasNlddSolver() const
Returns true if an NLDD solver exists.
Definition: BlackoilModel.hpp:338
void beginReportStep()
Definition: BlackoilModel.hpp:312
static constexpr int polymerConcentrationIdx
Definition: BlackoilModel.hpp:95
GetPropType< TypeTag, Properties::Simulator > Simulator
Definition: BlackoilModel.hpp:64
const SimulatorReportSingle & failureReport() const
return the statistics if the nonlinearIteration() method failed
Definition: BlackoilModel.hpp:286
Simulator & simulator()
Definition: BlackoilModel.hpp:282
static constexpr bool has_bioeffects_
Definition: BlackoilModel.hpp:352
void updateTUNINGDP(const TuningDp &tuning_dp)
Definition: BlackoilModel_impl.hpp:839
GetPropType< TypeTag, Properties::MaterialLaw > MaterialLaw
Definition: BlackoilModel.hpp:73
static constexpr int contiPolymerMWEqIdx
Definition: BlackoilModel.hpp:85
static constexpr bool has_micp_
Definition: BlackoilModel.hpp:353
const Simulator & simulator() const
Definition: BlackoilModel.hpp:279
CnvPvSplitData characteriseCnvPvSplit(const std::vector< Scalar > &B_avg, const double dt)
Compute pore-volume/cell count split among "converged", "relaxed converged", "unconverged" cells base...
Definition: BlackoilModel_impl.hpp:745
BlackoilWellModel< TypeTag > & well_model_
Definition: BlackoilModel.hpp:359
BlackoilWellModel< TypeTag > & wellModel()
return the StandardWells object
Definition: BlackoilModel.hpp:305
std::vector< std::vector< Scalar > > residual_norms_history_
Definition: BlackoilModel.hpp:366
typename SparseMatrixAdapter::IstlMatrix Mat
Definition: BlackoilModel.hpp:108
Dune::FieldVector< Scalar, numEq > VectorBlockType
Definition: BlackoilModel.hpp:106
static constexpr bool has_solvent_
Definition: BlackoilModel.hpp:345
const BlackoilWellModel< TypeTag > & wellModel() const
Definition: BlackoilModel.hpp:309
BVector dx_old_
Definition: BlackoilModel.hpp:368
Simulator & simulator_
Definition: BlackoilModel.hpp:343
static constexpr int calciteVolumeFractionIdx
Definition: BlackoilModel.hpp:104
static constexpr bool has_polymer_
Definition: BlackoilModel.hpp:347
void storeSolutionUpdate(const BVector &dx)
Definition: BlackoilModel_impl.hpp:563
static constexpr bool has_brine_
Definition: BlackoilModel.hpp:351
void solveJacobianSystem(BVector &x)
Definition: BlackoilModel_impl.hpp:445
static constexpr int oxygenConcentrationIdx
Definition: BlackoilModel.hpp:101
MaxSolutionUpdateData getMaxSolutionUpdate(const std::vector< unsigned > &ixCells)
Definition: BlackoilModel_impl.hpp:582
ComponentName compNames_
Definition: BlackoilModel.hpp:373
SimulatorReportSingle prepareStep(const SimulatorTimerInterface &timer)
Definition: BlackoilModel_impl.hpp:86
static constexpr int ureaConcentrationIdx
Definition: BlackoilModel.hpp:102
GetPropType< TypeTag, Properties::Grid > Grid
Definition: BlackoilModel.hpp:65
void convergencePerCell(const std::vector< Scalar > &B_avg, const double dt, const double tol_cnv, const double tol_cnv_energy, const int iteration)
Compute the number of Newtons required by each cell in order to satisfy the solution change convergen...
Definition: BlackoilModel_impl.hpp:1099
Dune::BlockVector< VectorBlockType > BVector
Definition: BlackoilModel.hpp:109
void prepareStoringSolutionUpdate()
Get solution update vector as a PrimaryVariable.
Definition: BlackoilModel_impl.hpp:542
GetPropType< TypeTag, Properties::Scalar > Scalar
Definition: BlackoilModel.hpp:75
static constexpr int saltConcentrationIdx
Definition: BlackoilModel.hpp:99
static constexpr int zFractionIdx
Definition: BlackoilModel.hpp:94
static constexpr int contiBrineEqIdx
Definition: BlackoilModel.hpp:87
Scalar relativeChange() const
Definition: BlackoilModel_impl.hpp:358
static constexpr int contiZfracEqIdx
Definition: BlackoilModel.hpp:82
static constexpr int contiSolventEqIdx
Definition: BlackoilModel.hpp:81
GetPropType< TypeTag, Properties::MaterialLawParams > MaterialLawParams
Definition: BlackoilModel.hpp:74
static constexpr int contiPolymerEqIdx
Definition: BlackoilModel.hpp:83
SolutionVector solUpd_
Definition: BlackoilModel.hpp:370
static constexpr bool enableSaltPrecipitation
Definition: BlackoilModel.hpp:78
const EclipseState & eclState() const
Definition: BlackoilModel.hpp:150
bool terminal_output_
Whether we print something to std::cout.
Definition: BlackoilModel.hpp:362
Class for handling the blackoil well model.
Definition: BlackoilWellModel.hpp:98
Definition: ComponentName.hpp:34
Definition: ConvergenceReport.hpp:38
Manages the initializing and running of time dependent problems.
Definition: simulator.hh:84
Interface class for SimulatorTimer objects, to be improved.
Definition: SimulatorTimerInterface.hpp:34
Dune::Communication< MPIComm > Communication
Definition: ParallelCommunication.hpp:30
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
Definition: BlackoilModel.hpp:114
std::vector< unsigned > ixCells
Definition: BlackoilModel.hpp:116
std::pair< std::vector< double >, std::vector< int > > cnvPvSplit
Definition: BlackoilModel.hpp:115
Definition: BlackoilModel.hpp:119
Scalar dRvMax
Definition: BlackoilModel.hpp:123
Scalar dPMax
Definition: BlackoilModel.hpp:120
Scalar dSMax
Definition: BlackoilModel.hpp:121
Scalar dRsMax
Definition: BlackoilModel.hpp:122
Solver parameters for the BlackoilModel.
Definition: BlackoilModelParameters.hpp:194
Scalar max_residual_allowed_
Absolute max limit for residuals.
Definition: BlackoilModelParameters.hpp:207
Definition: SimulatorReport.hpp:122
A struct for returning timing data from a simulator to its caller.
Definition: SimulatorReport.hpp:34