vtkcompositionmodule.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_COMPOSITION_MODULE_HH
28#define EWOMS_VTK_COMPOSITION_MODULE_HH
29
30#include "vtkmultiwriter.hh"
31#include "baseoutputmodule.hh"
32
35
36#include <opm/material/common/MathToolbox.hpp>
37
38namespace Opm::Properties {
39
40namespace TTag {
41
42// create new type tag for the VTK composition output
44
45} // namespace TTag
46
47// create the property tags needed for the composition module
48template<class TypeTag, class MyTypeTag>
50template<class TypeTag, class MyTypeTag>
52template<class TypeTag, class MyTypeTag>
54template<class TypeTag, class MyTypeTag>
56template<class TypeTag, class MyTypeTag>
58template<class TypeTag, class MyTypeTag>
60template<class TypeTag, class MyTypeTag>
62
63// set default values for what quantities to output
64template<class TypeTag>
65struct VtkWriteMassFractions<TypeTag, TTag::VtkComposition> { static constexpr bool value = false; };
66template<class TypeTag>
67struct VtkWriteMoleFractions<TypeTag, TTag::VtkComposition> { static constexpr bool value = true; };
68template<class TypeTag>
69struct VtkWriteTotalMassFractions<TypeTag, TTag::VtkComposition> { static constexpr bool value = false; };
70template<class TypeTag>
71struct VtkWriteTotalMoleFractions<TypeTag, TTag::VtkComposition> { static constexpr bool value = false; };
72template<class TypeTag>
73struct VtkWriteMolarities<TypeTag, TTag::VtkComposition> { static constexpr bool value = false; };
74template<class TypeTag>
75struct VtkWriteFugacities<TypeTag, TTag::VtkComposition> { static constexpr bool value = false; };
76template<class TypeTag>
77struct VtkWriteFugacityCoeffs<TypeTag, TTag::VtkComposition> { static constexpr bool value = false; };
78
79} // namespace Opm::Properties
80
81namespace Opm {
82
95template <class TypeTag>
97{
99
104
106
107 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
108 enum { numComponents = getPropValue<TypeTag, Properties::NumComponents>() };
109
110 static const int vtkFormat = getPropValue<TypeTag, Properties::VtkOutputFormat>();
112
115
116public:
117 VtkCompositionModule(const Simulator& simulator)
118 : ParentType(simulator)
119 { }
120
124 static void registerParameters()
125 {
126 Parameters::registerParam<TypeTag, Properties::VtkWriteMassFractions>
127 ("Include mass fractions in the VTK output files");
128 Parameters::registerParam<TypeTag, Properties::VtkWriteMoleFractions>
129 ("Include mole fractions in the VTK output files");
130 Parameters::registerParam<TypeTag, Properties::VtkWriteTotalMassFractions>
131 ("Include total mass fractions in the VTK output files");
132 Parameters::registerParam<TypeTag, Properties::VtkWriteTotalMoleFractions>
133 ("Include total mole fractions in the VTK output files");
134 Parameters::registerParam<TypeTag, Properties::VtkWriteMolarities>
135 ("Include component molarities in the VTK output files");
136 Parameters::registerParam<TypeTag, Properties::VtkWriteFugacities>
137 ("Include component fugacities in the VTK output files");
138 Parameters::registerParam<TypeTag, Properties::VtkWriteFugacityCoeffs>
139 ("Include component fugacity coefficients in the VTK output files");
140 }
141
147 {
148 if (moleFracOutput_())
149 this->resizePhaseComponentBuffer_(moleFrac_);
150 if (massFracOutput_())
151 this->resizePhaseComponentBuffer_(massFrac_);
152 if (totalMassFracOutput_())
153 this->resizeComponentBuffer_(totalMassFrac_);
154 if (totalMoleFracOutput_())
155 this->resizeComponentBuffer_(totalMoleFrac_);
156 if (molarityOutput_())
157 this->resizePhaseComponentBuffer_(molarity_);
158
159 if (fugacityOutput_())
160 this->resizeComponentBuffer_(fugacity_);
161 if (fugacityCoeffOutput_())
162 this->resizePhaseComponentBuffer_(fugacityCoeff_);
163 }
164
169 void processElement(const ElementContext& elemCtx)
170 {
171 using Toolbox = MathToolbox<Evaluation>;
172
173 if (!Parameters::get<TypeTag, Properties::EnableVtkOutput>())
174 return;
175
176 for (unsigned i = 0; i < elemCtx.numPrimaryDof(/*timeIdx=*/0); ++i) {
177 unsigned I = elemCtx.globalSpaceIndex(i, /*timeIdx=*/0);
178 const auto& intQuants = elemCtx.intensiveQuantities(i, /*timeIdx=*/0);
179 const auto& fs = intQuants.fluidState();
180
181 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
182 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
183 if (moleFracOutput_())
184 moleFrac_[phaseIdx][compIdx][I] = Toolbox::value(fs.moleFraction(phaseIdx, compIdx));
185 if (massFracOutput_())
186 massFrac_[phaseIdx][compIdx][I] = Toolbox::value(fs.massFraction(phaseIdx, compIdx));
187 if (molarityOutput_())
188 molarity_[phaseIdx][compIdx][I] = Toolbox::value(fs.molarity(phaseIdx, compIdx));
189
190 if (fugacityCoeffOutput_())
191 fugacityCoeff_[phaseIdx][compIdx][I] =
192 Toolbox::value(fs.fugacityCoefficient(phaseIdx, compIdx));
193 }
194 }
195
196 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
197 if (totalMassFracOutput_()) {
198 Scalar compMass = 0;
199 Scalar totalMass = 0;
200 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
201 totalMass += Toolbox::value(fs.density(phaseIdx)) * Toolbox::value(fs.saturation(phaseIdx));
202 compMass +=
203 Toolbox::value(fs.density(phaseIdx))
204 *Toolbox::value(fs.saturation(phaseIdx))
205 *Toolbox::value(fs.massFraction(phaseIdx, compIdx));
206 }
207 totalMassFrac_[compIdx][I] = compMass / totalMass;
208 }
209 if (totalMoleFracOutput_()) {
210 Scalar compMoles = 0;
211 Scalar totalMoles = 0;
212 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
213 totalMoles +=
214 Toolbox::value(fs.molarDensity(phaseIdx))
215 *Toolbox::value(fs.saturation(phaseIdx));
216 compMoles +=
217 Toolbox::value(fs.molarDensity(phaseIdx))
218 *Toolbox::value(fs.saturation(phaseIdx))
219 *Toolbox::value(fs.moleFraction(phaseIdx, compIdx));
220 }
221 totalMoleFrac_[compIdx][I] = compMoles / totalMoles;
222 }
223 if (fugacityOutput_())
224 fugacity_[compIdx][I] = Toolbox::value(intQuants.fluidState().fugacity(/*phaseIdx=*/0, compIdx));
225 }
226 }
227 }
228
233 {
234 VtkMultiWriter *vtkWriter = dynamic_cast<VtkMultiWriter*>(&baseWriter);
235 if (!vtkWriter) {
236 return;
237 }
238
239 if (moleFracOutput_())
240 this->commitPhaseComponentBuffer_(baseWriter, "moleFrac_%s^%s", moleFrac_);
241 if (massFracOutput_())
242 this->commitPhaseComponentBuffer_(baseWriter, "massFrac_%s^%s", massFrac_);
243 if (molarityOutput_())
244 this->commitPhaseComponentBuffer_(baseWriter, "molarity_%s^%s", molarity_);
245 if (totalMassFracOutput_())
246 this->commitComponentBuffer_(baseWriter, "totalMassFrac^%s", totalMassFrac_);
247 if (totalMoleFracOutput_())
248 this->commitComponentBuffer_(baseWriter, "totalMoleFrac^%s", totalMoleFrac_);
249
250 if (fugacityOutput_())
251 this->commitComponentBuffer_(baseWriter, "fugacity^%s", fugacity_);
252 if (fugacityCoeffOutput_())
253 this->commitPhaseComponentBuffer_(baseWriter, "fugacityCoeff_%s^%s", fugacityCoeff_);
254 }
255
256private:
257 static bool massFracOutput_()
258 {
259 static bool val = Parameters::get<TypeTag, Properties::VtkWriteMassFractions>();
260 return val;
261 }
262
263 static bool moleFracOutput_()
264 {
265 static bool val = Parameters::get<TypeTag, Properties::VtkWriteMoleFractions>();
266 return val;
267 }
268
269 static bool totalMassFracOutput_()
270 {
271 static bool val = Parameters::get<TypeTag, Properties::VtkWriteTotalMassFractions>();
272 return val;
273 }
274
275 static bool totalMoleFracOutput_()
276 {
277 static bool val = Parameters::get<TypeTag, Properties::VtkWriteTotalMoleFractions>();
278 return val;
279 }
280
281 static bool molarityOutput_()
282 {
283 static bool val = Parameters::get<TypeTag, Properties::VtkWriteMolarities>();
284 return val;
285 }
286
287 static bool fugacityOutput_()
288 {
289 static bool val = Parameters::get<TypeTag, Properties::VtkWriteFugacities>();
290 return val;
291 }
292
293 static bool fugacityCoeffOutput_()
294 {
295 static bool val = Parameters::get<TypeTag, Properties::VtkWriteFugacityCoeffs>();
296 return val;
297 }
298
299 PhaseComponentBuffer moleFrac_;
300 PhaseComponentBuffer massFrac_;
301 PhaseComponentBuffer molarity_;
302 ComponentBuffer totalMassFrac_;
303 ComponentBuffer totalMoleFrac_;
304
305 ComponentBuffer fugacity_;
306 PhaseComponentBuffer fugacityCoeff_;
307};
308
309} // namespace Opm
310
311#endif
The base class for writer modules.
Definition: baseoutputmodule.hh:73
void commitComponentBuffer_(BaseOutputWriter &baseWriter, const char *pattern, ComponentBuffer &buffer, BufferType bufferType=DofBuffer)
Add a component-specific buffer to the result file.
Definition: baseoutputmodule.hh:442
std::array< ScalarBuffer, numComponents > ComponentBuffer
Definition: baseoutputmodule.hh:97
void resizeComponentBuffer_(ComponentBuffer &buffer, BufferType bufferType=DofBuffer)
Allocate the space for a buffer storing a component specific quantity.
Definition: baseoutputmodule.hh:269
void resizePhaseComponentBuffer_(PhaseComponentBuffer &buffer, BufferType bufferType=DofBuffer)
Allocate the space for a buffer storing a phase and component specific buffer.
Definition: baseoutputmodule.hh:292
void commitPhaseComponentBuffer_(BaseOutputWriter &baseWriter, const char *pattern, PhaseComponentBuffer &buffer, BufferType bufferType=DofBuffer)
Add a phase and component specific quantities to the output.
Definition: baseoutputmodule.hh:465
std::array< std::array< ScalarBuffer, numComponents >, numPhases > PhaseComponentBuffer
Definition: baseoutputmodule.hh:98
The base class for all output writers.
Definition: baseoutputwriter.hh:44
VTK output module for the fluid composition.
Definition: vtkcompositionmodule.hh:97
void processElement(const ElementContext &elemCtx)
Modify the internal buffers according to the intensive quantities relevant for an element.
Definition: vtkcompositionmodule.hh:169
void commitBuffers(BaseOutputWriter &baseWriter)
Add all buffers to the VTK output writer.
Definition: vtkcompositionmodule.hh:232
static void registerParameters()
Register all run-time parameters for the Vtk output module.
Definition: vtkcompositionmodule.hh:124
void allocBuffers()
Allocate memory for the scalar fields we would like to write to the VTK file.
Definition: vtkcompositionmodule.hh:146
VtkCompositionModule(const Simulator &simulator)
Definition: vtkcompositionmodule.hh:117
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: vtkcompositionmodule.hh:43
a tag to mark properties as undefined
Definition: propertysystem.hh:40
Definition: vtkcompositionmodule.hh:59
Definition: vtkcompositionmodule.hh:61
Definition: vtkcompositionmodule.hh:49
Definition: vtkcompositionmodule.hh:57
Definition: vtkcompositionmodule.hh:51
Definition: vtkcompositionmodule.hh:53
Definition: vtkcompositionmodule.hh:55