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 "vtkmultiwriter.hh"
32#include "baseoutputmodule.hh"
33
36
37#include <opm/material/densead/Evaluation.hpp>
38#include <opm/material/densead/Math.hpp>
39
40namespace Opm::Properties {
41
42namespace TTag {
43
44// create new type tag for the VTK output of the quantities for molecular
45// diffusion
46struct VtkDiffusion {};
47
48} // namespace TTag
49
50// create the property tags needed for the diffusion module
51template<class TypeTag, class MyTypeTag>
53template<class TypeTag, class MyTypeTag>
55template<class TypeTag, class MyTypeTag>
57
58// set default values for what quantities to output
59template<class TypeTag>
60struct VtkWriteTortuosities<TypeTag, TTag::VtkDiffusion> { static constexpr bool value = false; };
61template<class TypeTag>
62struct VtkWriteDiffusionCoefficients<TypeTag, TTag::VtkDiffusion> { static constexpr bool value = false; };
63template<class TypeTag>
64struct VtkWriteEffectiveDiffusionCoefficients<TypeTag, TTag::VtkDiffusion> { static constexpr bool value = false; };
65
66} // namespace Opm::Properties
67
68namespace Opm {
80template <class TypeTag>
81class VtkDiffusionModule : public BaseOutputModule<TypeTag>
82{
84
89
90 using Toolbox = MathToolbox<Evaluation>;
91
93 using PhaseBuffer = typename ParentType::PhaseBuffer;
94
95 static const int vtkFormat = getPropValue<TypeTag, Properties::VtkOutputFormat>();
97
98 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
99 enum { numComponents = getPropValue<TypeTag, Properties::NumComponents>() };
100
101public:
102 VtkDiffusionModule(const Simulator& simulator)
103 : ParentType(simulator)
104 { }
105
109 static void registerParameters()
110 {
111 Parameters::registerParam<TypeTag, Properties::VtkWriteTortuosities>
112 ("Include the tortuosity for each phase in the VTK output files");
113 Parameters::registerParam<TypeTag, Properties::VtkWriteDiffusionCoefficients>
114 ("Include the molecular diffusion coefficients in "
115 "the VTK output files");
116 Parameters::registerParam<TypeTag, Properties::VtkWriteEffectiveDiffusionCoefficients>
117 ("Include the effective molecular diffusion "
118 "coefficients the medium in the VTK output files");
119 }
120
126 {
127 if (tortuosityOutput_())
128 this->resizePhaseBuffer_(tortuosity_);
129 if (diffusionCoefficientOutput_())
130 this->resizePhaseComponentBuffer_(diffusionCoefficient_);
131 if (effectiveDiffusionCoefficientOutput_())
132 this->resizePhaseComponentBuffer_(effectiveDiffusionCoefficient_);
133 }
134
139 void processElement(const ElementContext& elemCtx)
140 {
141 if (!Parameters::get<TypeTag, Properties::EnableVtkOutput>())
142 return;
143
144 for (unsigned i = 0; i < elemCtx.numPrimaryDof(/*timeIdx=*/0); ++i) {
145 unsigned I = elemCtx.globalSpaceIndex(i, /*timeIdx=*/0);
146 const auto& intQuants = elemCtx.intensiveQuantities(i, /*timeIdx=*/0);
147
148 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
149 if (tortuosityOutput_())
150 tortuosity_[phaseIdx][I] = Toolbox::value(intQuants.tortuosity(phaseIdx));
151 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
152 if (diffusionCoefficientOutput_())
153 diffusionCoefficient_[phaseIdx][compIdx][I] =
154 Toolbox::value(intQuants.diffusionCoefficient(phaseIdx, compIdx));
155 if (effectiveDiffusionCoefficientOutput_())
156 effectiveDiffusionCoefficient_[phaseIdx][compIdx][I] =
157 Toolbox::value(intQuants.effectiveDiffusionCoefficient(phaseIdx, compIdx));
158 }
159 }
160 }
161 }
162
167 {
168 VtkMultiWriter *vtkWriter = dynamic_cast<VtkMultiWriter*>(&baseWriter);
169 if (!vtkWriter) {
170 return;
171 }
172
173 if (tortuosityOutput_())
174 this->commitPhaseBuffer_(baseWriter, "tortuosity", tortuosity_);
175 if (diffusionCoefficientOutput_())
176 this->commitPhaseComponentBuffer_(baseWriter, "diffusionCoefficient",
177 diffusionCoefficient_);
178 if (effectiveDiffusionCoefficientOutput_())
179 this->commitPhaseComponentBuffer_(baseWriter,
180 "effectiveDiffusionCoefficient",
181 effectiveDiffusionCoefficient_);
182 }
183
184private:
185 static bool tortuosityOutput_()
186 {
187 static bool val = Parameters::get<TypeTag, Properties::VtkWriteTortuosities>();
188 return val;
189 }
190
191 static bool diffusionCoefficientOutput_()
192 {
193 static bool val = Parameters::get<TypeTag, Properties::VtkWriteDiffusionCoefficients>();
194 return val;
195 }
196
197 static bool effectiveDiffusionCoefficientOutput_()
198 {
199 static bool val = Parameters::get<TypeTag, Properties::VtkWriteEffectiveDiffusionCoefficients>();
200 return val;
201 }
202
203 PhaseBuffer tortuosity_;
204 PhaseComponentBuffer diffusionCoefficient_;
205 PhaseComponentBuffer effectiveDiffusionCoefficient_;
206};
207
208} // namespace Opm
209
210#endif
The base class for writer modules.
Definition: baseoutputmodule.hh:73
void commitPhaseBuffer_(BaseOutputWriter &baseWriter, const char *pattern, PhaseBuffer &buffer, BufferType bufferType=DofBuffer)
Add a phase-specific buffer to the result file.
Definition: baseoutputmodule.hh:419
void resizePhaseBuffer_(PhaseBuffer &buffer, BufferType bufferType=DofBuffer)
Allocate the space for a buffer storing a phase-specific quantity.
Definition: baseoutputmodule.hh:246
std::array< ScalarBuffer, numPhases > PhaseBuffer
Definition: baseoutputmodule.hh:96
void resizePhaseComponentBuffer_(PhaseComponentBuffer &buffer, BufferType bufferType=DofBuffer)
Allocate the space for a buffer storing a phase and component specific buffer.
Definition: baseoutputmodule.hh:292
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:465
std::array< std::array< ScalarBuffer, numComponents >, numPhases > PhaseComponentBuffer
Definition: baseoutputmodule.hh:98
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:82
void allocBuffers()
Allocate memory for the scalar fields we would like to write to the VTK file.
Definition: vtkdiffusionmodule.hh:125
void processElement(const ElementContext &elemCtx)
Modify the internal buffers according to the intensive quanties relevant for an element.
Definition: vtkdiffusionmodule.hh:139
VtkDiffusionModule(const Simulator &simulator)
Definition: vtkdiffusionmodule.hh:102
void commitBuffers(BaseOutputWriter &baseWriter)
Add all buffers to the VTK output writer.
Definition: vtkdiffusionmodule.hh:166
static void registerParameters()
Register all run-time parameters for the Vtk output module.
Definition: vtkdiffusionmodule.hh:109
Simplifies writing multi-file VTK datasets.
Definition: vtkmultiwriter.hh:66
Definition: blackoilmodel.hh:72
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:242
This file provides the infrastructure to retrieve run-time parameters.
The Opm property system, traits with inheritance.
Definition: vtkdiffusionmodule.hh:46
a tag to mark properties as undefined
Definition: propertysystem.hh:40
Definition: vtkdiffusionmodule.hh:54
Definition: vtkdiffusionmodule.hh:52