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 <dune/common/fvector.hh>
31
32#include <opm/material/densead/Math.hpp>
33
35
37
40
43
44namespace Opm::Parameters {
45
46// set default values for what quantities to output
47struct VtkWriteRockInternalEnergy { static constexpr bool value = true; };
48struct VtkWriteTotalThermalConductivity { static constexpr bool value = true; };
49struct VtkWriteFluidInternalEnergies { static constexpr bool value = true; };
50struct VtkWriteFluidEnthalpies { static constexpr bool value = true; };
51
52} // namespace Opm::Parameters
53
54namespace Opm {
60template <class TypeTag>
62{
64
71
72 static const int vtkFormat = getPropValue<TypeTag, Properties::VtkOutputFormat>();
74
75 enum { enableEnergy = getPropValue<TypeTag, Properties::EnableEnergy>() };
76 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
77
79 using PhaseBuffer = typename ParentType::PhaseBuffer;
80
81public:
82 VtkBlackOilEnergyModule(const Simulator& simulator)
83 : ParentType(simulator)
84 { }
85
90 static void registerParameters()
91 {
92 if (!enableEnergy)
93 return;
94
95 Parameters::Register<Parameters::VtkWriteRockInternalEnergy>
96 ("Include the volumetric internal energy of rock "
97 "in the VTK output files");
98 Parameters::Register<Parameters::VtkWriteTotalThermalConductivity>
99 ("Include the total thermal conductivity of the medium and the fluids "
100 "in the VTK output files");
101 Parameters::Register<Parameters::VtkWriteFluidInternalEnergies>
102 ("Include the internal energies of the fluids in the VTK output files");
103 Parameters::Register<Parameters::VtkWriteFluidEnthalpies>
104 ("Include the enthalpies of the fluids in the VTK output files");
105 }
106
112 {
113 if (!Parameters::Get<Parameters::EnableVtkOutput>()) {
114 return;
115 }
116
117 if (!enableEnergy)
118 return;
119
120 if (rockInternalEnergyOutput_())
121 this->resizeScalarBuffer_(rockInternalEnergy_);
122 if (totalThermalConductivityOutput_())
123 this->resizeScalarBuffer_(totalThermalConductivity_);
124 if (fluidInternalEnergiesOutput_())
125 this->resizePhaseBuffer_(fluidInternalEnergies_);
126 if (fluidEnthalpiesOutput_())
127 this->resizePhaseBuffer_(fluidEnthalpies_);
128 }
129
134 void processElement(const ElementContext& elemCtx)
135 {
136 if (!Parameters::Get<Parameters::EnableVtkOutput>()) {
137 return;
138 }
139
140 if (!enableEnergy)
141 return;
142
143 for (unsigned dofIdx = 0; dofIdx < elemCtx.numPrimaryDof(/*timeIdx=*/0); ++dofIdx) {
144 const auto& intQuants = elemCtx.intensiveQuantities(dofIdx, /*timeIdx=*/0);
145 unsigned globalDofIdx = elemCtx.globalSpaceIndex(dofIdx, /*timeIdx=*/0);
146
147 if (rockInternalEnergyOutput_())
148 rockInternalEnergy_[globalDofIdx] =
149 scalarValue(intQuants.rockInternalEnergy());
150
151 if (totalThermalConductivityOutput_())
152 totalThermalConductivity_[globalDofIdx] =
153 scalarValue(intQuants.totalThermalConductivity());
154
155 for (int phaseIdx = 0; phaseIdx < numPhases; ++ phaseIdx) {
156 if (FluidSystem::phaseIsActive(phaseIdx)) {
157 if (fluidInternalEnergiesOutput_())
158 fluidInternalEnergies_[phaseIdx][globalDofIdx] =
159 scalarValue(intQuants.fluidState().internalEnergy(phaseIdx));
160
161 if (fluidEnthalpiesOutput_())
162 fluidEnthalpies_[phaseIdx][globalDofIdx] =
163 scalarValue(intQuants.fluidState().enthalpy(phaseIdx));
164 }
165 }
166 }
167 }
168
173 {
174 VtkMultiWriter *vtkWriter = dynamic_cast<VtkMultiWriter*>(&baseWriter);
175 if (!vtkWriter)
176 return;
177
178 if (!enableEnergy)
179 return;
180
181 if (rockInternalEnergyOutput_())
182 this->commitScalarBuffer_(baseWriter, "volumetric internal energy rock", rockInternalEnergy_);
183
184 if (totalThermalConductivityOutput_())
185 this->commitScalarBuffer_(baseWriter, "total thermal conductivity", totalThermalConductivity_);
186
187 if (fluidInternalEnergiesOutput_())
188 this->commitPhaseBuffer_(baseWriter, "internal energy_%s", fluidInternalEnergies_);
189
190 if (fluidEnthalpiesOutput_())
191 this->commitPhaseBuffer_(baseWriter, "enthalpy_%s", fluidEnthalpies_);
192 }
193
194private:
195 static bool rockInternalEnergyOutput_()
196 {
197 static bool val = Parameters::Get<Parameters::VtkWriteRockInternalEnergy>();
198 return val;
199 }
200
201 static bool totalThermalConductivityOutput_()
202 {
203 static bool val = Parameters::Get<Parameters::VtkWriteTotalThermalConductivity>();
204 return val;
205 }
206
207 static bool fluidInternalEnergiesOutput_()
208 {
209 static bool val = Parameters::Get<Parameters::VtkWriteFluidInternalEnergies>();
210 return val;
211 }
212
213 static bool fluidEnthalpiesOutput_()
214 {
215 static bool val = Parameters::Get<Parameters::VtkWriteFluidEnthalpies>();
216 return val;
217 }
218
219 ScalarBuffer rockInternalEnergy_;
220 ScalarBuffer totalThermalConductivity_;
221 PhaseBuffer fluidInternalEnergies_;
222 PhaseBuffer fluidEnthalpies_;
223};
224} // namespace Opm
225
226#endif
Declares the properties required by the black oil model.
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 the black oil model's energy related quantities.
Definition: vtkblackoilenergymodule.hh:62
void allocBuffers()
Allocate memory for the scalar fields we would like to write to the VTK file.
Definition: vtkblackoilenergymodule.hh:111
VtkBlackOilEnergyModule(const Simulator &simulator)
Definition: vtkblackoilenergymodule.hh:82
void commitBuffers(BaseOutputWriter &baseWriter)
Add all buffers to the VTK output writer.
Definition: vtkblackoilenergymodule.hh:172
static void registerParameters()
Register all run-time parameters for the multi-phase VTK output module.
Definition: vtkblackoilenergymodule.hh:90
void processElement(const ElementContext &elemCtx)
Modify the internal buffers according to the intensive quantities relevant for an element.
Definition: vtkblackoilenergymodule.hh:134
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: vtkblackoilenergymodule.hh:50
static constexpr bool value
Definition: vtkblackoilenergymodule.hh:50
Definition: vtkblackoilenergymodule.hh:49
static constexpr bool value
Definition: vtkblackoilenergymodule.hh:49
Definition: vtkblackoilenergymodule.hh:47
static constexpr bool value
Definition: vtkblackoilenergymodule.hh:47
Definition: vtkblackoilenergymodule.hh:48
static constexpr bool value
Definition: vtkblackoilenergymodule.hh:48