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 <opm/material/common/MathToolbox.hpp>
31
34
36
39
40namespace Opm::Parameters {
41
42// set default values for what quantities to output
43struct VtkWriteSolidInternalEnergy { static constexpr bool value = false; };
44struct VtkWriteThermalConductivity { static constexpr bool value = false; };
45struct VtkWriteInternalEnergies { static constexpr bool value = false; };
46struct VtkWriteEnthalpies { static constexpr bool value = false; };
47
48} // namespace Opm::Parameters
49
50namespace Opm {
64template <class TypeTag>
65class VtkEnergyModule : public BaseOutputModule<TypeTag>
66{
68
74
76 using PhaseBuffer = typename ParentType::PhaseBuffer;
77
78 static const int vtkFormat = getPropValue<TypeTag, Properties::VtkOutputFormat>();
79 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
80
81 using Toolbox = typename Opm::MathToolbox<Evaluation>;
83
84public:
85 VtkEnergyModule(const Simulator& simulator)
86 : ParentType(simulator)
87 {
88 }
89
93 static void registerParameters()
94 {
95 Parameters::Register<Parameters::VtkWriteSolidInternalEnergy>
96 ("Include the volumetric internal energy of solid"
97 "matrix in the VTK output files");
98 Parameters::Register<Parameters::VtkWriteThermalConductivity>
99 ("Include the total thermal conductivity of the"
100 "medium in the VTK output files");
101 Parameters::Register<Parameters::VtkWriteEnthalpies>
102 ("Include the specific enthalpy of the phases in "
103 "the VTK output files");
104 Parameters::Register<Parameters::VtkWriteInternalEnergies>
105 ("Include the specific internal energy of the "
106 "phases in the VTK output files");
107 }
108
114 {
115 if (enthalpyOutput_())
116 this->resizePhaseBuffer_(enthalpy_);
117 if (internalEnergyOutput_())
118 this->resizePhaseBuffer_(internalEnergy_);
119
120 if (solidInternalEnergyOutput_())
121 this->resizeScalarBuffer_(solidInternalEnergy_);
122 if (thermalConductivityOutput_())
123 this->resizeScalarBuffer_(thermalConductivity_);
124 }
125
130 void processElement(const ElementContext& elemCtx)
131 {
132 if (!Parameters::Get<Parameters::EnableVtkOutput>()) {
133 return;
134 }
135
136 for (unsigned i = 0; i < elemCtx.numPrimaryDof(/*timeIdx=*/0); ++i) {
137 unsigned I = elemCtx.globalSpaceIndex(i, /*timeIdx=*/0);
138 const auto& intQuants = elemCtx.intensiveQuantities(i, /*timeIdx=*/0);
139 const auto& fs = intQuants.fluidState();
140
141 if (solidInternalEnergyOutput_())
142 solidInternalEnergy_[I] = Toolbox::value(intQuants.solidInternalEnergy());
143 if (thermalConductivityOutput_())
144 thermalConductivity_[I] = Toolbox::value(intQuants.thermalConductivity());
145
146 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
147 if (enthalpyOutput_())
148 enthalpy_[phaseIdx][I] = Toolbox::value(fs.enthalpy(phaseIdx));
149 if (internalEnergyOutput_())
150 internalEnergy_[phaseIdx][I] = Toolbox::value(fs.internalEnergy(phaseIdx));
151 }
152 }
153 }
154
159 {
160 VtkMultiWriter *vtkWriter = dynamic_cast<VtkMultiWriter*>(&baseWriter);
161 if (!vtkWriter) {
162 return;
163 }
164
165 if (solidInternalEnergyOutput_())
166 this->commitScalarBuffer_(baseWriter, "internalEnergySolid", solidInternalEnergy_);
167 if (thermalConductivityOutput_())
168 this->commitScalarBuffer_(baseWriter, "thermalConductivity", thermalConductivity_);
169
170 if (enthalpyOutput_())
171 this->commitPhaseBuffer_(baseWriter, "enthalpy_%s", enthalpy_);
172 if (internalEnergyOutput_())
173 this->commitPhaseBuffer_(baseWriter, "internalEnergy_%s", internalEnergy_);
174 }
175
176private:
177 static bool solidInternalEnergyOutput_()
178 {
179 static bool val = Parameters::Get<Parameters::VtkWriteSolidInternalEnergy>();
180 return val;
181 }
182
183 static bool thermalConductivityOutput_()
184 {
185 static bool val = Parameters::Get<Parameters::VtkWriteThermalConductivity>();
186 return val;
187 }
188
189 static bool enthalpyOutput_()
190 {
191 static bool val = Parameters::Get<Parameters::VtkWriteEnthalpies>();
192 return val;
193 }
194
195 static bool internalEnergyOutput_()
196 {
197 static bool val = Parameters::Get<Parameters::VtkWriteInternalEnergies>();
198 return val;
199 }
200
201 PhaseBuffer enthalpy_;
202 PhaseBuffer internalEnergy_;
203
204 ScalarBuffer thermalConductivity_;
205 ScalarBuffer solidInternalEnergy_;
206};
207
208} // namespace Opm
209
210#endif
The base class for writer modules.
Definition: baseoutputmodule.hh:67
BaseOutputWriter::ScalarBuffer ScalarBuffer
Definition: baseoutputmodule.hh:85
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 commitScalarBuffer_(BaseOutputWriter &baseWriter, const char *name, ScalarBuffer &buffer, BufferType bufferType=DofBuffer)
Add a buffer containing scalar quantities to the result file.
Definition: baseoutputmodule.hh:244
void resizeScalarBuffer_(ScalarBuffer &buffer, BufferType bufferType=DofBuffer)
Allocate the space for a buffer storing a scalar quantity.
Definition: baseoutputmodule.hh:156
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:66
void processElement(const ElementContext &elemCtx)
Modify the internal buffers according to the intensive quanties relevant for an element.
Definition: vtkenergymodule.hh:130
void allocBuffers()
Allocate memory for the scalar fields we would like to write to the VTK file.
Definition: vtkenergymodule.hh:113
void commitBuffers(BaseOutputWriter &baseWriter)
Add all buffers to the VTK output writer.
Definition: vtkenergymodule.hh:158
static void registerParameters()
Register all run-time parameters for the Vtk output module.
Definition: vtkenergymodule.hh:93
VtkEnergyModule(const Simulator &simulator)
Definition: vtkenergymodule.hh:85
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: vtkenergymodule.hh:46
static constexpr bool value
Definition: vtkenergymodule.hh:46
Definition: vtkenergymodule.hh:45
static constexpr bool value
Definition: vtkenergymodule.hh:45
Definition: vtkenergymodule.hh:43
static constexpr bool value
Definition: vtkenergymodule.hh:43
Definition: vtkenergymodule.hh:44
static constexpr bool value
Definition: vtkenergymodule.hh:44