vtkdiffusionmodule.hh
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 EWOMS_VTK_DIFFUSION_MODULE_HH
29#define EWOMS_VTK_DIFFUSION_MODULE_HH
30
31#include <opm/material/densead/Evaluation.hpp>
32#include <opm/material/densead/Math.hpp>
33
35
38
41
42namespace Opm::Parameters {
43
44// set default values for what quantities to output
45struct VtkWriteTortuosities { static constexpr bool value = false; };
46struct VtkWriteDiffusionCoefficients { static constexpr bool value = false; };
47struct VtkWriteEffectiveDiffusionCoefficients { static constexpr bool value = false; };
48
49} // namespace Opm::Parameters
50
51namespace Opm {
63template <class TypeTag>
64class VtkDiffusionModule : public BaseOutputModule<TypeTag>
65{
67
72
73 using Toolbox = MathToolbox<Evaluation>;
74
76 using PhaseBuffer = typename ParentType::PhaseBuffer;
77
78 static const int vtkFormat = getPropValue<TypeTag, Properties::VtkOutputFormat>();
80
81 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
82 enum { numComponents = getPropValue<TypeTag, Properties::NumComponents>() };
83
84public:
85 VtkDiffusionModule(const Simulator& simulator)
86 : ParentType(simulator)
87 { }
88
92 static void registerParameters()
93 {
94 Parameters::Register<Parameters::VtkWriteTortuosities>
95 ("Include the tortuosity for each phase in the VTK output files");
96 Parameters::Register<Parameters::VtkWriteDiffusionCoefficients>
97 ("Include the molecular diffusion coefficients in "
98 "the VTK output files");
99 Parameters::Register<Parameters::VtkWriteEffectiveDiffusionCoefficients>
100 ("Include the effective molecular diffusion "
101 "coefficients the medium in the VTK output files");
102 }
103
109 {
110 if (tortuosityOutput_())
111 this->resizePhaseBuffer_(tortuosity_);
112 if (diffusionCoefficientOutput_())
113 this->resizePhaseComponentBuffer_(diffusionCoefficient_);
114 if (effectiveDiffusionCoefficientOutput_())
115 this->resizePhaseComponentBuffer_(effectiveDiffusionCoefficient_);
116 }
117
122 void processElement(const ElementContext& elemCtx)
123 {
124 if (!Parameters::Get<Parameters::EnableVtkOutput>()) {
125 return;
126 }
127
128 for (unsigned i = 0; i < elemCtx.numPrimaryDof(/*timeIdx=*/0); ++i) {
129 unsigned I = elemCtx.globalSpaceIndex(i, /*timeIdx=*/0);
130 const auto& intQuants = elemCtx.intensiveQuantities(i, /*timeIdx=*/0);
131
132 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
133 if (tortuosityOutput_())
134 tortuosity_[phaseIdx][I] = Toolbox::value(intQuants.tortuosity(phaseIdx));
135 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
136 if (diffusionCoefficientOutput_())
137 diffusionCoefficient_[phaseIdx][compIdx][I] =
138 Toolbox::value(intQuants.diffusionCoefficient(phaseIdx, compIdx));
139 if (effectiveDiffusionCoefficientOutput_())
140 effectiveDiffusionCoefficient_[phaseIdx][compIdx][I] =
141 Toolbox::value(intQuants.effectiveDiffusionCoefficient(phaseIdx, compIdx));
142 }
143 }
144 }
145 }
146
151 {
152 VtkMultiWriter *vtkWriter = dynamic_cast<VtkMultiWriter*>(&baseWriter);
153 if (!vtkWriter) {
154 return;
155 }
156
157 if (tortuosityOutput_())
158 this->commitPhaseBuffer_(baseWriter, "tortuosity", tortuosity_);
159 if (diffusionCoefficientOutput_())
160 this->commitPhaseComponentBuffer_(baseWriter, "diffusionCoefficient",
161 diffusionCoefficient_);
162 if (effectiveDiffusionCoefficientOutput_())
163 this->commitPhaseComponentBuffer_(baseWriter,
164 "effectiveDiffusionCoefficient",
165 effectiveDiffusionCoefficient_);
166 }
167
168private:
169 static bool tortuosityOutput_()
170 {
171 static bool val = Parameters::Get<Parameters::VtkWriteTortuosities>();
172 return val;
173 }
174
175 static bool diffusionCoefficientOutput_()
176 {
177 static bool val = Parameters::Get<Parameters::VtkWriteDiffusionCoefficients>();
178 return val;
179 }
180
181 static bool effectiveDiffusionCoefficientOutput_()
182 {
183 static bool val = Parameters::Get<Parameters::VtkWriteEffectiveDiffusionCoefficients>();
184 return val;
185 }
186
187 PhaseBuffer tortuosity_;
188 PhaseComponentBuffer diffusionCoefficient_;
189 PhaseComponentBuffer effectiveDiffusionCoefficient_;
190};
191
192} // namespace Opm
193
194#endif
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.hh:65
void allocBuffers()
Allocate memory for the scalar fields we would like to write to the VTK file.
Definition: vtkdiffusionmodule.hh:108
void processElement(const ElementContext &elemCtx)
Modify the internal buffers according to the intensive quanties relevant for an element.
Definition: vtkdiffusionmodule.hh:122
VtkDiffusionModule(const Simulator &simulator)
Definition: vtkdiffusionmodule.hh:85
void commitBuffers(BaseOutputWriter &baseWriter)
Add all buffers to the VTK output writer.
Definition: vtkdiffusionmodule.hh:150
static void registerParameters()
Register all run-time parameters for the Vtk output module.
Definition: vtkdiffusionmodule.hh:92
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
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: vtkdiffusionmodule.hh:46
static constexpr bool value
Definition: vtkdiffusionmodule.hh:46
static constexpr bool value
Definition: vtkdiffusionmodule.hh:47
Definition: vtkdiffusionmodule.hh:45
static constexpr bool value
Definition: vtkdiffusionmodule.hh:45