vtkblackoilenergymodule.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_BLACK_OIL_ENERGY_MODULE_HPP
28#define OPM_VTK_BLACK_OIL_ENERGY_MODULE_HPP
29
30#include <dune/common/fvector.hh>
31
32#include <opm/material/densead/Math.hpp>
33
35
37
41
44
45namespace Opm {
46
52template <class TypeTag>
54{
56
63
64 static const int vtkFormat = getPropValue<TypeTag, Properties::VtkOutputFormat>();
66
67 enum { enableEnergy = getPropValue<TypeTag, Properties::EnableEnergy>() };
68 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
69
71 using PhaseBuffer = typename ParentType::PhaseBuffer;
72
73public:
74 VtkBlackOilEnergyModule(const Simulator& simulator)
75 : ParentType(simulator)
76 {
77 if constexpr (enableEnergy) {
78 params_.read();
79 }
80 }
81
86 static void registerParameters()
87 {
88 if constexpr (enableEnergy) {
90 }
91 }
92
98 {
99 if constexpr (enableEnergy) {
100 if (!Parameters::Get<Parameters::EnableVtkOutput>()) {
101 return;
102 }
103
104 if (params_.rockInternalEnergyOutput_) {
105 this->resizeScalarBuffer_(rockInternalEnergy_);
106 }
108 this->resizeScalarBuffer_(totalThermalConductivity_);
109 }
110 if (params_.fluidInternalEnergiesOutput_) {
111 this->resizePhaseBuffer_(fluidInternalEnergies_);
112 }
113 if (params_.fluidEnthalpiesOutput_) {
114 this->resizePhaseBuffer_(fluidEnthalpies_);
115 }
116 }
117 }
118
123 void processElement(const ElementContext& elemCtx)
124 {
125 if constexpr (enableEnergy) {
126 if (!Parameters::Get<Parameters::EnableVtkOutput>()) {
127 return;
128 }
129
130 for (unsigned dofIdx = 0; dofIdx < elemCtx.numPrimaryDof(/*timeIdx=*/0); ++dofIdx) {
131 const auto& intQuants = elemCtx.intensiveQuantities(dofIdx, /*timeIdx=*/0);
132 unsigned globalDofIdx = elemCtx.globalSpaceIndex(dofIdx, /*timeIdx=*/0);
133
134 if (params_.rockInternalEnergyOutput_) {
135 rockInternalEnergy_[globalDofIdx] =
136 scalarValue(intQuants.rockInternalEnergy());
137 }
138
140 totalThermalConductivity_[globalDofIdx] =
141 scalarValue(intQuants.totalThermalConductivity());
142 }
143
144 for (int phaseIdx = 0; phaseIdx < numPhases; ++ phaseIdx) {
145 if (FluidSystem::phaseIsActive(phaseIdx)) {
146 if (params_.fluidInternalEnergiesOutput_) {
147 fluidInternalEnergies_[phaseIdx][globalDofIdx] =
148 scalarValue(intQuants.fluidState().internalEnergy(phaseIdx));
149 }
150
151 if (params_.fluidEnthalpiesOutput_) {
152 fluidEnthalpies_[phaseIdx][globalDofIdx] =
153 scalarValue(intQuants.fluidState().enthalpy(phaseIdx));
154 }
155 }
156 }
157 }
158 }
159 }
160
165 {
166 if constexpr (enableEnergy) {
167 VtkMultiWriter* vtkWriter = dynamic_cast<VtkMultiWriter*>(&baseWriter);
168 if (!vtkWriter) {
169 return;
170 }
171
172 if (params_.rockInternalEnergyOutput_) {
173 this->commitScalarBuffer_(baseWriter, "volumetric internal energy rock", rockInternalEnergy_);
174 }
175
177 this->commitScalarBuffer_(baseWriter, "total thermal conductivity", totalThermalConductivity_);
178 }
179
180 if (params_.fluidInternalEnergiesOutput_) {
181 this->commitPhaseBuffer_(baseWriter, "internal energy_%s", fluidInternalEnergies_);
182 }
183
184 if (params_.fluidEnthalpiesOutput_) {
185 this->commitPhaseBuffer_(baseWriter, "enthalpy_%s", fluidEnthalpies_);
186 }
187 }
188 }
189
190private:
191 VtkBlackoilEnergyParams params_{};
192
193 ScalarBuffer rockInternalEnergy_{};
194 ScalarBuffer totalThermalConductivity_{};
195 PhaseBuffer fluidInternalEnergies_{};
196 PhaseBuffer fluidEnthalpies_{};
197};
198
199} // namespace Opm
200
201#endif // OPM_VTK_BLACKOIL_ENERGY_MODULE_HPP
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.hpp:54
void allocBuffers()
Allocate memory for the scalar fields we would like to write to the VTK file.
Definition: vtkblackoilenergymodule.hpp:97
VtkBlackOilEnergyModule(const Simulator &simulator)
Definition: vtkblackoilenergymodule.hpp:74
void commitBuffers(BaseOutputWriter &baseWriter)
Add all buffers to the VTK output writer.
Definition: vtkblackoilenergymodule.hpp:164
static void registerParameters()
Register all run-time parameters for the multi-phase VTK output module.
Definition: vtkblackoilenergymodule.hpp:86
void processElement(const ElementContext &elemCtx)
Modify the internal buffers according to the intensive quantities relevant for an element.
Definition: vtkblackoilenergymodule.hpp:123
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 VtkBlackoilEnergyOutputModule.
Definition: vtkblackoilenergyparams.hpp:46
bool fluidEnthalpiesOutput_
Definition: vtkblackoilenergyparams.hpp:56
bool rockInternalEnergyOutput_
Definition: vtkblackoilenergyparams.hpp:53
void read()
Reads the parameter values from the parameter system.
bool totalThermalConductivityOutput_
Definition: vtkblackoilenergyparams.hpp:54
static void registerParameters()
Registers the parameters in parameter system.
bool fluidInternalEnergiesOutput_
Definition: vtkblackoilenergyparams.hpp:55