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 <opm/material/common/MathToolbox.hpp>
31
33
36
39
40namespace Opm::Parameters {
41
42// set default values for what quantities to output
43struct VtkWriteMassFractions { static constexpr bool value = false; };
44struct VtkWriteMoleFractions { static constexpr bool value = true; };
45struct VtkWriteTotalMassFractions { static constexpr bool value = false; };
46struct VtkWriteTotalMoleFractions { static constexpr bool value = false; };
47struct VtkWriteMolarities { static constexpr bool value = false; };
48struct VtkWriteFugacities { static constexpr bool value = false; };
49struct VtkWriteFugacityCoeffs { static constexpr bool value = false; };
50
51} // namespace Opm::Properties
52
53namespace Opm {
54
67template <class TypeTag>
69{
71
76
78
79 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
80 enum { numComponents = getPropValue<TypeTag, Properties::NumComponents>() };
81
82 static const int vtkFormat = getPropValue<TypeTag, Properties::VtkOutputFormat>();
84
87
88public:
89 VtkCompositionModule(const Simulator& simulator)
90 : ParentType(simulator)
91 { }
92
96 static void registerParameters()
97 {
98 Parameters::Register<Parameters::VtkWriteMassFractions>
99 ("Include mass fractions in the VTK output files");
100 Parameters::Register<Parameters::VtkWriteMoleFractions>
101 ("Include mole fractions in the VTK output files");
102 Parameters::Register<Parameters::VtkWriteTotalMassFractions>
103 ("Include total mass fractions in the VTK output files");
104 Parameters::Register<Parameters::VtkWriteTotalMoleFractions>
105 ("Include total mole fractions in the VTK output files");
106 Parameters::Register<Parameters::VtkWriteMolarities>
107 ("Include component molarities in the VTK output files");
108 Parameters::Register<Parameters::VtkWriteFugacities>
109 ("Include component fugacities in the VTK output files");
110 Parameters::Register<Parameters::VtkWriteFugacityCoeffs>
111 ("Include component fugacity coefficients in the VTK output files");
112 }
113
119 {
120 if (moleFracOutput_())
121 this->resizePhaseComponentBuffer_(moleFrac_);
122 if (massFracOutput_())
123 this->resizePhaseComponentBuffer_(massFrac_);
124 if (totalMassFracOutput_())
125 this->resizeComponentBuffer_(totalMassFrac_);
126 if (totalMoleFracOutput_())
127 this->resizeComponentBuffer_(totalMoleFrac_);
128 if (molarityOutput_())
129 this->resizePhaseComponentBuffer_(molarity_);
130
131 if (fugacityOutput_())
132 this->resizeComponentBuffer_(fugacity_);
133 if (fugacityCoeffOutput_())
134 this->resizePhaseComponentBuffer_(fugacityCoeff_);
135 }
136
141 void processElement(const ElementContext& elemCtx)
142 {
143 using Toolbox = MathToolbox<Evaluation>;
144
145 if (!Parameters::Get<Parameters::EnableVtkOutput>()) {
146 return;
147 }
148
149 for (unsigned i = 0; i < elemCtx.numPrimaryDof(/*timeIdx=*/0); ++i) {
150 unsigned I = elemCtx.globalSpaceIndex(i, /*timeIdx=*/0);
151 const auto& intQuants = elemCtx.intensiveQuantities(i, /*timeIdx=*/0);
152 const auto& fs = intQuants.fluidState();
153
154 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
155 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
156 if (moleFracOutput_())
157 moleFrac_[phaseIdx][compIdx][I] = Toolbox::value(fs.moleFraction(phaseIdx, compIdx));
158 if (massFracOutput_())
159 massFrac_[phaseIdx][compIdx][I] = Toolbox::value(fs.massFraction(phaseIdx, compIdx));
160 if (molarityOutput_())
161 molarity_[phaseIdx][compIdx][I] = Toolbox::value(fs.molarity(phaseIdx, compIdx));
162
163 if (fugacityCoeffOutput_())
164 fugacityCoeff_[phaseIdx][compIdx][I] =
165 Toolbox::value(fs.fugacityCoefficient(phaseIdx, compIdx));
166 }
167 }
168
169 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
170 if (totalMassFracOutput_()) {
171 Scalar compMass = 0;
172 Scalar totalMass = 0;
173 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
174 totalMass += Toolbox::value(fs.density(phaseIdx)) * Toolbox::value(fs.saturation(phaseIdx));
175 compMass +=
176 Toolbox::value(fs.density(phaseIdx))
177 *Toolbox::value(fs.saturation(phaseIdx))
178 *Toolbox::value(fs.massFraction(phaseIdx, compIdx));
179 }
180 totalMassFrac_[compIdx][I] = compMass / totalMass;
181 }
182 if (totalMoleFracOutput_()) {
183 Scalar compMoles = 0;
184 Scalar totalMoles = 0;
185 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
186 totalMoles +=
187 Toolbox::value(fs.molarDensity(phaseIdx))
188 *Toolbox::value(fs.saturation(phaseIdx));
189 compMoles +=
190 Toolbox::value(fs.molarDensity(phaseIdx))
191 *Toolbox::value(fs.saturation(phaseIdx))
192 *Toolbox::value(fs.moleFraction(phaseIdx, compIdx));
193 }
194 totalMoleFrac_[compIdx][I] = compMoles / totalMoles;
195 }
196 if (fugacityOutput_())
197 fugacity_[compIdx][I] = Toolbox::value(intQuants.fluidState().fugacity(/*phaseIdx=*/0, compIdx));
198 }
199 }
200 }
201
206 {
207 VtkMultiWriter *vtkWriter = dynamic_cast<VtkMultiWriter*>(&baseWriter);
208 if (!vtkWriter) {
209 return;
210 }
211
212 if (moleFracOutput_())
213 this->commitPhaseComponentBuffer_(baseWriter, "moleFrac_%s^%s", moleFrac_);
214 if (massFracOutput_())
215 this->commitPhaseComponentBuffer_(baseWriter, "massFrac_%s^%s", massFrac_);
216 if (molarityOutput_())
217 this->commitPhaseComponentBuffer_(baseWriter, "molarity_%s^%s", molarity_);
218 if (totalMassFracOutput_())
219 this->commitComponentBuffer_(baseWriter, "totalMassFrac^%s", totalMassFrac_);
220 if (totalMoleFracOutput_())
221 this->commitComponentBuffer_(baseWriter, "totalMoleFrac^%s", totalMoleFrac_);
222
223 if (fugacityOutput_())
224 this->commitComponentBuffer_(baseWriter, "fugacity^%s", fugacity_);
225 if (fugacityCoeffOutput_())
226 this->commitPhaseComponentBuffer_(baseWriter, "fugacityCoeff_%s^%s", fugacityCoeff_);
227 }
228
229private:
230 static bool massFracOutput_()
231 {
232 static bool val = Parameters::Get<Parameters::VtkWriteMassFractions>();
233 return val;
234 }
235
236 static bool moleFracOutput_()
237 {
238 static bool val = Parameters::Get<Parameters::VtkWriteMoleFractions>();
239 return val;
240 }
241
242 static bool totalMassFracOutput_()
243 {
244 static bool val = Parameters::Get<Parameters::VtkWriteTotalMassFractions>();
245 return val;
246 }
247
248 static bool totalMoleFracOutput_()
249 {
250 static bool val = Parameters::Get<Parameters::VtkWriteTotalMoleFractions>();
251 return val;
252 }
253
254 static bool molarityOutput_()
255 {
256 static bool val = Parameters::Get<Parameters::VtkWriteMolarities>();
257 return val;
258 }
259
260 static bool fugacityOutput_()
261 {
262 static bool val = Parameters::Get<Parameters::VtkWriteFugacities>();
263 return val;
264 }
265
266 static bool fugacityCoeffOutput_()
267 {
268 static bool val = Parameters::Get<Parameters::VtkWriteFugacityCoeffs>();
269 return val;
270 }
271
272 PhaseComponentBuffer moleFrac_;
273 PhaseComponentBuffer massFrac_;
274 PhaseComponentBuffer molarity_;
275 ComponentBuffer totalMassFrac_;
276 ComponentBuffer totalMoleFrac_;
277
278 ComponentBuffer fugacity_;
279 PhaseComponentBuffer fugacityCoeff_;
280};
281
282} // namespace Opm
283
284#endif
The base class for writer modules.
Definition: baseoutputmodule.hh:67
void commitComponentBuffer_(BaseOutputWriter &baseWriter, const char *pattern, ComponentBuffer &buffer, BufferType bufferType=DofBuffer)
Add a component-specific buffer to the result file.
Definition: baseoutputmodule.hh:361
std::array< ScalarBuffer, numComponents > ComponentBuffer
Definition: baseoutputmodule.hh:91
void resizeComponentBuffer_(ComponentBuffer &buffer, BufferType bufferType=DofBuffer)
Allocate the space for a buffer storing a component specific quantity.
Definition: baseoutputmodule.hh:215
void resizePhaseComponentBuffer_(PhaseComponentBuffer &buffer, BufferType bufferType=DofBuffer)
Allocate the space for a buffer storing a phase and component specific buffer.
Definition: baseoutputmodule.hh:229
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:377
std::array< std::array< ScalarBuffer, numComponents >, numPhases > PhaseComponentBuffer
Definition: baseoutputmodule.hh:92
The base class for all output writers.
Definition: baseoutputwriter.hh:44
VTK output module for the fluid composition.
Definition: vtkcompositionmodule.hh:69
void processElement(const ElementContext &elemCtx)
Modify the internal buffers according to the intensive quantities relevant for an element.
Definition: vtkcompositionmodule.hh:141
void commitBuffers(BaseOutputWriter &baseWriter)
Add all buffers to the VTK output writer.
Definition: vtkcompositionmodule.hh:205
static void registerParameters()
Register all run-time parameters for the Vtk output module.
Definition: vtkcompositionmodule.hh:96
void allocBuffers()
Allocate memory for the scalar fields we would like to write to the VTK file.
Definition: vtkcompositionmodule.hh:118
VtkCompositionModule(const Simulator &simulator)
Definition: vtkcompositionmodule.hh:89
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: vtkcompositionmodule.hh:48
static constexpr bool value
Definition: vtkcompositionmodule.hh:48
Definition: vtkcompositionmodule.hh:49
static constexpr bool value
Definition: vtkcompositionmodule.hh:49
Definition: vtkcompositionmodule.hh:43
static constexpr bool value
Definition: vtkcompositionmodule.hh:43
Definition: vtkcompositionmodule.hh:47
static constexpr bool value
Definition: vtkcompositionmodule.hh:47
Definition: vtkcompositionmodule.hh:44
static constexpr bool value
Definition: vtkcompositionmodule.hh:44
Definition: vtkcompositionmodule.hh:45
static constexpr bool value
Definition: vtkcompositionmodule.hh:45
Definition: vtkcompositionmodule.hh:46
static constexpr bool value
Definition: vtkcompositionmodule.hh:46