27#ifndef EWOMS_VTK_BLACK_OIL_POLYMER_MODULE_HH
28#define EWOMS_VTK_BLACK_OIL_POLYMER_MODULE_HH
30#include <dune/common/fvector.hh>
32#include <opm/material/densead/Math.hpp>
69template <
class TypeTag>
80 static const int vtkFormat = getPropValue<TypeTag, Properties::VtkOutputFormat>();
83 enum { enablePolymer = getPropValue<TypeTag, Properties::EnablePolymer>() };
101 Parameters::Register<Parameters::VtkWritePolymerConcentration>
102 (
"Include the concentration of the polymer component in the water phase "
103 "in the VTK output files");
104 Parameters::Register<Parameters::VtkWritePolymerDeadPoreVolume>
105 (
"Include the fraction of the \"dead\" pore volume "
106 "in the VTK output files");
107 Parameters::Register<Parameters::VtkWritePolymerRockDensity>
108 (
"Include the amount of already adsorbed polymer component"
109 "in the VTK output files");
110 Parameters::Register<Parameters::VtkWritePolymerAdsorption>
111 (
"Include the adsorption rate of the polymer component"
112 "in the VTK output files");
113 Parameters::Register<Parameters::VtkWritePolymerViscosityCorrection>
114 (
"Include the viscosity correction of the polymer component "
115 "in the VTK output files");
116 Parameters::Register<Parameters::VtkWriteWaterViscosityCorrection>
117 (
"Include the viscosity correction of the water component "
118 "due to polymers in the VTK output files");
127 if (!Parameters::Get<Parameters::EnableVtkOutput>()) {
134 if (polymerConcentrationOutput_())
136 if (polymerDeadPoreVolumeOutput_())
138 if (polymerRockDensityOutput_())
140 if (polymerAdsorptionOutput_())
142 if (polymerViscosityCorrectionOutput_())
144 if (waterViscosityCorrectionOutput_())
154 if (!Parameters::Get<Parameters::EnableVtkOutput>()) {
161 for (
unsigned dofIdx = 0; dofIdx < elemCtx.numPrimaryDof(0); ++dofIdx) {
162 const auto& intQuants = elemCtx.intensiveQuantities(dofIdx, 0);
163 unsigned globalDofIdx = elemCtx.globalSpaceIndex(dofIdx, 0);
165 if (polymerConcentrationOutput_())
166 polymerConcentration_[globalDofIdx] =
167 scalarValue(intQuants.polymerConcentration());
169 if (polymerDeadPoreVolumeOutput_())
170 polymerDeadPoreVolume_[globalDofIdx] =
171 scalarValue(intQuants.polymerDeadPoreVolume());
173 if (polymerRockDensityOutput_())
174 polymerRockDensity_[globalDofIdx] =
175 scalarValue(intQuants.polymerRockDensity());
177 if (polymerAdsorptionOutput_())
178 polymerAdsorption_[globalDofIdx] =
179 scalarValue(intQuants.polymerAdsorption());
181 if (polymerViscosityCorrectionOutput_())
182 polymerViscosityCorrection_[globalDofIdx] =
183 scalarValue(intQuants.polymerViscosityCorrection());
185 if (waterViscosityCorrectionOutput_())
186 waterViscosityCorrection_[globalDofIdx] =
187 scalarValue(intQuants.waterViscosityCorrection());
203 if (polymerConcentrationOutput_())
206 if (polymerDeadPoreVolumeOutput_())
207 this->
commitScalarBuffer_(baseWriter,
"dead pore volume fraction", polymerDeadPoreVolume_);
209 if (polymerRockDensityOutput_())
212 if (polymerAdsorptionOutput_())
215 if (polymerViscosityCorrectionOutput_())
216 this->
commitScalarBuffer_(baseWriter,
"polymer viscosity correction", polymerViscosityCorrection_);
218 if (waterViscosityCorrectionOutput_())
219 this->
commitScalarBuffer_(baseWriter,
"water viscosity correction", waterViscosityCorrection_);
223 static bool polymerConcentrationOutput_()
225 static bool val = Parameters::Get<Parameters::VtkWritePolymerConcentration>();
229 static bool polymerDeadPoreVolumeOutput_()
231 static bool val = Parameters::Get<Parameters::VtkWritePolymerDeadPoreVolume>();
235 static bool polymerRockDensityOutput_()
237 static bool val = Parameters::Get<Parameters::VtkWritePolymerRockDensity>();
241 static bool polymerAdsorptionOutput_()
243 static bool val = Parameters::Get<Parameters::VtkWritePolymerAdsorption>();
247 static bool polymerViscosityCorrectionOutput_()
249 static bool val = Parameters::Get<Parameters::VtkWritePolymerViscosityCorrection>();
253 static bool waterViscosityCorrectionOutput_()
255 static bool val = Parameters::Get<Parameters::VtkWritePolymerViscosityCorrection>();
259 ScalarBuffer polymerConcentration_;
260 ScalarBuffer polymerDeadPoreVolume_;
261 ScalarBuffer polymerRockDensity_;
262 ScalarBuffer polymerAdsorption_;
263 ScalarBuffer polymerViscosityCorrection_;
264 ScalarBuffer waterViscosityCorrection_;
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 black oil model's polymer related quantities.
Definition: vtkblackoilpolymermodule.hh:71
void allocBuffers()
Allocate memory for the scalar fields we would like to write to the VTK file.
Definition: vtkblackoilpolymermodule.hh:125
static void registerParameters()
Register all run-time parameters for the multi-phase VTK output module.
Definition: vtkblackoilpolymermodule.hh:96
void processElement(const ElementContext &elemCtx)
Modify the internal buffers according to the intensive quantities relevant for an element.
Definition: vtkblackoilpolymermodule.hh:152
VtkBlackOilPolymerModule(const Simulator &simulator)
Definition: vtkblackoilpolymermodule.hh:88
void commitBuffers(BaseOutputWriter &baseWriter)
Add all buffers to the VTK output writer.
Definition: vtkblackoilpolymermodule.hh:194
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
The generic type tag for problems using the immiscible multi-phase model.
Definition: blackoilmodel.hh:74
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: vtkblackoilpolymermodule.hh:59
static constexpr bool value
Definition: vtkblackoilpolymermodule.hh:59
Definition: vtkblackoilpolymermodule.hh:54
static constexpr bool value
Definition: vtkblackoilpolymermodule.hh:54
Definition: vtkblackoilpolymermodule.hh:55
static constexpr bool value
Definition: vtkblackoilpolymermodule.hh:55
Definition: vtkblackoilpolymermodule.hh:58
static constexpr bool value
Definition: vtkblackoilpolymermodule.hh:58
Definition: vtkblackoilpolymermodule.hh:56
static constexpr bool value
Definition: vtkblackoilpolymermodule.hh:56
Definition: vtkblackoilpolymermodule.hh:57
static constexpr bool value
Definition: vtkblackoilpolymermodule.hh:57
Definition: vtkblackoilpolymermodule.hh:47