vtkblackoilpolymermodule.hh
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 EWOMS_VTK_BLACK_OIL_POLYMER_MODULE_HH
28#define EWOMS_VTK_BLACK_OIL_POLYMER_MODULE_HH
29
30#include <opm/material/densead/Math.hpp>
31
32#include "vtkmultiwriter.hh"
33#include "baseoutputmodule.hh"
34
38
39#include <dune/common/fvector.hh>
40
41#include <cstdio>
42
43namespace Opm::Properties {
44
45namespace TTag {
46
47// create new type tag for the VTK multi-phase output
49
50} // namespace TTag
51
52// create the property tags needed for the polymer output module
53template<class TypeTag, class MyTypeTag>
55template<class TypeTag, class MyTypeTag>
57template<class TypeTag, class MyTypeTag>
59template<class TypeTag, class MyTypeTag>
61template<class TypeTag, class MyTypeTag>
63template<class TypeTag, class MyTypeTag>
65
66// set default values for what quantities to output
67template<class TypeTag>
68struct VtkWritePolymerConcentration<TypeTag, TTag::VtkBlackOilPolymer> { static constexpr bool value = true; };
69template<class TypeTag>
70struct VtkWritePolymerDeadPoreVolume<TypeTag, TTag::VtkBlackOilPolymer> { static constexpr bool value = true; };
71template<class TypeTag>
72struct VtkWritePolymerViscosityCorrection<TypeTag, TTag::VtkBlackOilPolymer> { static constexpr bool value = true; };
73template<class TypeTag>
74struct VtkWriteWaterViscosityCorrection<TypeTag, TTag::VtkBlackOilPolymer> { static constexpr bool value = true; };
75template<class TypeTag>
76struct VtkWritePolymerRockDensity<TypeTag, TTag::VtkBlackOilPolymer> { static constexpr bool value = true; };
77template<class TypeTag>
78struct VtkWritePolymerAdsorption<TypeTag, TTag::VtkBlackOilPolymer> { static constexpr bool value = true; };
79
80} // namespace Opm::Properties
81
82namespace Opm {
88template <class TypeTag>
90{
92
98
99 static const int vtkFormat = getPropValue<TypeTag, Properties::VtkOutputFormat>();
101
102 enum { enablePolymer = getPropValue<TypeTag, Properties::EnablePolymer>() };
103
104 using ScalarBuffer = typename ParentType::ScalarBuffer;
105
106public:
107 VtkBlackOilPolymerModule(const Simulator& simulator)
108 : ParentType(simulator)
109 { }
110
115 static void registerParameters()
116 {
117 if (!enablePolymer)
118 return;
119
120 Parameters::registerParam<TypeTag, Properties::VtkWritePolymerConcentration>
121 ("Include the concentration of the polymer component in the water phase "
122 "in the VTK output files");
123 Parameters::registerParam<TypeTag, Properties::VtkWritePolymerDeadPoreVolume>
124 ("Include the fraction of the \"dead\" pore volume "
125 "in the VTK output files");
126 Parameters::registerParam<TypeTag, Properties::VtkWritePolymerRockDensity>
127 ("Include the amount of already adsorbed polymer component"
128 "in the VTK output files");
129 Parameters::registerParam<TypeTag, Properties::VtkWritePolymerAdsorption>
130 ("Include the adsorption rate of the polymer component"
131 "in the VTK output files");
132 Parameters::registerParam<TypeTag, Properties::VtkWritePolymerViscosityCorrection>
133 ("Include the viscosity correction of the polymer component "
134 "in the VTK output files");
135 Parameters::registerParam<TypeTag, Properties::VtkWriteWaterViscosityCorrection>
136 ("Include the viscosity correction of the water component "
137 "due to polymers in the VTK output files");
138 }
139
145 {
146 if (!Parameters::get<TypeTag, Properties::EnableVtkOutput>())
147 return;
148
149 if (!enablePolymer)
150 return;
151
152 if (polymerConcentrationOutput_())
153 this->resizeScalarBuffer_(polymerConcentration_);
154 if (polymerDeadPoreVolumeOutput_())
155 this->resizeScalarBuffer_(polymerDeadPoreVolume_);
156 if (polymerRockDensityOutput_())
157 this->resizeScalarBuffer_(polymerRockDensity_);
158 if (polymerAdsorptionOutput_())
159 this->resizeScalarBuffer_(polymerAdsorption_);
160 if (polymerViscosityCorrectionOutput_())
161 this->resizeScalarBuffer_(polymerViscosityCorrection_);
162 if (waterViscosityCorrectionOutput_())
163 this->resizeScalarBuffer_(waterViscosityCorrection_);
164 }
165
170 void processElement(const ElementContext& elemCtx)
171 {
172 if (!Parameters::get<TypeTag, Properties::EnableVtkOutput>())
173 return;
174
175 if (!enablePolymer)
176 return;
177
178 for (unsigned dofIdx = 0; dofIdx < elemCtx.numPrimaryDof(/*timeIdx=*/0); ++dofIdx) {
179 const auto& intQuants = elemCtx.intensiveQuantities(dofIdx, /*timeIdx=*/0);
180 unsigned globalDofIdx = elemCtx.globalSpaceIndex(dofIdx, /*timeIdx=*/0);
181
182 if (polymerConcentrationOutput_())
183 polymerConcentration_[globalDofIdx] =
184 scalarValue(intQuants.polymerConcentration());
185
186 if (polymerDeadPoreVolumeOutput_())
187 polymerDeadPoreVolume_[globalDofIdx] =
188 scalarValue(intQuants.polymerDeadPoreVolume());
189
190 if (polymerRockDensityOutput_())
191 polymerRockDensity_[globalDofIdx] =
192 scalarValue(intQuants.polymerRockDensity());
193
194 if (polymerAdsorptionOutput_())
195 polymerAdsorption_[globalDofIdx] =
196 scalarValue(intQuants.polymerAdsorption());
197
198 if (polymerViscosityCorrectionOutput_())
199 polymerViscosityCorrection_[globalDofIdx] =
200 scalarValue(intQuants.polymerViscosityCorrection());
201
202 if (waterViscosityCorrectionOutput_())
203 waterViscosityCorrection_[globalDofIdx] =
204 scalarValue(intQuants.waterViscosityCorrection());
205 }
206 }
207
212 {
213 VtkMultiWriter *vtkWriter = dynamic_cast<VtkMultiWriter*>(&baseWriter);
214 if (!vtkWriter)
215 return;
216
217 if (!enablePolymer)
218 return;
219
220 if (polymerConcentrationOutput_())
221 this->commitScalarBuffer_(baseWriter, "polymer concentration", polymerConcentration_);
222
223 if (polymerDeadPoreVolumeOutput_())
224 this->commitScalarBuffer_(baseWriter, "dead pore volume fraction", polymerDeadPoreVolume_);
225
226 if (polymerRockDensityOutput_())
227 this->commitScalarBuffer_(baseWriter, "polymer rock density", polymerRockDensity_);
228
229 if (polymerAdsorptionOutput_())
230 this->commitScalarBuffer_(baseWriter, "polymer adsorption", polymerAdsorption_);
231
232 if (polymerViscosityCorrectionOutput_())
233 this->commitScalarBuffer_(baseWriter, "polymer viscosity correction", polymerViscosityCorrection_);
234
235 if (waterViscosityCorrectionOutput_())
236 this->commitScalarBuffer_(baseWriter, "water viscosity correction", waterViscosityCorrection_);
237 }
238
239private:
240 static bool polymerConcentrationOutput_()
241 {
242 static bool val = Parameters::get<TypeTag, Properties::VtkWritePolymerConcentration>();
243 return val;
244 }
245
246 static bool polymerDeadPoreVolumeOutput_()
247 {
248 static bool val = Parameters::get<TypeTag, Properties::VtkWritePolymerDeadPoreVolume>();
249 return val;
250 }
251
252 static bool polymerRockDensityOutput_()
253 {
254 static bool val = Parameters::get<TypeTag, Properties::VtkWritePolymerRockDensity>();
255 return val;
256 }
257
258 static bool polymerAdsorptionOutput_()
259 {
260 static bool val = Parameters::get<TypeTag, Properties::VtkWritePolymerAdsorption>();
261 return val;
262 }
263
264 static bool polymerViscosityCorrectionOutput_()
265 {
266 static bool val = Parameters::get<TypeTag, Properties::VtkWritePolymerViscosityCorrection>();
267 return val;
268 }
269
270 static bool waterViscosityCorrectionOutput_()
271 {
272 static bool val = Parameters::get<TypeTag, Properties::VtkWritePolymerViscosityCorrection>();
273 return val;
274 }
275
276 ScalarBuffer polymerConcentration_;
277 ScalarBuffer polymerDeadPoreVolume_;
278 ScalarBuffer polymerRockDensity_;
279 ScalarBuffer polymerAdsorption_;
280 ScalarBuffer polymerViscosityCorrection_;
281 ScalarBuffer waterViscosityCorrection_;
282};
283} // namespace Opm
284
285#endif
Declares the properties required by the black oil model.
The base class for writer modules.
Definition: baseoutputmodule.hh:73
BaseOutputWriter::ScalarBuffer ScalarBuffer
Definition: baseoutputmodule.hh:91
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:316
void resizeScalarBuffer_(ScalarBuffer &buffer, BufferType bufferType=DofBuffer)
Allocate the space for a buffer storing a scalar quantity.
Definition: baseoutputmodule.hh:162
The base class for all output writers.
Definition: baseoutputwriter.hh:44
VTK output module for the black oil model's polymer related quantities.
Definition: vtkblackoilpolymermodule.hh:90
void allocBuffers()
Allocate memory for the scalar fields we would like to write to the VTK file.
Definition: vtkblackoilpolymermodule.hh:144
static void registerParameters()
Register all run-time parameters for the multi-phase VTK output module.
Definition: vtkblackoilpolymermodule.hh:115
void processElement(const ElementContext &elemCtx)
Modify the internal buffers according to the intensive quantities relevant for an element.
Definition: vtkblackoilpolymermodule.hh:170
VtkBlackOilPolymerModule(const Simulator &simulator)
Definition: vtkblackoilpolymermodule.hh:107
void commitBuffers(BaseOutputWriter &baseWriter)
Add all buffers to the VTK output writer.
Definition: vtkblackoilpolymermodule.hh:211
Simplifies writing multi-file VTK datasets.
Definition: vtkmultiwriter.hh:66
Definition: blackoilmodel.hh:72
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:242
This file provides the infrastructure to retrieve run-time parameters.
The Opm property system, traits with inheritance.
Definition: vtkblackoilpolymermodule.hh:48
a tag to mark properties as undefined
Definition: propertysystem.hh:40
Definition: vtkblackoilpolymermodule.hh:58
Definition: vtkblackoilpolymermodule.hh:54
Definition: vtkblackoilpolymermodule.hh:56
Definition: vtkblackoilpolymermodule.hh:60
Definition: vtkblackoilpolymermodule.hh:62
Definition: vtkblackoilpolymermodule.hh:64