vtkdiffusionmodule.hpp
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3/*
4 This file is part of the Open Porous Media project (OPM).
5
6 OPM is free software: you can redistribute it and/or modify
7 it under the terms of the GNU General Public License as published by
8 the Free Software Foundation, either version 2 of the License, or
9 (at your option) any later version.
10
11 OPM is distributed in the hope that it will be useful,
12 but WITHOUT ANY WARRANTY; without even the implied warranty of
13 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14 GNU General Public License for more details.
15
16 You should have received a copy of the GNU General Public License
17 along with OPM. If not, see <http://www.gnu.org/licenses/>.
18
19 Consult the COPYING file in the top-level source directory of this
20 module for the precise wording of the license and the list of
21 copyright holders.
22*/
28#ifndef OPM_VTK_DIFFUSION_MODULE_HPP
29#define OPM_VTK_DIFFUSION_MODULE_HPP
30
31#include <opm/material/densead/Evaluation.hpp>
32#include <opm/material/densead/Math.hpp>
33
35
39
42
43namespace Opm {
44
56template <class TypeTag>
57class VtkDiffusionModule : public BaseOutputModule<TypeTag>
58{
60
65
66 using Toolbox = MathToolbox<Evaluation>;
67
69 using PhaseBuffer = typename ParentType::PhaseBuffer;
70
71 static const int vtkFormat = getPropValue<TypeTag, Properties::VtkOutputFormat>();
73
74 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
75 enum { numComponents = getPropValue<TypeTag, Properties::NumComponents>() };
76
77public:
78 VtkDiffusionModule(const Simulator& simulator)
79 : ParentType(simulator)
80 {
81 params_.read();
82 }
83
87 static void registerParameters()
88 {
90 }
91
97 {
98 if (params_.tortuosityOutput_) {
99 this->resizePhaseBuffer_(tortuosity_);
100 }
101 if (params_.diffusionCoefficientOutput_) {
102 this->resizePhaseComponentBuffer_(diffusionCoefficient_);
103 }
105 this->resizePhaseComponentBuffer_(effectiveDiffusionCoefficient_);
106 }
107 }
108
113 void processElement(const ElementContext& elemCtx)
114 {
115 if (!Parameters::Get<Parameters::EnableVtkOutput>()) {
116 return;
117 }
118
119 for (unsigned i = 0; i < elemCtx.numPrimaryDof(/*timeIdx=*/0); ++i) {
120 unsigned I = elemCtx.globalSpaceIndex(i, /*timeIdx=*/0);
121 const auto& intQuants = elemCtx.intensiveQuantities(i, /*timeIdx=*/0);
122
123 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
124 if (params_.tortuosityOutput_) {
125 tortuosity_[phaseIdx][I] = Toolbox::value(intQuants.tortuosity(phaseIdx));
126 }
127 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
128 if (params_.diffusionCoefficientOutput_) {
129 diffusionCoefficient_[phaseIdx][compIdx][I] =
130 Toolbox::value(intQuants.diffusionCoefficient(phaseIdx, compIdx));
131 }
133 effectiveDiffusionCoefficient_[phaseIdx][compIdx][I] =
134 Toolbox::value(intQuants.effectiveDiffusionCoefficient(phaseIdx, compIdx));
135 }
136 }
137 }
138 }
139 }
140
145 {
146 VtkMultiWriter* vtkWriter = dynamic_cast<VtkMultiWriter*>(&baseWriter);
147 if (!vtkWriter) {
148 return;
149 }
150
151 if (params_.tortuosityOutput_) {
152 this->commitPhaseBuffer_(baseWriter, "tortuosity", tortuosity_);
153 }
154 if (params_.diffusionCoefficientOutput_) {
155 this->commitPhaseComponentBuffer_(baseWriter, "diffusionCoefficient",
156 diffusionCoefficient_);
157 }
159 this->commitPhaseComponentBuffer_(baseWriter,
160 "effectiveDiffusionCoefficient",
161 effectiveDiffusionCoefficient_);
162 }
163 }
164
165private:
166 VtkDiffusionParams params_{};
167 PhaseBuffer tortuosity_{};
168 PhaseComponentBuffer diffusionCoefficient_{};
169 PhaseComponentBuffer effectiveDiffusionCoefficient_{};
170};
171
172} // namespace Opm
173
174#endif // OPM_VTK_DIFFUSION_MODULE_HPP
The base class for writer modules.
Definition: baseoutputmodule.hh:67
void commitPhaseBuffer_(BaseOutputWriter &baseWriter, const char *pattern, PhaseBuffer &buffer, BufferType bufferType=DofBuffer)
Add a phase-specific buffer to the result file.
Definition: baseoutputmodule.hh:345
void resizePhaseBuffer_(PhaseBuffer &buffer, BufferType bufferType=DofBuffer)
Allocate the space for a buffer storing a phase-specific quantity.
Definition: baseoutputmodule.hh:201
std::array< ScalarBuffer, numPhases > PhaseBuffer
Definition: baseoutputmodule.hh:90
void resizePhaseComponentBuffer_(PhaseComponentBuffer &buffer, BufferType bufferType=DofBuffer)
Allocate the space for a buffer storing a phase and component specific buffer.
Definition: baseoutputmodule.hh:229
void commitPhaseComponentBuffer_(BaseOutputWriter &baseWriter, const char *pattern, PhaseComponentBuffer &buffer, BufferType bufferType=DofBuffer)
Add a phase and component specific quantities to the output.
Definition: baseoutputmodule.hh:377
std::array< std::array< ScalarBuffer, numComponents >, numPhases > PhaseComponentBuffer
Definition: baseoutputmodule.hh:92
The base class for all output writers.
Definition: baseoutputwriter.hh:44
VTK output module for quantities which make sense for models which incorperate molecular diffusion.
Definition: vtkdiffusionmodule.hpp:58
void allocBuffers()
Allocate memory for the scalar fields we would like to write to the VTK file.
Definition: vtkdiffusionmodule.hpp:96
void processElement(const ElementContext &elemCtx)
Modify the internal buffers according to the intensive quanties relevant for an element.
Definition: vtkdiffusionmodule.hpp:113
VtkDiffusionModule(const Simulator &simulator)
Definition: vtkdiffusionmodule.hpp:78
void commitBuffers(BaseOutputWriter &baseWriter)
Add all buffers to the VTK output writer.
Definition: vtkdiffusionmodule.hpp:144
static void registerParameters()
Register all run-time parameters for the Vtk output module.
Definition: vtkdiffusionmodule.hpp:87
Simplifies writing multi-file VTK datasets.
Definition: vtkmultiwriter.hh:66
Declare the properties used by the infrastructure code of the finite volume discretizations.
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.
Struct holding the parameters for VtkDiffusionModule.
Definition: vtkdiffusionparams.hpp:46
static void registerParameters()
Registers the parameters in parameter system.
bool diffusionCoefficientOutput_
Definition: vtkdiffusionparams.hpp:54
bool tortuosityOutput_
Definition: vtkdiffusionparams.hpp:53
bool effectiveDiffusionCoefficientOutput_
Definition: vtkdiffusionparams.hpp:55
void read()
Reads the parameter values from the parameter system.