vtkenergymodule.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*/
27#ifndef OPM_VTK_ENERGY_MODULE_HPP
28#define OPM_VTK_ENERGY_MODULE_HPP
29
30#include <opm/material/common/MathToolbox.hpp>
31
35
37
40
41namespace Opm {
42
56template <class TypeTag>
57class VtkEnergyModule : public BaseOutputModule<TypeTag>
58{
60
66
68 using PhaseBuffer = typename ParentType::PhaseBuffer;
69
70 static const int vtkFormat = getPropValue<TypeTag, Properties::VtkOutputFormat>();
71 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
72
73 using Toolbox = typename Opm::MathToolbox<Evaluation>;
75
76public:
77 VtkEnergyModule(const Simulator& simulator)
78 : ParentType(simulator)
79 {
80 params_.read();
81 }
82
86 static void registerParameters()
87 {
89 }
90
96 {
97 if (params_.enthalpyOutput_) {
98 this->resizePhaseBuffer_(enthalpy_);
99 }
100 if (params_.internalEnergyOutput_) {
101 this->resizePhaseBuffer_(internalEnergy_);
102 }
103
104 if (params_.solidInternalEnergyOutput_) {
105 this->resizeScalarBuffer_(solidInternalEnergy_);
106 }
107 if (params_.thermalConductivityOutput_) {
108 this->resizeScalarBuffer_(thermalConductivity_);
109 }
110 }
111
116 void processElement(const ElementContext& elemCtx)
117 {
118 if (!Parameters::Get<Parameters::EnableVtkOutput>()) {
119 return;
120 }
121
122 for (unsigned i = 0; i < elemCtx.numPrimaryDof(/*timeIdx=*/0); ++i) {
123 unsigned I = elemCtx.globalSpaceIndex(i, /*timeIdx=*/0);
124 const auto& intQuants = elemCtx.intensiveQuantities(i, /*timeIdx=*/0);
125 const auto& fs = intQuants.fluidState();
126
127 if (params_.solidInternalEnergyOutput_) {
128 solidInternalEnergy_[I] = Toolbox::value(intQuants.solidInternalEnergy());
129 }
130 if (params_.thermalConductivityOutput_) {
131 thermalConductivity_[I] = Toolbox::value(intQuants.thermalConductivity());
132 }
133
134 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
135 if (params_.enthalpyOutput_) {
136 enthalpy_[phaseIdx][I] = Toolbox::value(fs.enthalpy(phaseIdx));
137 }
138 if (params_.internalEnergyOutput_) {
139 internalEnergy_[phaseIdx][I] = Toolbox::value(fs.internalEnergy(phaseIdx));
140 }
141 }
142 }
143 }
144
149 {
150 VtkMultiWriter* vtkWriter = dynamic_cast<VtkMultiWriter*>(&baseWriter);
151 if (!vtkWriter) {
152 return;
153 }
154
155 if (params_.solidInternalEnergyOutput_) {
156 this->commitScalarBuffer_(baseWriter, "internalEnergySolid", solidInternalEnergy_);
157 }
158 if (params_.thermalConductivityOutput_) {
159 this->commitScalarBuffer_(baseWriter, "thermalConductivity", thermalConductivity_);
160 }
161
162 if (params_.enthalpyOutput_) {
163 this->commitPhaseBuffer_(baseWriter, "enthalpy_%s", enthalpy_);
164 }
165 if (params_.internalEnergyOutput_) {
166 this->commitPhaseBuffer_(baseWriter, "internalEnergy_%s", internalEnergy_);
167 }
168 }
169
170private:
171 VtkEnergyParams params_{};
172 PhaseBuffer enthalpy_{};
173 PhaseBuffer internalEnergy_{};
174
175 ScalarBuffer thermalConductivity_{};
176 ScalarBuffer solidInternalEnergy_{};
177};
178
179} // namespace Opm
180
181#endif // OPM_VTK_ENERGY_MODULE_HPP
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.hpp:58
void processElement(const ElementContext &elemCtx)
Modify the internal buffers according to the intensive quanties relevant for an element.
Definition: vtkenergymodule.hpp:116
void allocBuffers()
Allocate memory for the scalar fields we would like to write to the VTK file.
Definition: vtkenergymodule.hpp:95
void commitBuffers(BaseOutputWriter &baseWriter)
Add all buffers to the VTK output writer.
Definition: vtkenergymodule.hpp:148
static void registerParameters()
Register all run-time parameters for the Vtk output module.
Definition: vtkenergymodule.hpp:86
VtkEnergyModule(const Simulator &simulator)
Definition: vtkenergymodule.hpp:77
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 VtkEnergyModule.
Definition: vtkenergyparams.hpp:46
static void registerParameters()
Registers the parameters in parameter system.
bool solidInternalEnergyOutput_
Definition: vtkenergyparams.hpp:53
bool thermalConductivityOutput_
Definition: vtkenergyparams.hpp:54
bool internalEnergyOutput_
Definition: vtkenergyparams.hpp:56
bool enthalpyOutput_
Definition: vtkenergyparams.hpp:55
void read()
Reads the parameter values from the parameter system.