vtkenergymodule.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*/
27#ifndef EWOMS_VTK_ENERGY_MODULE_HH
28#define EWOMS_VTK_ENERGY_MODULE_HH
29
30#include "vtkmultiwriter.hh"
31#include "baseoutputmodule.hh"
32
35
36#include <opm/material/common/MathToolbox.hpp>
37
38namespace Opm::Properties {
39
40namespace TTag {
41
42// create new type tag for the VTK energy output
43struct VtkEnergy {};
44
45} // namespace TTag
46
47// create the property tags needed for the energy module
48template<class TypeTag, class MyTypeTag>
50template<class TypeTag, class MyTypeTag>
52template<class TypeTag, class MyTypeTag>
54template<class TypeTag, class MyTypeTag>
56
57// set default values for what quantities to output
58template<class TypeTag>
59struct VtkWriteSolidInternalEnergy<TypeTag, TTag::VtkEnergy> { static constexpr bool value = false; };
60template<class TypeTag>
61struct VtkWriteThermalConductivity<TypeTag, TTag::VtkEnergy> { static constexpr bool value = false; };
62template<class TypeTag>
63struct VtkWriteInternalEnergies<TypeTag, TTag::VtkEnergy> { static constexpr bool value = false; };
64template<class TypeTag>
65struct VtkWriteEnthalpies<TypeTag, TTag::VtkEnergy> { static constexpr bool value = false; };
66
67} // namespace Opm::Properties
68
69namespace Opm {
83template <class TypeTag>
84class VtkEnergyModule : public BaseOutputModule<TypeTag>
85{
87
93
95 using PhaseBuffer = typename ParentType::PhaseBuffer;
96
97 static const int vtkFormat = getPropValue<TypeTag, Properties::VtkOutputFormat>();
98 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
99
100 using Toolbox = typename Opm::MathToolbox<Evaluation>;
102
103public:
104 VtkEnergyModule(const Simulator& simulator)
105 : ParentType(simulator)
106 {
107 }
108
112 static void registerParameters()
113 {
114 Parameters::registerParam<TypeTag, Properties::VtkWriteSolidInternalEnergy>
115 ("Include the volumetric internal energy of solid"
116 "matrix in the VTK output files");
117 Parameters::registerParam<TypeTag, Properties::VtkWriteThermalConductivity>
118 ("Include the total thermal conductivity of the"
119 "medium in the VTK output files");
120 Parameters::registerParam<TypeTag, Properties::VtkWriteEnthalpies>
121 ("Include the specific enthalpy of the phases in "
122 "the VTK output files");
123 Parameters::registerParam<TypeTag, Properties::VtkWriteInternalEnergies>
124 ("Include the specific internal energy of the "
125 "phases in the VTK output files");
126 }
127
133 {
134 if (enthalpyOutput_())
135 this->resizePhaseBuffer_(enthalpy_);
136 if (internalEnergyOutput_())
137 this->resizePhaseBuffer_(internalEnergy_);
138
139 if (solidInternalEnergyOutput_())
140 this->resizeScalarBuffer_(solidInternalEnergy_);
141 if (thermalConductivityOutput_())
142 this->resizeScalarBuffer_(thermalConductivity_);
143 }
144
149 void processElement(const ElementContext& elemCtx)
150 {
151 if (!Parameters::get<TypeTag, Properties::EnableVtkOutput>())
152 return;
153
154 for (unsigned i = 0; i < elemCtx.numPrimaryDof(/*timeIdx=*/0); ++i) {
155 unsigned I = elemCtx.globalSpaceIndex(i, /*timeIdx=*/0);
156 const auto& intQuants = elemCtx.intensiveQuantities(i, /*timeIdx=*/0);
157 const auto& fs = intQuants.fluidState();
158
159 if (solidInternalEnergyOutput_())
160 solidInternalEnergy_[I] = Toolbox::value(intQuants.solidInternalEnergy());
161 if (thermalConductivityOutput_())
162 thermalConductivity_[I] = Toolbox::value(intQuants.thermalConductivity());
163
164 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
165 if (enthalpyOutput_())
166 enthalpy_[phaseIdx][I] = Toolbox::value(fs.enthalpy(phaseIdx));
167 if (internalEnergyOutput_())
168 internalEnergy_[phaseIdx][I] = Toolbox::value(fs.internalEnergy(phaseIdx));
169 }
170 }
171 }
172
177 {
178 VtkMultiWriter *vtkWriter = dynamic_cast<VtkMultiWriter*>(&baseWriter);
179 if (!vtkWriter) {
180 return;
181 }
182
183 if (solidInternalEnergyOutput_())
184 this->commitScalarBuffer_(baseWriter, "internalEnergySolid", solidInternalEnergy_);
185 if (thermalConductivityOutput_())
186 this->commitScalarBuffer_(baseWriter, "thermalConductivity", thermalConductivity_);
187
188 if (enthalpyOutput_())
189 this->commitPhaseBuffer_(baseWriter, "enthalpy_%s", enthalpy_);
190 if (internalEnergyOutput_())
191 this->commitPhaseBuffer_(baseWriter, "internalEnergy_%s", internalEnergy_);
192 }
193
194private:
195 static bool solidInternalEnergyOutput_()
196 {
197 static bool val = Parameters::get<TypeTag, Properties::VtkWriteSolidInternalEnergy>();
198 return val;
199 }
200
201 static bool thermalConductivityOutput_()
202 {
203 static bool val = Parameters::get<TypeTag, Properties::VtkWriteThermalConductivity>();
204 return val;
205 }
206
207 static bool enthalpyOutput_()
208 {
209 static bool val = Parameters::get<TypeTag, Properties::VtkWriteEnthalpies>();
210 return val;
211 }
212
213 static bool internalEnergyOutput_()
214 {
215 static bool val = Parameters::get<TypeTag, Properties::VtkWriteInternalEnergies>();
216 return val;
217 }
218
219 PhaseBuffer enthalpy_;
220 PhaseBuffer internalEnergy_;
221
222 ScalarBuffer thermalConductivity_;
223 ScalarBuffer solidInternalEnergy_;
224};
225
226} // namespace Opm
227
228#endif
The base class for writer modules.
Definition: baseoutputmodule.hh:73
BaseOutputWriter::ScalarBuffer ScalarBuffer
Definition: baseoutputmodule.hh:91
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 commitScalarBuffer_(BaseOutputWriter &baseWriter, const char *name, ScalarBuffer &buffer, BufferType bufferType=DofBuffer)
Add a buffer containing scalar quantities to the result file.
Definition: baseoutputmodule.hh:316
void resizeScalarBuffer_(ScalarBuffer &buffer, BufferType bufferType=DofBuffer)
Allocate the space for a buffer storing a scalar quantity.
Definition: baseoutputmodule.hh:162
The base class for all output writers.
Definition: baseoutputwriter.hh:44
VTK output module for quantities which make sense for models which assume thermal equilibrium.
Definition: vtkenergymodule.hh:85
void processElement(const ElementContext &elemCtx)
Modify the internal buffers according to the intensive quanties relevant for an element.
Definition: vtkenergymodule.hh:149
void allocBuffers()
Allocate memory for the scalar fields we would like to write to the VTK file.
Definition: vtkenergymodule.hh:132
void commitBuffers(BaseOutputWriter &baseWriter)
Add all buffers to the VTK output writer.
Definition: vtkenergymodule.hh:176
static void registerParameters()
Register all run-time parameters for the Vtk output module.
Definition: vtkenergymodule.hh:112
VtkEnergyModule(const Simulator &simulator)
Definition: vtkenergymodule.hh:104
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: vtkenergymodule.hh:43
a tag to mark properties as undefined
Definition: propertysystem.hh:40
Definition: vtkenergymodule.hh:55
Definition: vtkenergymodule.hh:53
Definition: vtkenergymodule.hh:49
Definition: vtkenergymodule.hh:51