27#ifndef OPM_VTK_COMPOSITION_MODULE_HPP
28#define OPM_VTK_COMPOSITION_MODULE_HPP
30#include <opm/material/common/MathToolbox.hpp>
55template <
class TypeTag>
56class VtkCompositionModule :
public BaseOutputModule<TypeTag>
58 using ParentType = BaseOutputModule<TypeTag>;
74 using ComponentBuffer =
typename ParentType::ComponentBuffer;
75 using PhaseComponentBuffer =
typename ParentType::PhaseComponentBuffer;
78 explicit VtkCompositionModule(
const Simulator& simulator)
79 : ParentType(simulator)
98 if (params_.moleFracOutput_) {
101 if (params_.massFracOutput_) {
104 if (params_.totalMassFracOutput_) {
107 if (params_.totalMoleFracOutput_) {
110 if (params_.molarityOutput_) {
114 if (params_.fugacityOutput_) {
117 if (params_.fugacityCoeffOutput_) {
128 using Toolbox = MathToolbox<Evaluation>;
134 for (
unsigned i = 0; i < elemCtx.numPrimaryDof(0); ++i) {
135 const unsigned I = elemCtx.globalSpaceIndex(i, 0);
136 const auto& intQuants = elemCtx.intensiveQuantities(i, 0);
137 const auto& fs = intQuants.fluidState();
139 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
140 for (
unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
141 if (params_.moleFracOutput_) {
142 moleFrac_[phaseIdx][compIdx][I] = Toolbox::value(fs.moleFraction(phaseIdx, compIdx));
144 if (params_.massFracOutput_) {
145 massFrac_[phaseIdx][compIdx][I] = Toolbox::value(fs.massFraction(phaseIdx, compIdx));
147 if (params_.molarityOutput_) {
148 molarity_[phaseIdx][compIdx][I] = Toolbox::value(fs.molarity(phaseIdx, compIdx));
151 if (params_.fugacityCoeffOutput_) {
152 fugacityCoeff_[phaseIdx][compIdx][I] =
153 Toolbox::value(fs.fugacityCoefficient(phaseIdx, compIdx));
158 for (
unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
159 if (params_.totalMassFracOutput_) {
161 Scalar totalMass = 0;
162 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
163 totalMass += Toolbox::value(fs.density(phaseIdx)) * Toolbox::value(fs.saturation(phaseIdx));
165 Toolbox::value(fs.density(phaseIdx)) *
166 Toolbox::value(fs.saturation(phaseIdx)) *
167 Toolbox::value(fs.massFraction(phaseIdx, compIdx));
169 totalMassFrac_[compIdx][I] = compMass / totalMass;
171 if (params_.totalMoleFracOutput_) {
172 Scalar compMoles = 0;
173 Scalar totalMoles = 0;
174 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
176 Toolbox::value(fs.molarDensity(phaseIdx)) *
177 Toolbox::value(fs.saturation(phaseIdx));
179 Toolbox::value(fs.molarDensity(phaseIdx)) *
180 Toolbox::value(fs.saturation(phaseIdx)) *
181 Toolbox::value(fs.moleFraction(phaseIdx, compIdx));
183 totalMoleFrac_[compIdx][I] = compMoles / totalMoles;
185 if (params_.fugacityOutput_) {
186 fugacity_[compIdx][I] = Toolbox::value(intQuants.fluidState().fugacity(0, compIdx));
197 if (!
dynamic_cast<VtkMultiWriter*
>(&baseWriter)) {
201 if (params_.moleFracOutput_) {
203 moleFrac_, BufferType::Dof);
205 if (params_.massFracOutput_) {
207 massFrac_, BufferType::Dof);
209 if (params_.molarityOutput_) {
211 molarity_, BufferType::Dof);
213 if (params_.totalMassFracOutput_) {
215 totalMassFrac_, BufferType::Dof);
217 if (params_.totalMoleFracOutput_) {
219 totalMoleFrac_, BufferType::Dof);
222 if (params_.fugacityOutput_) {
224 fugacity_, BufferType::Dof);
226 if (params_.fugacityCoeffOutput_) {
228 fugacityCoeff_, BufferType::Dof);
234 PhaseComponentBuffer moleFrac_{};
235 PhaseComponentBuffer massFrac_{};
236 PhaseComponentBuffer molarity_{};
237 ComponentBuffer totalMassFrac_{};
238 ComponentBuffer totalMoleFrac_{};
240 ComponentBuffer fugacity_{};
241 PhaseComponentBuffer fugacityCoeff_{};
The base class for writer modules.
void commitPhaseComponentBuffer_(BaseOutputWriter &baseWriter, const char *pattern, PhaseComponentBuffer &buffer, BufferType bufferType)
Add a phase and component specific quantities to the output.
Definition baseoutputmodule.hh:369
void commitComponentBuffer_(BaseOutputWriter &baseWriter, const char *pattern, ComponentBuffer &buffer, BufferType bufferType)
Add a component-specific buffer to the result file.
Definition baseoutputmodule.hh:353
BufferType
Definition baseoutputmodule.hh:143
void resizePhaseComponentBuffer_(PhaseComponentBuffer &buffer, BufferType bufferType)
Allocate the space for a buffer storing a phase and component specific buffer.
Definition baseoutputmodule.hh:224
void resizeComponentBuffer_(ComponentBuffer &buffer, BufferType bufferType)
Allocate the space for a buffer storing a component specific quantity.
Definition baseoutputmodule.hh:211
The base class for all output writers.
Definition baseoutputwriter.hh:46
void processElement(const ElementContext &elemCtx) override
Modify the internal buffers according to the intensive quantities relevant for an element.
Definition vtkcompositionmodule.hpp:126
void commitBuffers(BaseOutputWriter &baseWriter) override
Add all buffers to the VTK output writer.
Definition vtkcompositionmodule.hpp:195
void allocBuffers() override
Allocate memory for the scalar fields we would like to write to the VTK file.
Definition vtkcompositionmodule.hpp:96
static void registerParameters()
Register all run-time parameters for the Vtk output module.
Definition vtkcompositionmodule.hpp:87
Simplifies writing multi-file VTK datasets.
Definition vtkmultiwriter.hh:65
Declare the properties used by the infrastructure code of the finite volume discretizations.
This file contains a set of helper functions used by VFPProd / VFPInj.
Definition blackoilbioeffectsmodules.hh:45
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:233
constexpr auto getPropValue()
get the value data member of a property
Definition propertysystem.hh:240
This file provides the infrastructure to retrieve run-time parameters.
auto Get(bool errorIfNotRegistered=true)
Retrieve a runtime parameter.
Definition parametersystem.hpp:187
The Opm property system, traits with inheritance.
Struct holding the parameters for VtkCompositionModule.
Definition vtkcompositionparams.hpp:49
static void registerParameters()
Registers the parameters in parameter system.
Definition vtkcompositionparams.cpp:31
VTK output module for the fluid composition.
Simplifies writing multi-file VTK datasets.