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 <dune/common/fvector.hh>
31
32#include <opm/material/densead/Math.hpp>
33
35
37
40
43
44namespace Opm::Properties::TTag {
45
46// create new type tag for the VTK multi-phase output
48
49} // namespace Opm::Properties::TTag
50
51namespace Opm::Parameters {
52
53// set default values for what quantities to output
54struct VtkWritePolymerConcentration { static constexpr bool value = true; };
55struct VtkWritePolymerDeadPoreVolume { static constexpr bool value = true; };
56struct VtkWritePolymerViscosityCorrection { static constexpr bool value = true; };
57struct VtkWriteWaterViscosityCorrection { static constexpr bool value = true; };
58struct VtkWritePolymerRockDensity { static constexpr bool value = true; };
59struct VtkWritePolymerAdsorption { static constexpr bool value = true; };
60
61} // namespace Opm::Parameters
62
63namespace Opm {
69template <class TypeTag>
71{
73
79
80 static const int vtkFormat = getPropValue<TypeTag, Properties::VtkOutputFormat>();
82
83 enum { enablePolymer = getPropValue<TypeTag, Properties::EnablePolymer>() };
84
86
87public:
88 VtkBlackOilPolymerModule(const Simulator& simulator)
89 : ParentType(simulator)
90 { }
91
96 static void registerParameters()
97 {
98 if (!enablePolymer)
99 return;
100
101 Parameters::Register<Parameters::VtkWritePolymerConcentration>
102 ("Include the concentration of the polymer component in the water phase "
103 "in the VTK output files");
104 Parameters::Register<Parameters::VtkWritePolymerDeadPoreVolume>
105 ("Include the fraction of the \"dead\" pore volume "
106 "in the VTK output files");
107 Parameters::Register<Parameters::VtkWritePolymerRockDensity>
108 ("Include the amount of already adsorbed polymer component"
109 "in the VTK output files");
110 Parameters::Register<Parameters::VtkWritePolymerAdsorption>
111 ("Include the adsorption rate of the polymer component"
112 "in the VTK output files");
113 Parameters::Register<Parameters::VtkWritePolymerViscosityCorrection>
114 ("Include the viscosity correction of the polymer component "
115 "in the VTK output files");
116 Parameters::Register<Parameters::VtkWriteWaterViscosityCorrection>
117 ("Include the viscosity correction of the water component "
118 "due to polymers in the VTK output files");
119 }
120
126 {
127 if (!Parameters::Get<Parameters::EnableVtkOutput>()) {
128 return;
129 }
130
131 if (!enablePolymer)
132 return;
133
134 if (polymerConcentrationOutput_())
135 this->resizeScalarBuffer_(polymerConcentration_);
136 if (polymerDeadPoreVolumeOutput_())
137 this->resizeScalarBuffer_(polymerDeadPoreVolume_);
138 if (polymerRockDensityOutput_())
139 this->resizeScalarBuffer_(polymerRockDensity_);
140 if (polymerAdsorptionOutput_())
141 this->resizeScalarBuffer_(polymerAdsorption_);
142 if (polymerViscosityCorrectionOutput_())
143 this->resizeScalarBuffer_(polymerViscosityCorrection_);
144 if (waterViscosityCorrectionOutput_())
145 this->resizeScalarBuffer_(waterViscosityCorrection_);
146 }
147
152 void processElement(const ElementContext& elemCtx)
153 {
154 if (!Parameters::Get<Parameters::EnableVtkOutput>()) {
155 return;
156 }
157
158 if (!enablePolymer)
159 return;
160
161 for (unsigned dofIdx = 0; dofIdx < elemCtx.numPrimaryDof(/*timeIdx=*/0); ++dofIdx) {
162 const auto& intQuants = elemCtx.intensiveQuantities(dofIdx, /*timeIdx=*/0);
163 unsigned globalDofIdx = elemCtx.globalSpaceIndex(dofIdx, /*timeIdx=*/0);
164
165 if (polymerConcentrationOutput_())
166 polymerConcentration_[globalDofIdx] =
167 scalarValue(intQuants.polymerConcentration());
168
169 if (polymerDeadPoreVolumeOutput_())
170 polymerDeadPoreVolume_[globalDofIdx] =
171 scalarValue(intQuants.polymerDeadPoreVolume());
172
173 if (polymerRockDensityOutput_())
174 polymerRockDensity_[globalDofIdx] =
175 scalarValue(intQuants.polymerRockDensity());
176
177 if (polymerAdsorptionOutput_())
178 polymerAdsorption_[globalDofIdx] =
179 scalarValue(intQuants.polymerAdsorption());
180
181 if (polymerViscosityCorrectionOutput_())
182 polymerViscosityCorrection_[globalDofIdx] =
183 scalarValue(intQuants.polymerViscosityCorrection());
184
185 if (waterViscosityCorrectionOutput_())
186 waterViscosityCorrection_[globalDofIdx] =
187 scalarValue(intQuants.waterViscosityCorrection());
188 }
189 }
190
195 {
196 VtkMultiWriter *vtkWriter = dynamic_cast<VtkMultiWriter*>(&baseWriter);
197 if (!vtkWriter)
198 return;
199
200 if (!enablePolymer)
201 return;
202
203 if (polymerConcentrationOutput_())
204 this->commitScalarBuffer_(baseWriter, "polymer concentration", polymerConcentration_);
205
206 if (polymerDeadPoreVolumeOutput_())
207 this->commitScalarBuffer_(baseWriter, "dead pore volume fraction", polymerDeadPoreVolume_);
208
209 if (polymerRockDensityOutput_())
210 this->commitScalarBuffer_(baseWriter, "polymer rock density", polymerRockDensity_);
211
212 if (polymerAdsorptionOutput_())
213 this->commitScalarBuffer_(baseWriter, "polymer adsorption", polymerAdsorption_);
214
215 if (polymerViscosityCorrectionOutput_())
216 this->commitScalarBuffer_(baseWriter, "polymer viscosity correction", polymerViscosityCorrection_);
217
218 if (waterViscosityCorrectionOutput_())
219 this->commitScalarBuffer_(baseWriter, "water viscosity correction", waterViscosityCorrection_);
220 }
221
222private:
223 static bool polymerConcentrationOutput_()
224 {
225 static bool val = Parameters::Get<Parameters::VtkWritePolymerConcentration>();
226 return val;
227 }
228
229 static bool polymerDeadPoreVolumeOutput_()
230 {
231 static bool val = Parameters::Get<Parameters::VtkWritePolymerDeadPoreVolume>();
232 return val;
233 }
234
235 static bool polymerRockDensityOutput_()
236 {
237 static bool val = Parameters::Get<Parameters::VtkWritePolymerRockDensity>();
238 return val;
239 }
240
241 static bool polymerAdsorptionOutput_()
242 {
243 static bool val = Parameters::Get<Parameters::VtkWritePolymerAdsorption>();
244 return val;
245 }
246
247 static bool polymerViscosityCorrectionOutput_()
248 {
249 static bool val = Parameters::Get<Parameters::VtkWritePolymerViscosityCorrection>();
250 return val;
251 }
252
253 static bool waterViscosityCorrectionOutput_()
254 {
255 static bool val = Parameters::Get<Parameters::VtkWritePolymerViscosityCorrection>();
256 return val;
257 }
258
259 ScalarBuffer polymerConcentration_;
260 ScalarBuffer polymerDeadPoreVolume_;
261 ScalarBuffer polymerRockDensity_;
262 ScalarBuffer polymerAdsorption_;
263 ScalarBuffer polymerViscosityCorrection_;
264 ScalarBuffer waterViscosityCorrection_;
265};
266} // namespace Opm
267
268#endif
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 polymer related quantities.
Definition: vtkblackoilpolymermodule.hh:71
void allocBuffers()
Allocate memory for the scalar fields we would like to write to the VTK file.
Definition: vtkblackoilpolymermodule.hh:125
static void registerParameters()
Register all run-time parameters for the multi-phase VTK output module.
Definition: vtkblackoilpolymermodule.hh:96
void processElement(const ElementContext &elemCtx)
Modify the internal buffers according to the intensive quantities relevant for an element.
Definition: vtkblackoilpolymermodule.hh:152
VtkBlackOilPolymerModule(const Simulator &simulator)
Definition: vtkblackoilpolymermodule.hh:88
void commitBuffers(BaseOutputWriter &baseWriter)
Add all buffers to the VTK output writer.
Definition: vtkblackoilpolymermodule.hh:194
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
The generic type tag for problems using the immiscible multi-phase model.
Definition: blackoilmodel.hh:74
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: vtkblackoilpolymermodule.hh:59
static constexpr bool value
Definition: vtkblackoilpolymermodule.hh:59
Definition: vtkblackoilpolymermodule.hh:54
static constexpr bool value
Definition: vtkblackoilpolymermodule.hh:54
Definition: vtkblackoilpolymermodule.hh:55
static constexpr bool value
Definition: vtkblackoilpolymermodule.hh:55
Definition: vtkblackoilpolymermodule.hh:58
static constexpr bool value
Definition: vtkblackoilpolymermodule.hh:58
Definition: vtkblackoilpolymermodule.hh:56
static constexpr bool value
Definition: vtkblackoilpolymermodule.hh:56
Definition: vtkblackoilpolymermodule.hh:57
static constexpr bool value
Definition: vtkblackoilpolymermodule.hh:57
Definition: vtkblackoilpolymermodule.hh:47