27#ifndef EWOMS_VTK_BLACK_OIL_MODULE_HH
28#define EWOMS_VTK_BLACK_OIL_MODULE_HH
33#include <dune/common/fvector.hh>
35#include <opm/material/densead/Math.hpp>
67template <
class TypeTag>
80 static const int vtkFormat = getPropValue<TypeTag, Properties::VtkOutputFormat>();
83 enum { oilPhaseIdx = FluidSystem::oilPhaseIdx };
84 enum { gasPhaseIdx = FluidSystem::gasPhaseIdx };
85 enum { waterPhaseIdx = FluidSystem::waterPhaseIdx };
87 enum { gasCompIdx = FluidSystem::gasCompIdx };
88 enum { oilCompIdx = FluidSystem::oilCompIdx };
89 enum { waterCompIdx = FluidSystem::waterCompIdx };
104 Parameters::Register<Parameters::VtkWriteGasDissolutionFactor>
105 (
"Include the gas dissolution factor (R_s) of the observed oil "
106 "in the VTK output files");
107 Parameters::Register<Parameters::VtkWriteOilVaporizationFactor>
108 (
"Include the oil vaporization factor (R_v) of the observed gas "
109 "in the VTK output files");
110 Parameters::Register<Parameters::VtkWriteOilFormationVolumeFactor>
111 (
"Include the oil formation volume factor (B_o) in the VTK output files");
112 Parameters::Register<Parameters::VtkWriteGasFormationVolumeFactor>
113 (
"Include the gas formation volume factor (B_g) in the "
115 Parameters::Register<Parameters::VtkWriteWaterFormationVolumeFactor>
116 (
"Include the water formation volume factor (B_w) in the "
118 Parameters::Register<Parameters::VtkWriteOilSaturationPressure>
119 (
"Include the saturation pressure of oil (p_o,sat) in the "
121 Parameters::Register<Parameters::VtkWriteGasSaturationPressure>
122 (
"Include the saturation pressure of gas (p_g,sat) in the "
124 Parameters::Register<Parameters::VtkWriteSaturatedOilGasDissolutionFactor>
125 (
"Include the gas dissolution factor (R_s,sat) of gas saturated "
126 "oil in the VTK output files");
127 Parameters::Register<Parameters::VtkWriteSaturatedGasOilVaporizationFactor>
128 (
"Include the oil vaporization factor (R_v,sat) of oil saturated "
129 "gas in the VTK output files");
130 Parameters::Register<Parameters::VtkWriteSaturationRatios>
131 (
"Write the ratio of the actually and maximum dissolved component of "
133 Parameters::Register<Parameters::VtkWritePrimaryVarsMeaning>
134 (
"Include how the primary variables should be interpreted to the "
144 if (gasDissolutionFactorOutput_())
146 if (oilVaporizationFactorOutput_())
148 if (oilFormationVolumeFactorOutput_())
150 if (gasFormationVolumeFactorOutput_())
152 if (waterFormationVolumeFactorOutput_())
154 if (oilSaturationPressureOutput_())
156 if (gasSaturationPressureOutput_())
158 if (saturatedOilGasDissolutionFactorOutput_())
160 if (saturatedGasOilVaporizationFactorOutput_())
162 if (saturationRatiosOutput_()) {
166 if (primaryVarsMeaningOutput_()) {
179 if (!Parameters::Get<Parameters::EnableVtkOutput>()) {
183 for (
unsigned dofIdx = 0; dofIdx < elemCtx.numPrimaryDof(0); ++dofIdx) {
184 const auto& fs = elemCtx.intensiveQuantities(dofIdx, 0).fluidState();
185 using FluidState =
typename std::remove_const<
typename std::remove_reference<
decltype(fs)>::type>::type;
186 unsigned globalDofIdx = elemCtx.globalSpaceIndex(dofIdx, 0);
188 const auto& primaryVars = elemCtx.primaryVars(dofIdx, 0);
190 unsigned pvtRegionIdx = elemCtx.primaryVars(dofIdx, 0).pvtRegionIndex();
192 if (FluidSystem::phaseIsActive(oilPhaseIdx))
193 SoMax = std::max(getValue(fs.saturation(oilPhaseIdx)),
194 elemCtx.problem().maxOilSaturation(globalDofIdx));
196 if (FluidSystem::phaseIsActive(gasPhaseIdx) && FluidSystem::phaseIsActive(oilPhaseIdx)) {
197 Scalar x_oG = getValue(fs.moleFraction(oilPhaseIdx, gasCompIdx));
198 Scalar x_gO = getValue(fs.moleFraction(gasPhaseIdx, oilCompIdx));
199 Scalar X_oG = getValue(fs.massFraction(oilPhaseIdx, gasCompIdx));
200 Scalar X_gO = getValue(fs.massFraction(gasPhaseIdx, oilCompIdx));
201 Scalar Rs = FluidSystem::convertXoGToRs(X_oG, pvtRegionIdx);
202 Scalar Rv = FluidSystem::convertXgOToRv(X_gO, pvtRegionIdx);
205 FluidSystem::template saturatedDissolutionFactor<FluidState, Scalar>(fs,
209 Scalar X_oG_sat = FluidSystem::convertRsToXoG(RsSat, pvtRegionIdx);
210 Scalar x_oG_sat = FluidSystem::convertXoGToxoG(X_oG_sat, pvtRegionIdx);
213 FluidSystem::template saturatedDissolutionFactor<FluidState, Scalar>(fs,
217 Scalar X_gO_sat = FluidSystem::convertRvToXgO(RvSat, pvtRegionIdx);
218 Scalar x_gO_sat = FluidSystem::convertXgOToxgO(X_gO_sat, pvtRegionIdx);
219 if (gasDissolutionFactorOutput_())
220 gasDissolutionFactor_[globalDofIdx] = Rs;
221 if (oilVaporizationFactorOutput_())
222 oilVaporizationFactor_[globalDofIdx] = Rv;
223 if (oilSaturationPressureOutput_())
224 oilSaturationPressure_[globalDofIdx] =
225 FluidSystem::template saturationPressure<FluidState, Scalar>(fs, oilPhaseIdx, pvtRegionIdx);
226 if (gasSaturationPressureOutput_())
227 gasSaturationPressure_[globalDofIdx] =
228 FluidSystem::template saturationPressure<FluidState, Scalar>(fs, gasPhaseIdx, pvtRegionIdx);
229 if (saturatedOilGasDissolutionFactorOutput_())
230 saturatedOilGasDissolutionFactor_[globalDofIdx] = RsSat;
231 if (saturatedGasOilVaporizationFactorOutput_())
232 saturatedGasOilVaporizationFactor_[globalDofIdx] = RvSat;
233 if (saturationRatiosOutput_()) {
235 oilSaturationRatio_[globalDofIdx] = 1.0;
237 oilSaturationRatio_[globalDofIdx] = x_oG / x_oG_sat;
240 gasSaturationRatio_[globalDofIdx] = 1.0;
242 gasSaturationRatio_[globalDofIdx] = x_gO / x_gO_sat;
245 if (oilFormationVolumeFactorOutput_())
246 oilFormationVolumeFactor_[globalDofIdx] =
247 1.0/FluidSystem::template inverseFormationVolumeFactor<FluidState, Scalar>(fs, oilPhaseIdx, pvtRegionIdx);
248 if (gasFormationVolumeFactorOutput_())
249 gasFormationVolumeFactor_[globalDofIdx] =
250 1.0/FluidSystem::template inverseFormationVolumeFactor<FluidState, Scalar>(fs, gasPhaseIdx, pvtRegionIdx);
251 if (waterFormationVolumeFactorOutput_())
252 waterFormationVolumeFactor_[globalDofIdx] =
253 1.0/FluidSystem::template inverseFormationVolumeFactor<FluidState, Scalar>(fs, waterPhaseIdx, pvtRegionIdx);
255 if (primaryVarsMeaningOutput_()) {
256 primaryVarsMeaningWater_[globalDofIdx] =
257 static_cast<int>(primaryVars.primaryVarsMeaningWater());
258 primaryVarsMeaningGas_[globalDofIdx] =
259 static_cast<int>(primaryVars.primaryVarsMeaningGas());
260 primaryVarsMeaningPressure_[globalDofIdx] =
261 static_cast<int>(primaryVars.primaryVarsMeaningPressure());
275 if (gasDissolutionFactorOutput_())
277 if (oilVaporizationFactorOutput_())
279 if (oilFormationVolumeFactorOutput_())
281 if (gasFormationVolumeFactorOutput_())
283 if (waterFormationVolumeFactorOutput_())
285 if (oilSaturationPressureOutput_())
287 if (gasSaturationPressureOutput_())
289 if (saturatedOilGasDissolutionFactorOutput_())
291 if (saturatedGasOilVaporizationFactorOutput_())
293 if (saturationRatiosOutput_()) {
298 if (primaryVarsMeaningOutput_()) {
299 this->
commitScalarBuffer_(baseWriter,
"primary vars meaning water", primaryVarsMeaningWater_);
301 this->
commitScalarBuffer_(baseWriter,
"primary vars meaning pressure", primaryVarsMeaningPressure_);
306 static bool gasDissolutionFactorOutput_()
308 static bool val = Parameters::Get<Parameters::VtkWriteGasDissolutionFactor>();
312 static bool oilVaporizationFactorOutput_()
314 static bool val = Parameters::Get<Parameters::VtkWriteOilVaporizationFactor>();
318 static bool oilFormationVolumeFactorOutput_()
320 static bool val = Parameters::Get<Parameters::VtkWriteOilFormationVolumeFactor>();
324 static bool gasFormationVolumeFactorOutput_()
326 static bool val = Parameters::Get<Parameters::VtkWriteGasFormationVolumeFactor>();
330 static bool waterFormationVolumeFactorOutput_()
332 static bool val = Parameters::Get<Parameters::VtkWriteWaterFormationVolumeFactor>();
336 static bool oilSaturationPressureOutput_()
338 static bool val = Parameters::Get<Parameters::VtkWriteOilSaturationPressure>();
342 static bool gasSaturationPressureOutput_()
344 static bool val = Parameters::Get<Parameters::VtkWriteGasSaturationPressure>();
348 static bool saturatedOilGasDissolutionFactorOutput_()
350 static bool val = Parameters::Get<Parameters::VtkWriteSaturatedOilGasDissolutionFactor>();
354 static bool saturatedGasOilVaporizationFactorOutput_()
356 static bool val = Parameters::Get<Parameters::VtkWriteSaturatedGasOilVaporizationFactor>();
360 static bool saturationRatiosOutput_()
362 static bool val = Parameters::Get<Parameters::VtkWriteSaturationRatios>();
366 static bool primaryVarsMeaningOutput_()
368 static bool val = Parameters::Get<Parameters::VtkWritePrimaryVarsMeaning>();
372 ScalarBuffer gasDissolutionFactor_;
373 ScalarBuffer oilVaporizationFactor_;
374 ScalarBuffer oilFormationVolumeFactor_;
375 ScalarBuffer gasFormationVolumeFactor_;
376 ScalarBuffer waterFormationVolumeFactor_;
377 ScalarBuffer oilSaturationPressure_;
378 ScalarBuffer gasSaturationPressure_;
380 ScalarBuffer saturatedOilGasDissolutionFactor_;
381 ScalarBuffer saturatedGasOilVaporizationFactor_;
382 ScalarBuffer oilSaturationRatio_;
383 ScalarBuffer gasSaturationRatio_;
385 ScalarBuffer primaryVarsMeaningPressure_;
386 ScalarBuffer primaryVarsMeaningWater_;
387 ScalarBuffer primaryVarsMeaningGas_;
Declares the properties required by the black oil model.
The base class for writer modules.
Definition: baseoutputmodule.hh:67
BaseOutputWriter::ScalarBuffer ScalarBuffer
Definition: baseoutputmodule.hh:85
void commitScalarBuffer_(BaseOutputWriter &baseWriter, const char *name, ScalarBuffer &buffer, BufferType bufferType=DofBuffer)
Add a buffer containing scalar quantities to the result file.
Definition: baseoutputmodule.hh:244
void resizeScalarBuffer_(ScalarBuffer &buffer, BufferType bufferType=DofBuffer)
Allocate the space for a buffer storing a scalar quantity.
Definition: baseoutputmodule.hh:156
The base class for all output writers.
Definition: baseoutputwriter.hh:44
VTK output module for the black oil model's parameters.
Definition: vtkblackoilmodule.hh:69
void processElement(const ElementContext &elemCtx)
Modify the internal buffers according to the intensive quantities relevant for an element.
Definition: vtkblackoilmodule.hh:177
static void registerParameters()
Register all run-time parameters for the multi-phase VTK output module.
Definition: vtkblackoilmodule.hh:102
void commitBuffers(BaseOutputWriter &baseWriter)
Add all buffers to the VTK output writer.
Definition: vtkblackoilmodule.hh:269
VtkBlackOilModule(const Simulator &simulator)
Definition: vtkblackoilmodule.hh:94
void allocBuffers()
Allocate memory for the scalar fields we would like to write to the VTK file.
Definition: vtkblackoilmodule.hh:142
Simplifies writing multi-file VTK datasets.
Definition: vtkmultiwriter.hh:66
Declare the properties used by the infrastructure code of the finite volume discretizations.
Definition: blackoilnewtonmethodparameters.hh:31
Definition: blackoilboundaryratevector.hh:37
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:235
This file provides the infrastructure to retrieve run-time parameters.
The Opm property system, traits with inheritance.
Definition: vtkblackoilmodule.hh:47
static constexpr bool value
Definition: vtkblackoilmodule.hh:47
Definition: vtkblackoilmodule.hh:53
static constexpr bool value
Definition: vtkblackoilmodule.hh:53
Definition: vtkblackoilmodule.hh:52
static constexpr bool value
Definition: vtkblackoilmodule.hh:52
Definition: vtkblackoilmodule.hh:48
static constexpr bool value
Definition: vtkblackoilmodule.hh:48
Definition: vtkblackoilmodule.hh:57
static constexpr bool value
Definition: vtkblackoilmodule.hh:57
Definition: vtkblackoilmodule.hh:56
static constexpr bool value
Definition: vtkblackoilmodule.hh:56
Definition: vtkblackoilmodule.hh:55
static constexpr bool value
Definition: vtkblackoilmodule.hh:55
Definition: vtkblackoilmodule.hh:54
static constexpr bool value
Definition: vtkblackoilmodule.hh:54