vtkcompositionmodule.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_COMPOSITION_MODULE_HPP
28#define OPM_VTK_COMPOSITION_MODULE_HPP
29
30#include <opm/material/common/MathToolbox.hpp>
31
33
37
40
41namespace Opm {
42
55template <class TypeTag>
57{
59
64
66
67 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
68 enum { numComponents = getPropValue<TypeTag, Properties::NumComponents>() };
69
70 static const int vtkFormat = getPropValue<TypeTag, Properties::VtkOutputFormat>();
72
75
76public:
77 VtkCompositionModule(const Simulator& simulator)
78 : ParentType(simulator)
79 {
80 params_.read();
81 }
82
86 static void registerParameters()
87 {
89 }
90
96 {
97 if (params_.moleFracOutput_) {
98 this->resizePhaseComponentBuffer_(moleFrac_);
99 }
100 if (params_.massFracOutput_) {
101 this->resizePhaseComponentBuffer_(massFrac_);
102 }
103 if (params_.totalMassFracOutput_) {
104 this->resizeComponentBuffer_(totalMassFrac_);
105 }
106 if (params_.totalMoleFracOutput_) {
107 this->resizeComponentBuffer_(totalMoleFrac_);
108 }
109 if (params_.molarityOutput_) {
110 this->resizePhaseComponentBuffer_(molarity_);
111 }
112
113 if (params_.fugacityOutput_) {
114 this->resizeComponentBuffer_(fugacity_);
115 }
116 if (params_.fugacityCoeffOutput_) {
117 this->resizePhaseComponentBuffer_(fugacityCoeff_);
118 }
119 }
120
125 void processElement(const ElementContext& elemCtx)
126 {
127 using Toolbox = MathToolbox<Evaluation>;
128
129 if (!Parameters::Get<Parameters::EnableVtkOutput>()) {
130 return;
131 }
132
133 for (unsigned i = 0; i < elemCtx.numPrimaryDof(/*timeIdx=*/0); ++i) {
134 unsigned I = elemCtx.globalSpaceIndex(i, /*timeIdx=*/0);
135 const auto& intQuants = elemCtx.intensiveQuantities(i, /*timeIdx=*/0);
136 const auto& fs = intQuants.fluidState();
137
138 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
139 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
140 if (params_.moleFracOutput_) {
141 moleFrac_[phaseIdx][compIdx][I] = Toolbox::value(fs.moleFraction(phaseIdx, compIdx));
142 }
143 if (params_.massFracOutput_) {
144 massFrac_[phaseIdx][compIdx][I] = Toolbox::value(fs.massFraction(phaseIdx, compIdx));
145 }
146 if (params_.molarityOutput_) {
147 molarity_[phaseIdx][compIdx][I] = Toolbox::value(fs.molarity(phaseIdx, compIdx));
148 }
149
150 if (params_.fugacityCoeffOutput_) {
151 fugacityCoeff_[phaseIdx][compIdx][I] =
152 Toolbox::value(fs.fugacityCoefficient(phaseIdx, compIdx));
153 }
154 }
155 }
156
157 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
158 if (params_.totalMassFracOutput_) {
159 Scalar compMass = 0;
160 Scalar totalMass = 0;
161 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
162 totalMass += Toolbox::value(fs.density(phaseIdx)) * Toolbox::value(fs.saturation(phaseIdx));
163 compMass +=
164 Toolbox::value(fs.density(phaseIdx))
165 *Toolbox::value(fs.saturation(phaseIdx))
166 *Toolbox::value(fs.massFraction(phaseIdx, compIdx));
167 }
168 totalMassFrac_[compIdx][I] = compMass / totalMass;
169 }
170 if (params_.totalMoleFracOutput_) {
171 Scalar compMoles = 0;
172 Scalar totalMoles = 0;
173 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
174 totalMoles +=
175 Toolbox::value(fs.molarDensity(phaseIdx))
176 *Toolbox::value(fs.saturation(phaseIdx));
177 compMoles +=
178 Toolbox::value(fs.molarDensity(phaseIdx))
179 *Toolbox::value(fs.saturation(phaseIdx))
180 *Toolbox::value(fs.moleFraction(phaseIdx, compIdx));
181 }
182 totalMoleFrac_[compIdx][I] = compMoles / totalMoles;
183 }
184 if (params_.fugacityOutput_) {
185 fugacity_[compIdx][I] = Toolbox::value(intQuants.fluidState().fugacity(/*phaseIdx=*/0, compIdx));
186 }
187 }
188 }
189 }
190
195 {
196 VtkMultiWriter* vtkWriter = dynamic_cast<VtkMultiWriter*>(&baseWriter);
197 if (!vtkWriter) {
198 return;
199 }
200
201 if (params_.moleFracOutput_) {
202 this->commitPhaseComponentBuffer_(baseWriter, "moleFrac_%s^%s", moleFrac_);
203 }
204 if (params_.massFracOutput_) {
205 this->commitPhaseComponentBuffer_(baseWriter, "massFrac_%s^%s", massFrac_);
206 }
207 if (params_.molarityOutput_) {
208 this->commitPhaseComponentBuffer_(baseWriter, "molarity_%s^%s", molarity_);
209 }
210 if (params_.totalMassFracOutput_) {
211 this->commitComponentBuffer_(baseWriter, "totalMassFrac^%s", totalMassFrac_);
212 }
213 if (params_.totalMoleFracOutput_) {
214 this->commitComponentBuffer_(baseWriter, "totalMoleFrac^%s", totalMoleFrac_);
215 }
216
217 if (params_.fugacityOutput_) {
218 this->commitComponentBuffer_(baseWriter, "fugacity^%s", fugacity_);
219 }
220 if (params_.fugacityCoeffOutput_) {
221 this->commitPhaseComponentBuffer_(baseWriter, "fugacityCoeff_%s^%s", fugacityCoeff_);
222 }
223 }
224
225private:
226 VtkCompositionParams params_{};
227 PhaseComponentBuffer moleFrac_{};
228 PhaseComponentBuffer massFrac_{};
229 PhaseComponentBuffer molarity_{};
230 ComponentBuffer totalMassFrac_{};
231 ComponentBuffer totalMoleFrac_{};
232
233 ComponentBuffer fugacity_{};
234 PhaseComponentBuffer fugacityCoeff_{};
235};
236
237} // namespace Opm
238
239#endif // OPM_VTK_COMPOSITION_MODULE_HPP
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.hpp:57
void processElement(const ElementContext &elemCtx)
Modify the internal buffers according to the intensive quantities relevant for an element.
Definition: vtkcompositionmodule.hpp:125
void commitBuffers(BaseOutputWriter &baseWriter)
Add all buffers to the VTK output writer.
Definition: vtkcompositionmodule.hpp:194
static void registerParameters()
Register all run-time parameters for the Vtk output module.
Definition: vtkcompositionmodule.hpp:86
void allocBuffers()
Allocate memory for the scalar fields we would like to write to the VTK file.
Definition: vtkcompositionmodule.hpp:95
VtkCompositionModule(const Simulator &simulator)
Definition: vtkcompositionmodule.hpp:77
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 VtkCompositionModule.
Definition: vtkcompositionparams.hpp:49
bool fugacityCoeffOutput_
Definition: vtkcompositionparams.hpp:62
bool totalMoleFracOutput_
Definition: vtkcompositionparams.hpp:59
bool moleFracOutput_
Definition: vtkcompositionparams.hpp:57
bool molarityOutput_
Definition: vtkcompositionparams.hpp:60
bool massFracOutput_
Definition: vtkcompositionparams.hpp:56
bool totalMassFracOutput_
Definition: vtkcompositionparams.hpp:58
bool fugacityOutput_
Definition: vtkcompositionparams.hpp:61
static void registerParameters()
Registers the parameters in parameter system.
void read()
Reads the parameter values from the parameter system.