27#ifndef EWOMS_VTK_BLACK_OIL_MICP_MODULE_HH
28#define EWOMS_VTK_BLACK_OIL_MICP_MODULE_HH
30#include <dune/common/fvector.hh>
32#include <opm/material/densead/Math.hpp>
61template <
class TypeTag>
72 static const int vtkFormat = getPropValue<TypeTag, Properties::VtkOutputFormat>();
75 enum { enableMICP = getPropValue<TypeTag, Properties::EnableMICP>() };
93 Parameters::Register<Parameters::VtkWriteMicrobialConcentration>
94 (
"Include the concentration of the microbial component in the water phase "
95 "in the VTK output files");
96 Parameters::Register<Parameters::VtkWriteOxygenConcentration>
97 (
"Include the concentration of the oxygen component in the water phase "
98 "in the VTK output files");
99 Parameters::Register<Parameters::VtkWriteUreaConcentration>
100 (
"Include the concentration of the urea component in the water phase "
101 "in the VTK output files");
102 Parameters::Register<Parameters::VtkWriteBiofilmConcentration>
103 (
"Include the biofilm volume fraction in the VTK output files");
104 Parameters::Register<Parameters::VtkWriteCalciteConcentration>
105 (
"Include the calcite volume fraction in the VTK output files");
114 if (!Parameters::Get<Parameters::EnableVtkOutput>()) {
121 if (microbialConcentrationOutput_())
123 if (oxygenConcentrationOutput_())
125 if (ureaConcentrationOutput_())
127 if (biofilmConcentrationOutput_())
129 if (calciteConcentrationOutput_())
139 if (!Parameters::Get<Parameters::EnableVtkOutput>()) {
146 for (
unsigned dofIdx = 0; dofIdx < elemCtx.numPrimaryDof(0); ++dofIdx) {
147 const auto& intQuants = elemCtx.intensiveQuantities(dofIdx, 0);
148 unsigned globalDofIdx = elemCtx.globalSpaceIndex(dofIdx, 0);
150 if (microbialConcentrationOutput_())
151 microbialConcentration_[globalDofIdx] =
152 scalarValue(intQuants.microbialConcentration());
154 if (oxygenConcentrationOutput_())
155 oxygenConcentration_[globalDofIdx] =
156 scalarValue(intQuants.oxygenConcentration());
158 if (ureaConcentrationOutput_())
159 ureaConcentration_[globalDofIdx] =
160 10 * scalarValue(intQuants.ureaConcentration());
162 if (biofilmConcentrationOutput_())
163 biofilmConcentration_[globalDofIdx] =
164 scalarValue(intQuants.biofilmConcentration());
166 if (calciteConcentrationOutput_())
167 calciteConcentration_[globalDofIdx] =
168 scalarValue(intQuants.calciteConcentration());
185 if (microbialConcentrationOutput_())
188 if (oxygenConcentrationOutput_())
191 if (ureaConcentrationOutput_())
194 if (biofilmConcentrationOutput_())
197 if (calciteConcentrationOutput_())
203 static bool microbialConcentrationOutput_()
205 static bool val = Parameters::Get<Parameters::VtkWriteMicrobialConcentration>();
209 static bool oxygenConcentrationOutput_()
211 static bool val = Parameters::Get<Parameters::VtkWriteOxygenConcentration>();
215 static bool ureaConcentrationOutput_()
217 static bool val = Parameters::Get<Parameters::VtkWriteUreaConcentration>();
221 static bool biofilmConcentrationOutput_()
223 static bool val = Parameters::Get<Parameters::VtkWriteBiofilmConcentration>();
227 static bool calciteConcentrationOutput_()
229 static bool val = Parameters::Get<Parameters::VtkWriteCalciteConcentration>();
233 ScalarBuffer microbialConcentration_;
234 ScalarBuffer oxygenConcentration_;
235 ScalarBuffer ureaConcentration_;
236 ScalarBuffer biofilmConcentration_;
237 ScalarBuffer calciteConcentration_;
Declares the properties required by the black oil model.
The base class for writer modules.
Definition: baseoutputmodule.hh:67
BaseOutputWriter::ScalarBuffer ScalarBuffer
Definition: baseoutputmodule.hh:85
void commitScalarBuffer_(BaseOutputWriter &baseWriter, const char *name, ScalarBuffer &buffer, BufferType bufferType=DofBuffer)
Add a buffer containing scalar quantities to the result file.
Definition: baseoutputmodule.hh:244
void resizeScalarBuffer_(ScalarBuffer &buffer, BufferType bufferType=DofBuffer)
Allocate the space for a buffer storing a scalar quantity.
Definition: baseoutputmodule.hh:156
The base class for all output writers.
Definition: baseoutputwriter.hh:44
VTK output module for the MICP model's related quantities.
Definition: vtkblackoilmicpmodule.hh:63
void processElement(const ElementContext &elemCtx)
Modify the internal buffers according to the intensive quantities relevant for an element.
Definition: vtkblackoilmicpmodule.hh:137
void commitBuffers(BaseOutputWriter &baseWriter)
Add all buffers to the VTK output writer.
Definition: vtkblackoilmicpmodule.hh:176
VtkBlackOilMICPModule(const Simulator &simulator)
Definition: vtkblackoilmicpmodule.hh:80
void allocBuffers()
Allocate memory for the scalar fields we would like to write to the VTK file.
Definition: vtkblackoilmicpmodule.hh:112
static void registerParameters()
Register all run-time parameters for the multi-phase VTK output module.
Definition: vtkblackoilmicpmodule.hh:88
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: vtkblackoilmicpmodule.hh:50
static constexpr bool value
Definition: vtkblackoilmicpmodule.hh:50
Definition: vtkblackoilmicpmodule.hh:51
static constexpr bool value
Definition: vtkblackoilmicpmodule.hh:51
Definition: vtkblackoilmicpmodule.hh:47
static constexpr bool value
Definition: vtkblackoilmicpmodule.hh:47
Definition: vtkblackoilmicpmodule.hh:48
static constexpr bool value
Definition: vtkblackoilmicpmodule.hh:48
Definition: vtkblackoilmicpmodule.hh:49
static constexpr bool value
Definition: vtkblackoilmicpmodule.hh:49