27#ifndef EWOMS_VTK_COMPOSITION_MODULE_HH
28#define EWOMS_VTK_COMPOSITION_MODULE_HH
30#include <opm/material/common/MathToolbox.hpp>
67template <
class TypeTag>
79 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
80 enum { numComponents = getPropValue<TypeTag, Properties::NumComponents>() };
82 static const int vtkFormat = getPropValue<TypeTag, Properties::VtkOutputFormat>();
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");
120 if (moleFracOutput_())
122 if (massFracOutput_())
124 if (totalMassFracOutput_())
126 if (totalMoleFracOutput_())
128 if (molarityOutput_())
131 if (fugacityOutput_())
133 if (fugacityCoeffOutput_())
143 using Toolbox = MathToolbox<Evaluation>;
145 if (!Parameters::Get<Parameters::EnableVtkOutput>()) {
149 for (
unsigned i = 0; i < elemCtx.numPrimaryDof(0); ++i) {
150 unsigned I = elemCtx.globalSpaceIndex(i, 0);
151 const auto& intQuants = elemCtx.intensiveQuantities(i, 0);
152 const auto& fs = intQuants.fluidState();
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));
163 if (fugacityCoeffOutput_())
164 fugacityCoeff_[phaseIdx][compIdx][I] =
165 Toolbox::value(fs.fugacityCoefficient(phaseIdx, compIdx));
169 for (
unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
170 if (totalMassFracOutput_()) {
172 Scalar totalMass = 0;
173 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
174 totalMass += Toolbox::value(fs.density(phaseIdx)) * Toolbox::value(fs.saturation(phaseIdx));
176 Toolbox::value(fs.density(phaseIdx))
177 *Toolbox::value(fs.saturation(phaseIdx))
178 *Toolbox::value(fs.massFraction(phaseIdx, compIdx));
180 totalMassFrac_[compIdx][I] = compMass / totalMass;
182 if (totalMoleFracOutput_()) {
183 Scalar compMoles = 0;
184 Scalar totalMoles = 0;
185 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
187 Toolbox::value(fs.molarDensity(phaseIdx))
188 *Toolbox::value(fs.saturation(phaseIdx));
190 Toolbox::value(fs.molarDensity(phaseIdx))
191 *Toolbox::value(fs.saturation(phaseIdx))
192 *Toolbox::value(fs.moleFraction(phaseIdx, compIdx));
194 totalMoleFrac_[compIdx][I] = compMoles / totalMoles;
196 if (fugacityOutput_())
197 fugacity_[compIdx][I] = Toolbox::value(intQuants.fluidState().fugacity(0, compIdx));
212 if (moleFracOutput_())
214 if (massFracOutput_())
216 if (molarityOutput_())
218 if (totalMassFracOutput_())
220 if (totalMoleFracOutput_())
223 if (fugacityOutput_())
225 if (fugacityCoeffOutput_())
230 static bool massFracOutput_()
232 static bool val = Parameters::Get<Parameters::VtkWriteMassFractions>();
236 static bool moleFracOutput_()
238 static bool val = Parameters::Get<Parameters::VtkWriteMoleFractions>();
242 static bool totalMassFracOutput_()
244 static bool val = Parameters::Get<Parameters::VtkWriteTotalMassFractions>();
248 static bool totalMoleFracOutput_()
250 static bool val = Parameters::Get<Parameters::VtkWriteTotalMoleFractions>();
254 static bool molarityOutput_()
256 static bool val = Parameters::Get<Parameters::VtkWriteMolarities>();
260 static bool fugacityOutput_()
262 static bool val = Parameters::Get<Parameters::VtkWriteFugacities>();
266 static bool fugacityCoeffOutput_()
268 static bool val = Parameters::Get<Parameters::VtkWriteFugacityCoeffs>();
272 PhaseComponentBuffer moleFrac_;
273 PhaseComponentBuffer massFrac_;
274 PhaseComponentBuffer molarity_;
275 ComponentBuffer totalMassFrac_;
276 ComponentBuffer totalMoleFrac_;
278 ComponentBuffer fugacity_;
279 PhaseComponentBuffer fugacityCoeff_;
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