27#ifndef OPM_VTK_MULTI_PHASE_MODULE_HPP
28#define OPM_VTK_MULTI_PHASE_MODULE_HPP
30#include <dune/common/fvector.hh>
32#include <opm/material/common/MathToolbox.hpp>
33#include <opm/material/common/Valgrind.hpp>
71template<
class TypeTag>
87 enum { dimWorld = GridView::dimensionworld };
91 using ScalarBuffer =
typename ParentType::ScalarBuffer;
92 using VectorBuffer =
typename ParentType::VectorBuffer;
93 using TensorBuffer =
typename ParentType::TensorBuffer;
94 using PhaseBuffer =
typename ParentType::PhaseBuffer;
96 using DimVector = Dune::FieldVector<Scalar, dimWorld>;
98 using PhaseVectorBuffer = std::array<VectorBuffer, numPhases>;
121 if (params_.extrusionFactorOutput_) {
124 if (params_.pressureOutput_) {
127 if (params_.densityOutput_) {
130 if (params_.saturationOutput_) {
133 if (params_.mobilityOutput_) {
136 if (params_.relativePermeabilityOutput_) {
139 if (params_.viscosityOutput_) {
142 if (params_.averageMolarMassOutput_) {
146 if (params_.porosityOutput_) {
149 if (params_.intrinsicPermeabilityOutput_) {
153 if (params_.velocityOutput_) {
154 const std::size_t
nDof = this->simulator_.model().numGridDof();
157 for (
unsigned dofIdx = 0; dofIdx <
nDof; ++ dofIdx) {
158 velocity_[
phaseIdx][dofIdx].resize(dimWorld);
165 if (params_.potentialGradientOutput_) {
166 const std::size_t
nDof = this->simulator_.model().numGridDof();
169 for (
unsigned dofIdx = 0; dofIdx <
nDof; ++ dofIdx) {
170 potentialGradient_[
phaseIdx][dofIdx].resize(dimWorld);
171 potentialGradient_[
phaseIdx][dofIdx] = 0.0;
185 if (!Parameters::Get<Parameters::EnableVtkOutput>()) {
189 const auto& problem = elemCtx.problem();
190 for (
unsigned i = 0; i < elemCtx.numPrimaryDof(0); ++i) {
191 const unsigned I = elemCtx.globalSpaceIndex(i, 0);
192 const auto& intQuants = elemCtx.intensiveQuantities(i, 0);
193 const auto& fs = intQuants.fluidState();
195 if (params_.extrusionFactorOutput_) {
196 extrusionFactor_[I] = intQuants.extrusionFactor();
198 if (params_.porosityOutput_) {
199 porosity_[I] =
getValue(intQuants.porosity());
202 if (params_.intrinsicPermeabilityOutput_) {
203 const auto& K = problem.intrinsicPermeability(elemCtx, i, 0);
212 if (!FluidSystem::phaseIsActive(
phaseIdx)) {
215 if (params_.pressureOutput_) {
218 if (params_.densityOutput_) {
221 if (params_.saturationOutput_) {
224 if (params_.mobilityOutput_) {
227 if (params_.relativePermeabilityOutput_) {
230 if (params_.viscosityOutput_) {
233 if (params_.averageMolarMassOutput_) {
239 if (params_.potentialGradientOutput_) {
241 for (
unsigned faceIdx = 0; faceIdx < elemCtx.numInteriorFaces(0); ++faceIdx) {
242 const auto&
extQuants = elemCtx.extensiveQuantities(faceIdx, 0);
244 const unsigned i =
extQuants.interiorIndex();
245 const unsigned I = elemCtx.globalSpaceIndex(i, 0);
248 const Scalar weight =
extQuants.extrusionFactor();
250 potentialWeight_[
phaseIdx][I] += weight;
262 if (params_.velocityOutput_) {
264 for (
unsigned faceIdx = 0; faceIdx < elemCtx.numInteriorFaces(0); ++faceIdx) {
265 const auto&
extQuants = elemCtx.extensiveQuantities(faceIdx, 0);
267 const unsigned i =
extQuants.interiorIndex();
268 const unsigned I = elemCtx.globalSpaceIndex(i, 0);
270 const unsigned j =
extQuants.exteriorIndex();
271 const unsigned J = elemCtx.globalSpaceIndex(j, 0);
274 Scalar weight = std::max(Scalar{1
e-16},
276 Valgrind::CheckDefined(
extQuants.extrusionFactor());
282 for (
unsigned k = 0;
k < dimWorld; ++
k) {
285 if (
v.two_norm() > 1
e-20) {
286 weight /=
v.two_norm();
293 velocityWeight_[
phaseIdx][I] += weight;
294 velocityWeight_[
phaseIdx][J] += weight;
309 if (params_.extrusionFactorOutput_) {
311 extrusionFactor_, BufferType::Dof);
313 if (params_.pressureOutput_) {
316 if (params_.densityOutput_) {
319 if (params_.saturationOutput_) {
322 if (params_.mobilityOutput_) {
325 if (params_.relativePermeabilityOutput_) {
327 relativePermeability_, BufferType::Dof);
329 if (params_.viscosityOutput_) {
332 if (params_.averageMolarMassOutput_) {
334 averageMolarMass_, BufferType::Dof);
337 if (params_.porosityOutput_) {
340 if (params_.intrinsicPermeabilityOutput_) {
342 intrinsicPermeability_, BufferType::Dof);
345 if (params_.velocityOutput_) {
346 const std::size_t numDof = this->simulator_.model().numGridDof();
350 for (
unsigned i = 0; i < numDof; ++i) {
355 snprintf(name, 512,
"filterVelocity_%s", FluidSystem::phaseName(
phaseIdx).data());
361 if (params_.potentialGradientOutput_) {
362 const std::size_t numDof = this->simulator_.model().numGridDof();
366 for (
unsigned i = 0; i < numDof; ++i) {
373 DiscBaseOutputModule::attachVectorDofData_(
baseWriter,
390 return params_.velocityOutput_ || params_.potentialGradientOutput_;
395 ScalarBuffer extrusionFactor_{};
396 PhaseBuffer pressure_{};
397 PhaseBuffer density_{};
398 PhaseBuffer saturation_{};
399 PhaseBuffer mobility_{};
400 PhaseBuffer relativePermeability_{};
401 PhaseBuffer viscosity_{};
402 PhaseBuffer averageMolarMass_{};
404 ScalarBuffer porosity_{};
405 TensorBuffer intrinsicPermeability_{};
407 PhaseVectorBuffer velocity_{};
408 PhaseBuffer velocityWeight_{};
410 PhaseVectorBuffer potentialGradient_{};
411 PhaseBuffer potentialWeight_{};
The base class for writer modules.
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 commitTensorBuffer_(BaseOutputWriter &baseWriter, const char *name, TensorBuffer &buffer, BufferType bufferType)
Add a buffer containing tensorial quantities to the result file.
Definition baseoutputmodule.hh:282
void resizeTensorBuffer_(TensorBuffer &buffer, BufferType bufferType)
Allocate the space for a buffer storing a tensorial quantity.
Definition baseoutputmodule.hh:166
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 commitScalarBuffer_(BaseOutputWriter &baseWriter, const char *name, ScalarBuffer &buffer, BufferType bufferType)
Add a buffer containing scalar quantities to the result file.
Definition baseoutputmodule.hh:238
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
VTK output module for quantities which make sense for all models which deal with multiple fluid phase...
Definition vtkmultiphasemodule.hpp:73
static void registerParameters()
Register all run-time parameters for the multi-phase VTK output module.
Definition vtkmultiphasemodule.hpp:110
bool needExtensiveQuantities() const override
Returns true iff the module needs to access the extensive quantities of a context to do its job.
Definition vtkmultiphasemodule.hpp:388
void commitBuffers(BaseOutputWriter &baseWriter) override
Add all buffers to the VTK output writer.
Definition vtkmultiphasemodule.hpp:303
void allocBuffers() override
Allocate memory for the scalar fields we would like to write to the VTK file.
Definition vtkmultiphasemodule.hpp:119
void processElement(const ElementContext &elemCtx) override
Modify the internal buffers according to the intensive quantities seen on an element.
Definition vtkmultiphasemodule.hpp:183
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: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.
Struct holding the parameters for VtkMultiPhaseModule.
Definition vtkmultiphaseparams.hpp:54
static void registerParameters()
Registers the parameters in parameter system.
Definition vtkmultiphaseparams.cpp:31
void read()
Reads the parameter values from the parameter system.
Definition vtkmultiphaseparams.cpp:59
VTK output module for quantities which make sense for all models which deal with multiple fluid phase...
Simplifies writing multi-file VTK datasets.