27#ifndef OPM_VTK_BLACK_OIL_POLYMER_MODULE_HPP
28#define OPM_VTK_BLACK_OIL_POLYMER_MODULE_HPP
30#include <dune/common/fvector.hh>
32#include <opm/material/densead/Math.hpp>
59template <
class TypeTag>
76 using ScalarBuffer =
typename ParentType::ScalarBuffer;
82 if constexpr (enablePolymer) {
93 if constexpr (enablePolymer) {
104 if constexpr (enablePolymer) {
105 if (!Parameters::Get<Parameters::EnableVtkOutput>()) {
109 if (params_.polymerConcentrationOutput_) {
112 if (params_.polymerDeadPoreVolumeOutput_) {
115 if (params_.polymerRockDensityOutput_) {
118 if (params_.polymerAdsorptionOutput_) {
121 if (params_.polymerViscosityCorrectionOutput_) {
124 if (params_.waterViscosityCorrectionOutput_) {
136 if constexpr (enablePolymer) {
137 if (!Parameters::Get<Parameters::EnableVtkOutput>()) {
141 for (
unsigned dofIdx = 0; dofIdx < elemCtx.numPrimaryDof(0); ++dofIdx) {
142 const auto& intQuants = elemCtx.intensiveQuantities(dofIdx, 0);
143 const unsigned globalDofIdx = elemCtx.globalSpaceIndex(dofIdx, 0);
145 if (params_.polymerConcentrationOutput_) {
146 polymerConcentration_[globalDofIdx] =
150 if (params_.polymerDeadPoreVolumeOutput_) {
151 polymerDeadPoreVolume_[globalDofIdx] =
155 if (params_.polymerRockDensityOutput_) {
156 polymerRockDensity_[globalDofIdx] =
160 if (params_.polymerAdsorptionOutput_) {
161 polymerAdsorption_[globalDofIdx] =
165 if (params_.polymerViscosityCorrectionOutput_) {
166 polymerViscosityCorrection_[globalDofIdx] =
167 scalarValue(intQuants.polymerViscosityCorrection());
170 if (params_.waterViscosityCorrectionOutput_) {
171 waterViscosityCorrection_[globalDofIdx] =
183 if constexpr (enablePolymer) {
188 if (params_.polymerConcentrationOutput_) {
190 polymerConcentration_, BufferType::Dof);
193 if (params_.polymerDeadPoreVolumeOutput_) {
195 polymerDeadPoreVolume_, BufferType::Dof);
198 if (params_.polymerRockDensityOutput_) {
200 polymerRockDensity_, BufferType::Dof);
203 if (params_.polymerAdsorptionOutput_) {
205 polymerAdsorption_, BufferType::Dof);
208 if (params_.polymerViscosityCorrectionOutput_) {
210 polymerViscosityCorrection_, BufferType::Dof);
213 if (params_.waterViscosityCorrectionOutput_) {
215 waterViscosityCorrection_, BufferType::Dof);
222 ScalarBuffer polymerConcentration_{};
223 ScalarBuffer polymerDeadPoreVolume_{};
224 ScalarBuffer polymerRockDensity_{};
225 ScalarBuffer polymerAdsorption_{};
226 ScalarBuffer polymerViscosityCorrection_{};
227 ScalarBuffer waterViscosityCorrection_{};
The base class for writer modules.
Declares the properties required by the black oil model.
The base class for writer modules.
Definition baseoutputmodule.hh:68
void resizeScalarBuffer_(ScalarBuffer &buffer, BufferType bufferType)
Allocate the space for a buffer storing a scalar quantity.
Definition baseoutputmodule.hh:157
BufferType
Definition baseoutputmodule.hh:143
void commitScalarBuffer_(BaseOutputWriter &baseWriter, const char *name, ScalarBuffer &buffer, BufferType bufferType)
Add a buffer containing scalar quantities to the result file.
Definition baseoutputmodule.hh:238
The base class for all output writers.
Definition baseoutputwriter.hh:46
VTK output module for the black oil model's polymer related quantities.
Definition vtkblackoilpolymermodule.hpp:61
static void registerParameters()
Register all run-time parameters for the multi-phase VTK output module.
Definition vtkblackoilpolymermodule.hpp:91
void processElement(const ElementContext &elemCtx) override
Modify the internal buffers according to the intensive quantities relevant for an element.
Definition vtkblackoilpolymermodule.hpp:134
void allocBuffers() override
Allocate memory for the scalar fields we would like to write to the VTK file.
Definition vtkblackoilpolymermodule.hpp:102
void commitBuffers(BaseOutputWriter &baseWriter) override
Add all buffers to the VTK output writer.
Definition vtkblackoilpolymermodule.hpp:181
Simplifies writing multi-file VTK datasets.
Definition vtkmultiwriter.hh:65
Declare the properties used by the infrastructure code of the finite volume discretizations.
The generic type tag for problems using the immiscible multi-phase model.
Definition blackoilmodel.hh:81
This file contains a set of helper functions used by VFPProd / VFPInj.
Definition blackoilbioeffectsmodules.hh:43
constexpr auto getPropValue()
get the value data member of a property
Definition propertysystem.hh:240
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
This file provides the infrastructure to retrieve run-time parameters.
The Opm property system, traits with inheritance.
Definition vtkblackoilpolymermodule.hpp:48
Struct holding the parameters for VtkBlackoilPolymerModule.
Definition vtkblackoilpolymerparams.hpp:48
static void registerParameters()
Registers the parameters in parameter system.
Definition vtkblackoilpolymerparams.cpp:31
void read()
Reads the parameter values from the parameter system.
Definition vtkblackoilpolymerparams.cpp:53
VTK output module for the black oil model's polymer related quantities.
Simplifies writing multi-file VTK datasets.