vtkblackoilmodule.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*/
27#ifndef OPM_VTK_BLACK_OIL_MODULE_HPP
28#define OPM_VTK_BLACK_OIL_MODULE_HPP
29
30#include <dune/common/fvector.hh>
31
32#include <opm/material/densead/Math.hpp>
33
35
37
41
44
45namespace Opm {
51template <class TypeTag>
52class VtkBlackOilModule : public BaseOutputModule<TypeTag>
53{
55
60
63
64 static const int vtkFormat = getPropValue<TypeTag, Properties::VtkOutputFormat>();
66
67 enum { oilPhaseIdx = FluidSystem::oilPhaseIdx };
68 enum { gasPhaseIdx = FluidSystem::gasPhaseIdx };
69 enum { waterPhaseIdx = FluidSystem::waterPhaseIdx };
70
71 enum { gasCompIdx = FluidSystem::gasCompIdx };
72 enum { oilCompIdx = FluidSystem::oilCompIdx };
73 enum { waterCompIdx = FluidSystem::waterCompIdx };
74
76
77public:
78 VtkBlackOilModule(const Simulator& simulator)
79 : ParentType(simulator)
80 {
81 params_.read();
82 }
83
88 static void registerParameters()
89 {
91 }
92
98 {
99 if (params_.gasDissolutionFactorOutput_) {
100 this->resizeScalarBuffer_(gasDissolutionFactor_);
101 }
102 if (params_.oilVaporizationFactorOutput_) {
103 this->resizeScalarBuffer_(oilVaporizationFactor_);
104 }
106 this->resizeScalarBuffer_(oilFormationVolumeFactor_);
107 }
109 this->resizeScalarBuffer_(gasFormationVolumeFactor_);
110 }
112 this->resizeScalarBuffer_(waterFormationVolumeFactor_);
113 }
114 if (params_.oilSaturationPressureOutput_) {
115 this->resizeScalarBuffer_(oilSaturationPressure_);
116 }
117 if (params_.gasSaturationPressureOutput_) {
118 this->resizeScalarBuffer_(gasSaturationPressure_);
119 }
121 this->resizeScalarBuffer_(saturatedOilGasDissolutionFactor_);
122 }
124 this->resizeScalarBuffer_(saturatedGasOilVaporizationFactor_);
125 }
126 if (params_.saturationRatiosOutput_) {
127 this->resizeScalarBuffer_(oilSaturationRatio_);
128 this->resizeScalarBuffer_(gasSaturationRatio_);
129 }
130 if (params_.primaryVarsMeaningOutput_) {
131 this->resizeScalarBuffer_(primaryVarsMeaningPressure_);
132 this->resizeScalarBuffer_(primaryVarsMeaningWater_);
133 this->resizeScalarBuffer_(primaryVarsMeaningGas_);
134 }
135 }
136
141 void processElement(const ElementContext& elemCtx)
142 {
143 if (!Parameters::Get<Parameters::EnableVtkOutput>()) {
144 return;
145 }
146
147 for (unsigned dofIdx = 0; dofIdx < elemCtx.numPrimaryDof(/*timeIdx=*/0); ++dofIdx) {
148 const auto& fs = elemCtx.intensiveQuantities(dofIdx, /*timeIdx=*/0).fluidState();
149 using FluidState = typename std::remove_const<typename std::remove_reference<decltype(fs)>::type>::type;
150 unsigned globalDofIdx = elemCtx.globalSpaceIndex(dofIdx, /*timeIdx=*/0);
151
152 const auto& primaryVars = elemCtx.primaryVars(dofIdx, /*timeIdx=*/0);
153
154 unsigned pvtRegionIdx = elemCtx.primaryVars(dofIdx, /*timeIdx=*/0).pvtRegionIndex();
155 Scalar SoMax = 0.0;
156 if (FluidSystem::phaseIsActive(oilPhaseIdx))
157 SoMax = std::max(getValue(fs.saturation(oilPhaseIdx)),
158 elemCtx.problem().maxOilSaturation(globalDofIdx));
159
160 if (FluidSystem::phaseIsActive(gasPhaseIdx) && FluidSystem::phaseIsActive(oilPhaseIdx)) {
161 Scalar x_oG = getValue(fs.moleFraction(oilPhaseIdx, gasCompIdx));
162 Scalar x_gO = getValue(fs.moleFraction(gasPhaseIdx, oilCompIdx));
163 Scalar X_oG = getValue(fs.massFraction(oilPhaseIdx, gasCompIdx));
164 Scalar X_gO = getValue(fs.massFraction(gasPhaseIdx, oilCompIdx));
165 Scalar Rs = FluidSystem::convertXoGToRs(X_oG, pvtRegionIdx);
166 Scalar Rv = FluidSystem::convertXgOToRv(X_gO, pvtRegionIdx);
167
168 Scalar RsSat =
169 FluidSystem::template saturatedDissolutionFactor<FluidState, Scalar>(fs,
170 oilPhaseIdx,
171 pvtRegionIdx,
172 SoMax);
173 Scalar X_oG_sat = FluidSystem::convertRsToXoG(RsSat, pvtRegionIdx);
174 Scalar x_oG_sat = FluidSystem::convertXoGToxoG(X_oG_sat, pvtRegionIdx);
175
176 Scalar RvSat =
177 FluidSystem::template saturatedDissolutionFactor<FluidState, Scalar>(fs,
178 gasPhaseIdx,
179 pvtRegionIdx,
180 SoMax);
181 Scalar X_gO_sat = FluidSystem::convertRvToXgO(RvSat, pvtRegionIdx);
182 Scalar x_gO_sat = FluidSystem::convertXgOToxgO(X_gO_sat, pvtRegionIdx);
183 if (params_.gasDissolutionFactorOutput_) {
184 gasDissolutionFactor_[globalDofIdx] = Rs;
185 }
186 if (params_.oilVaporizationFactorOutput_) {
187 oilVaporizationFactor_[globalDofIdx] = Rv;
188 }
189 if (params_.oilSaturationPressureOutput_) {
190 oilSaturationPressure_[globalDofIdx] =
191 FluidSystem::template saturationPressure<FluidState, Scalar>(fs, oilPhaseIdx, pvtRegionIdx);
192 }
193 if (params_.gasSaturationPressureOutput_) {
194 gasSaturationPressure_[globalDofIdx] =
195 FluidSystem::template saturationPressure<FluidState, Scalar>(fs, gasPhaseIdx, pvtRegionIdx);
196 }
198 saturatedOilGasDissolutionFactor_[globalDofIdx] = RsSat;
199 }
201 saturatedGasOilVaporizationFactor_[globalDofIdx] = RvSat;
202 }
203 if (params_.saturationRatiosOutput_) {
204 if (x_oG_sat <= 0.0) {
205 oilSaturationRatio_[globalDofIdx] = 1.0;
206 }
207 else {
208 oilSaturationRatio_[globalDofIdx] = x_oG / x_oG_sat;
209 }
210
211 if (x_gO_sat <= 0.0) {
212 gasSaturationRatio_[globalDofIdx] = 1.0;
213 }
214 else {
215 gasSaturationRatio_[globalDofIdx] = x_gO / x_gO_sat;
216 }
217 }
218 }
220 oilFormationVolumeFactor_[globalDofIdx] =
221 1.0 / FluidSystem::template inverseFormationVolumeFactor<FluidState, Scalar>(fs, oilPhaseIdx, pvtRegionIdx);
222 }
224 gasFormationVolumeFactor_[globalDofIdx] =
225 1.0 / FluidSystem::template inverseFormationVolumeFactor<FluidState, Scalar>(fs, gasPhaseIdx, pvtRegionIdx);
226 }
228 waterFormationVolumeFactor_[globalDofIdx] =
229 1.0 / FluidSystem::template inverseFormationVolumeFactor<FluidState, Scalar>(fs, waterPhaseIdx, pvtRegionIdx);
230 }
231
232 if (params_.primaryVarsMeaningOutput_) {
233 primaryVarsMeaningWater_[globalDofIdx] =
234 static_cast<int>(primaryVars.primaryVarsMeaningWater());
235 primaryVarsMeaningGas_[globalDofIdx] =
236 static_cast<int>(primaryVars.primaryVarsMeaningGas());
237 primaryVarsMeaningPressure_[globalDofIdx] =
238 static_cast<int>(primaryVars.primaryVarsMeaningPressure());
239 }
240 }
241 }
242
247 {
248 VtkMultiWriter *vtkWriter = dynamic_cast<VtkMultiWriter*>(&baseWriter);
249 if (!vtkWriter) {
250 return;
251 }
252
253 if (params_.gasDissolutionFactorOutput_) {
254 this->commitScalarBuffer_(baseWriter, "R_s", gasDissolutionFactor_);
255 }
256 if (params_.oilVaporizationFactorOutput_) {
257 this->commitScalarBuffer_(baseWriter, "R_v", oilVaporizationFactor_);
258 }
260 this->commitScalarBuffer_(baseWriter, "B_o", oilFormationVolumeFactor_);
261 }
263 this->commitScalarBuffer_(baseWriter, "B_g", gasFormationVolumeFactor_);
264 }
266 this->commitScalarBuffer_(baseWriter, "B_w", waterFormationVolumeFactor_);
267 }
268 if (params_.oilSaturationPressureOutput_) {
269 this->commitScalarBuffer_(baseWriter, "p_o,sat", oilSaturationPressure_);
270 }
271 if (params_.gasSaturationPressureOutput_) {
272 this->commitScalarBuffer_(baseWriter, "p_g,sat", gasSaturationPressure_);
273 }
275 this->commitScalarBuffer_(baseWriter, "R_s,sat", saturatedOilGasDissolutionFactor_);
276 }
278 this->commitScalarBuffer_(baseWriter, "R_v,sat", saturatedGasOilVaporizationFactor_);
279 }
280 if (params_.saturationRatiosOutput_) {
281 this->commitScalarBuffer_(baseWriter, "saturation ratio_oil", oilSaturationRatio_);
282 this->commitScalarBuffer_(baseWriter, "saturation ratio_gas", gasSaturationRatio_);
283 }
284
285 if (params_.primaryVarsMeaningOutput_) {
286 this->commitScalarBuffer_(baseWriter, "primary vars meaning water", primaryVarsMeaningWater_);
287 this->commitScalarBuffer_(baseWriter, "primary vars meaning gas", primaryVarsMeaningGas_);
288 this->commitScalarBuffer_(baseWriter, "primary vars meaning pressure", primaryVarsMeaningPressure_);
289 }
290 }
291
292private:
293 VtkBlackoilParams params_{};
294 ScalarBuffer gasDissolutionFactor_{};
295 ScalarBuffer oilVaporizationFactor_{};
296 ScalarBuffer oilFormationVolumeFactor_{};
297 ScalarBuffer gasFormationVolumeFactor_{};
298 ScalarBuffer waterFormationVolumeFactor_{};
299 ScalarBuffer oilSaturationPressure_{};
300 ScalarBuffer gasSaturationPressure_{};
301
302 ScalarBuffer saturatedOilGasDissolutionFactor_{};
303 ScalarBuffer saturatedGasOilVaporizationFactor_{};
304 ScalarBuffer oilSaturationRatio_{};
305 ScalarBuffer gasSaturationRatio_{};
306
307 ScalarBuffer primaryVarsMeaningPressure_{};
308 ScalarBuffer primaryVarsMeaningWater_{};
309 ScalarBuffer primaryVarsMeaningGas_{};
310};
311
312} // namespace Opm
313
314#endif // OPM_VTK_BLACK_OIL_MODULE_HPP
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.hpp:53
void processElement(const ElementContext &elemCtx)
Modify the internal buffers according to the intensive quantities relevant for an element.
Definition: vtkblackoilmodule.hpp:141
static void registerParameters()
Register all run-time parameters for the multi-phase VTK output module.
Definition: vtkblackoilmodule.hpp:88
void commitBuffers(BaseOutputWriter &baseWriter)
Add all buffers to the VTK output writer.
Definition: vtkblackoilmodule.hpp:246
VtkBlackOilModule(const Simulator &simulator)
Definition: vtkblackoilmodule.hpp:78
void allocBuffers()
Allocate memory for the scalar fields we would like to write to the VTK file.
Definition: vtkblackoilmodule.hpp:97
Simplifies writing multi-file VTK datasets.
Definition: vtkmultiwriter.hh:66
Declare the properties used by the infrastructure code of the finite volume discretizations.
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.
Struct holding the parameters for VtkBlackoilOutputModule.
Definition: vtkblackoilparams.hpp:53
bool gasSaturationPressureOutput_
Definition: vtkblackoilparams.hpp:66
bool oilFormationVolumeFactorOutput_
Definition: vtkblackoilparams.hpp:62
static void registerParameters()
Registers the parameters in parameter system.
bool gasFormationVolumeFactorOutput_
Definition: vtkblackoilparams.hpp:63
bool oilVaporizationFactorOutput_
Definition: vtkblackoilparams.hpp:61
bool saturatedGasOilVaporizationFactorOutput_
Definition: vtkblackoilparams.hpp:68
bool primaryVarsMeaningOutput_
Definition: vtkblackoilparams.hpp:70
bool gasDissolutionFactorOutput_
Definition: vtkblackoilparams.hpp:60
bool waterFormationVolumeFactorOutput_
Definition: vtkblackoilparams.hpp:64
bool saturatedOilGasDissolutionFactorOutput_
Definition: vtkblackoilparams.hpp:67
bool saturationRatiosOutput_
Definition: vtkblackoilparams.hpp:69
bool oilSaturationPressureOutput_
Definition: vtkblackoilparams.hpp:65
void read()
Reads the parameter values from the parameter system.