28#ifndef OPM_VTK_DIFFUSION_MODULE_HPP
29#define OPM_VTK_DIFFUSION_MODULE_HPP
31#include <opm/material/densead/Evaluation.hpp>
32#include <opm/material/densead/Math.hpp>
56template <
class TypeTag>
57class VtkDiffusionModule :
public BaseOutputModule<TypeTag>
59 using ParentType = BaseOutputModule<TypeTag>;
66 using Toolbox = MathToolbox<Evaluation>;
69 using PhaseBuffer =
typename ParentType::PhaseBuffer;
70 using PhaseComponentBuffer =
typename ParentType::PhaseComponentBuffer;
79 explicit VtkDiffusionModule(
const Simulator& simulator)
80 : ParentType(simulator)
99 if (params_.tortuosityOutput_) {
102 if (params_.diffusionCoefficientOutput_) {
105 if (params_.effectiveDiffusionCoefficientOutput_) {
120 for (
unsigned i = 0; i < elemCtx.numPrimaryDof(0); ++i) {
121 const unsigned I = elemCtx.globalSpaceIndex(i, 0);
122 const auto& intQuants = elemCtx.intensiveQuantities(i, 0);
124 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
125 if (params_.tortuosityOutput_) {
126 tortuosity_[phaseIdx][I] = Toolbox::value(intQuants.tortuosity(phaseIdx));
128 for (
unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
129 if (params_.diffusionCoefficientOutput_) {
130 diffusionCoefficient_[phaseIdx][compIdx][I] =
131 Toolbox::value(intQuants.diffusionCoefficient(phaseIdx, compIdx));
133 if (params_.effectiveDiffusionCoefficientOutput_) {
134 effectiveDiffusionCoefficient_[phaseIdx][compIdx][I] =
135 Toolbox::value(intQuants.effectiveDiffusionCoefficient(phaseIdx, compIdx));
147 if (!
dynamic_cast<VtkMultiWriter*
>(&baseWriter)) {
151 if (params_.tortuosityOutput_) {
154 if (params_.diffusionCoefficientOutput_) {
156 diffusionCoefficient_, BufferType::Dof);
158 if (params_.effectiveDiffusionCoefficientOutput_) {
160 "effectiveDiffusionCoefficient",
161 effectiveDiffusionCoefficient_,
168 PhaseBuffer tortuosity_{};
169 PhaseComponentBuffer diffusionCoefficient_{};
170 PhaseComponentBuffer effectiveDiffusionCoefficient_{};
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
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 commitPhaseBuffer_(BaseOutputWriter &baseWriter, const char *pattern, PhaseBuffer &buffer, BufferType bufferType)
Add a phase-specific buffer to the result file.
Definition baseoutputmodule.hh:337
void resizePhaseBuffer_(PhaseBuffer &buffer, BufferType bufferType)
Allocate the space for a buffer storing a phase-specific quantity.
Definition baseoutputmodule.hh:198
The base class for all output writers.
Definition baseoutputwriter.hh:46
void allocBuffers() override
Allocate memory for the scalar fields we would like to write to the VTK file.
Definition vtkdiffusionmodule.hpp:97
void commitBuffers(BaseOutputWriter &baseWriter) override
Add all buffers to the VTK output writer.
Definition vtkdiffusionmodule.hpp:145
void processElement(const ElementContext &elemCtx) override
Modify the internal buffers according to the intensive quanties relevant for an element.
Definition vtkdiffusionmodule.hpp:114
static void registerParameters()
Register all run-time parameters for the Vtk output module.
Definition vtkdiffusionmodule.hpp:88
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 VtkDiffusionModule.
Definition vtkdiffusionparams.hpp:46
static void registerParameters()
Registers the parameters in parameter system.
Definition vtkdiffusionparams.cpp:31
VTK output module for quantities which make sense for models which incorperate molecular diffusion.
Simplifies writing multi-file VTK datasets.