opm-simulators
Loading...
Searching...
No Matches
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
37
41
44
45namespace Opm {
46
52template <class TypeTag>
53class VtkBlackOilEnergyModule : public BaseOutputModule<TypeTag>
54{
55 using ParentType = BaseOutputModule<TypeTag>;
56
63
64 static constexpr auto vtkFormat = getPropValue<TypeTag, Properties::VtkOutputFormat>();
65 using VtkMultiWriter = ::Opm::VtkMultiWriter<GridView, vtkFormat>;
66
68
69 using BufferType = typename ParentType::BufferType;
70 using PhaseBuffer = typename ParentType::PhaseBuffer;
71 using ScalarBuffer = typename ParentType::ScalarBuffer;
72
73public:
74 explicit VtkBlackOilEnergyModule(const Simulator& simulator)
75 : ParentType(simulator)
76 {
77 params_.read();
78 }
79
88
93 void allocBuffers() override
94 {
96 return;
97 }
98
99 if (params_.rockInternalEnergyOutput_) {
100 this->resizeScalarBuffer_(rockInternalEnergy_, BufferType::Dof);
101 }
102 if (params_.totalThermalConductivityOutput_) {
103 this->resizeScalarBuffer_(totalThermalConductivity_, BufferType::Dof);
104 }
105 if (params_.fluidInternalEnergiesOutput_) {
106 this->resizePhaseBuffer_(fluidInternalEnergies_, BufferType::Dof);
107 }
108 if (params_.fluidEnthalpiesOutput_) {
109 this->resizePhaseBuffer_(fluidEnthalpies_, BufferType::Dof);
110 }
111 }
112
117 void processElement(const ElementContext& elemCtx) override
118 {
120 return;
121 }
122
123 for (unsigned dofIdx = 0; dofIdx < elemCtx.numPrimaryDof(/*timeIdx=*/0); ++dofIdx) {
124 const auto& intQuants = elemCtx.intensiveQuantities(dofIdx, /*timeIdx=*/0);
125 const unsigned globalDofIdx = elemCtx.globalSpaceIndex(dofIdx, /*timeIdx=*/0);
126
127 if (params_.rockInternalEnergyOutput_) {
128 rockInternalEnergy_[globalDofIdx] =
129 scalarValue(intQuants.rockInternalEnergy());
130 }
131
132 if (params_.totalThermalConductivityOutput_) {
133 totalThermalConductivity_[globalDofIdx] =
134 scalarValue(intQuants.totalThermalConductivity());
135 }
136
137 for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
138 if (FluidSystem::phaseIsActive(phaseIdx)) {
139 if (params_.fluidInternalEnergiesOutput_) {
140 fluidInternalEnergies_[phaseIdx][globalDofIdx] =
141 scalarValue(intQuants.fluidState().internalEnergy(phaseIdx));
142 }
143
144 if (params_.fluidEnthalpiesOutput_) {
145 fluidEnthalpies_[phaseIdx][globalDofIdx] =
146 scalarValue(intQuants.fluidState().enthalpy(phaseIdx));
147 }
148 }
149 }
150 }
151 }
152
156 void commitBuffers(BaseOutputWriter& baseWriter) override
157 {
158 if (!dynamic_cast<VtkMultiWriter*>(&baseWriter)) {
159 return;
160 }
161
162 if (params_.rockInternalEnergyOutput_) {
163 this->commitScalarBuffer_(baseWriter, "volumetric internal energy rock",
164 rockInternalEnergy_, BufferType::Dof);
165 }
166
167 if (params_.totalThermalConductivityOutput_) {
168 this->commitScalarBuffer_(baseWriter, "total thermal conductivity",
169 totalThermalConductivity_, BufferType::Dof);
170 }
171
172 if (params_.fluidInternalEnergiesOutput_) {
173 this->commitPhaseBuffer_(baseWriter, "internal energy_%s",
174 fluidInternalEnergies_, BufferType::Dof);
175 }
176
177 if (params_.fluidEnthalpiesOutput_) {
178 this->commitPhaseBuffer_(baseWriter, "enthalpy_%s",
179 fluidEnthalpies_, BufferType::Dof);
180 }
181 }
182
183private:
184 VtkBlackoilEnergyParams params_{};
185
186 ScalarBuffer rockInternalEnergy_{};
187 ScalarBuffer totalThermalConductivity_{};
188 PhaseBuffer fluidInternalEnergies_{};
189 PhaseBuffer fluidEnthalpies_{};
190};
191
192} // namespace Opm
193
194#endif // OPM_VTK_BLACKOIL_ENERGY_MODULE_HPP
The base class for writer modules.
Contains the classes required to extend the black-oil model by energy.
Declares the properties required by the black oil model.
void resizeScalarBuffer_(ScalarBuffer &buffer, BufferType bufferType)
Allocate the space for a buffer storing a scalar quantity.
Definition baseoutputmodule.hh:157
BufferType
Definition baseoutputmodule.hh:143
void commitPhaseBuffer_(BaseOutputWriter &baseWriter, const char *pattern, PhaseBuffer &buffer, BufferType bufferType)
Add a phase-specific buffer to the result file.
Definition baseoutputmodule.hh:337
void commitScalarBuffer_(BaseOutputWriter &baseWriter, const char *name, ScalarBuffer &buffer, BufferType bufferType)
Add a buffer containing scalar quantities to the result file.
Definition baseoutputmodule.hh:238
void resizePhaseBuffer_(PhaseBuffer &buffer, BufferType bufferType)
Allocate the space for a buffer storing a phase-specific quantity.
Definition baseoutputmodule.hh:198
The base class for all output writers.
Definition baseoutputwriter.hh:46
void processElement(const ElementContext &elemCtx) override
Modify the internal buffers according to the intensive quantities relevant for an element.
Definition vtkblackoilenergymodule.hpp:117
void allocBuffers() override
Allocate memory for the scalar fields we would like to write to the VTK file.
Definition vtkblackoilenergymodule.hpp:93
void commitBuffers(BaseOutputWriter &baseWriter) override
Add all buffers to the VTK output writer.
Definition vtkblackoilenergymodule.hpp:156
static void registerParameters()
Register all run-time parameters for the multi-phase VTK output module.
Definition vtkblackoilenergymodule.hpp:84
Simplifies writing multi-file VTK datasets.
Definition vtkmultiwriter.hh:65
Declare the properties used by the infrastructure code of the finite volume discretizations.
This file contains a set of helper functions used by VFPProd / VFPInj.
Definition blackoilbioeffectsmodules.hh:45
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:233
constexpr auto getPropValue()
get the value data member of a property
Definition propertysystem.hh:240
This file provides the infrastructure to retrieve run-time parameters.
auto Get(bool errorIfNotRegistered=true)
Retrieve a runtime parameter.
Definition parametersystem.hpp:187
The Opm property system, traits with inheritance.
Struct holding the parameters for VtkBlackoilEnergyOutputModule.
Definition vtkblackoilenergyparams.hpp:46
static void registerParameters()
Registers the parameters in parameter system.
Definition vtkblackoilenergyparams.cpp:31
VTK output module for the black oil model's energy related quantities.
Simplifies writing multi-file VTK datasets.