vtkblackoilenergymodule.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_BLACK_OIL_ENERGY_MODULE_HH
28#define EWOMS_VTK_BLACK_OIL_ENERGY_MODULE_HH
29
30#include <opm/material/densead/Math.hpp>
31
32#include "vtkmultiwriter.hh"
33#include "baseoutputmodule.hh"
34
38
39#include <dune/common/fvector.hh>
40
41#include <cstdio>
42
43namespace Opm::Properties {
44
45namespace TTag {
46
47// create new type tag for the VTK multi-phase output
49
50} // namespace TTag
51
52// create the property tags needed for the energy module
53template<class TypeTag, class MyTypeTag>
55template<class TypeTag, class MyTypeTag>
57template<class TypeTag, class MyTypeTag>
59template<class TypeTag, class MyTypeTag>
61
62// set default values for what quantities to output
63template<class TypeTag>
64struct VtkWriteRockInternalEnergy<TypeTag, TTag::VtkBlackOilEnergy> { static constexpr bool value = true; };
65template<class TypeTag>
66struct VtkWriteTotalThermalConductivity<TypeTag, TTag::VtkBlackOilEnergy> { static constexpr bool value = true; };
67template<class TypeTag>
68struct VtkWriteFluidInternalEnergies<TypeTag, TTag::VtkBlackOilEnergy> { static constexpr bool value = true; };
69template<class TypeTag>
70struct VtkWriteFluidEnthalpies<TypeTag, TTag::VtkBlackOilEnergy> { static constexpr bool value = true; };
71
72} // namespace Opm::Properties
73
74namespace Opm {
80template <class TypeTag>
82{
84
91
92 static const int vtkFormat = getPropValue<TypeTag, Properties::VtkOutputFormat>();
94
95 enum { enableEnergy = getPropValue<TypeTag, Properties::EnableEnergy>() };
96 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
97
99 using PhaseBuffer = typename ParentType::PhaseBuffer;
100
101public:
102 VtkBlackOilEnergyModule(const Simulator& simulator)
103 : ParentType(simulator)
104 { }
105
110 static void registerParameters()
111 {
112 if (!enableEnergy)
113 return;
114
115 Parameters::registerParam<TypeTag, Properties::VtkWriteRockInternalEnergy>
116 ("Include the volumetric internal energy of rock "
117 "in the VTK output files");
118 Parameters::registerParam<TypeTag, Properties::VtkWriteTotalThermalConductivity>
119 ("Include the total thermal conductivity of the medium and the fluids "
120 "in the VTK output files");
121 Parameters::registerParam<TypeTag, Properties::VtkWriteFluidInternalEnergies>
122 ("Include the internal energies of the fluids in the VTK output files");
123 Parameters::registerParam<TypeTag, Properties::VtkWriteFluidEnthalpies>
124 ("Include the enthalpies of the fluids in the VTK output files");
125 }
126
132 {
133 if (!Parameters::get<TypeTag, Properties::EnableVtkOutput>())
134 return;
135
136 if (!enableEnergy)
137 return;
138
139 if (rockInternalEnergyOutput_())
140 this->resizeScalarBuffer_(rockInternalEnergy_);
141 if (totalThermalConductivityOutput_())
142 this->resizeScalarBuffer_(totalThermalConductivity_);
143 if (fluidInternalEnergiesOutput_())
144 this->resizePhaseBuffer_(fluidInternalEnergies_);
145 if (fluidEnthalpiesOutput_())
146 this->resizePhaseBuffer_(fluidEnthalpies_);
147 }
148
153 void processElement(const ElementContext& elemCtx)
154 {
155 if (!Parameters::get<TypeTag, Properties::EnableVtkOutput>())
156 return;
157
158 if (!enableEnergy)
159 return;
160
161 for (unsigned dofIdx = 0; dofIdx < elemCtx.numPrimaryDof(/*timeIdx=*/0); ++dofIdx) {
162 const auto& intQuants = elemCtx.intensiveQuantities(dofIdx, /*timeIdx=*/0);
163 unsigned globalDofIdx = elemCtx.globalSpaceIndex(dofIdx, /*timeIdx=*/0);
164
165 if (rockInternalEnergyOutput_())
166 rockInternalEnergy_[globalDofIdx] =
167 scalarValue(intQuants.rockInternalEnergy());
168
169 if (totalThermalConductivityOutput_())
170 totalThermalConductivity_[globalDofIdx] =
171 scalarValue(intQuants.totalThermalConductivity());
172
173 for (int phaseIdx = 0; phaseIdx < numPhases; ++ phaseIdx) {
174 if (FluidSystem::phaseIsActive(phaseIdx)) {
175 if (fluidInternalEnergiesOutput_())
176 fluidInternalEnergies_[phaseIdx][globalDofIdx] =
177 scalarValue(intQuants.fluidState().internalEnergy(phaseIdx));
178
179 if (fluidEnthalpiesOutput_())
180 fluidEnthalpies_[phaseIdx][globalDofIdx] =
181 scalarValue(intQuants.fluidState().enthalpy(phaseIdx));
182 }
183 }
184 }
185 }
186
191 {
192 VtkMultiWriter *vtkWriter = dynamic_cast<VtkMultiWriter*>(&baseWriter);
193 if (!vtkWriter)
194 return;
195
196 if (!enableEnergy)
197 return;
198
199 if (rockInternalEnergyOutput_())
200 this->commitScalarBuffer_(baseWriter, "volumetric internal energy rock", rockInternalEnergy_);
201
202 if (totalThermalConductivityOutput_())
203 this->commitScalarBuffer_(baseWriter, "total thermal conductivity", totalThermalConductivity_);
204
205 if (fluidInternalEnergiesOutput_())
206 this->commitPhaseBuffer_(baseWriter, "internal energy_%s", fluidInternalEnergies_);
207
208 if (fluidEnthalpiesOutput_())
209 this->commitPhaseBuffer_(baseWriter, "enthalpy_%s", fluidEnthalpies_);
210 }
211
212private:
213 static bool rockInternalEnergyOutput_()
214 {
215 static bool val = Parameters::get<TypeTag, Properties::VtkWriteRockInternalEnergy>();
216 return val;
217 }
218
219 static bool totalThermalConductivityOutput_()
220 {
221 static bool val = Parameters::get<TypeTag, Properties::VtkWriteTotalThermalConductivity>();
222 return val;
223 }
224
225 static bool fluidInternalEnergiesOutput_()
226 {
227 static bool val = Parameters::get<TypeTag, Properties::VtkWriteFluidInternalEnergies>();
228 return val;
229 }
230
231 static bool fluidEnthalpiesOutput_()
232 {
233 static bool val = Parameters::get<TypeTag, Properties::VtkWriteFluidEnthalpies>();
234 return val;
235 }
236
237 ScalarBuffer rockInternalEnergy_;
238 ScalarBuffer totalThermalConductivity_;
239 PhaseBuffer fluidInternalEnergies_;
240 PhaseBuffer fluidEnthalpies_;
241};
242} // namespace Opm
243
244#endif
Declares the properties required by the black oil model.
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 the black oil model's energy related quantities.
Definition: vtkblackoilenergymodule.hh:82
void allocBuffers()
Allocate memory for the scalar fields we would like to write to the VTK file.
Definition: vtkblackoilenergymodule.hh:131
VtkBlackOilEnergyModule(const Simulator &simulator)
Definition: vtkblackoilenergymodule.hh:102
void commitBuffers(BaseOutputWriter &baseWriter)
Add all buffers to the VTK output writer.
Definition: vtkblackoilenergymodule.hh:190
static void registerParameters()
Register all run-time parameters for the multi-phase VTK output module.
Definition: vtkblackoilenergymodule.hh:110
void processElement(const ElementContext &elemCtx)
Modify the internal buffers according to the intensive quantities relevant for an element.
Definition: vtkblackoilenergymodule.hh:153
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: vtkblackoilenergymodule.hh:48
a tag to mark properties as undefined
Definition: propertysystem.hh:40
Definition: vtkblackoilenergymodule.hh:60
Definition: vtkblackoilenergymodule.hh:58
Definition: vtkblackoilenergymodule.hh:54
Definition: vtkblackoilenergymodule.hh:56