27#ifndef EWOMS_VTK_BLACK_OIL_ENERGY_MODULE_HH
28#define EWOMS_VTK_BLACK_OIL_ENERGY_MODULE_HH
30#include <dune/common/fvector.hh>
32#include <opm/material/densead/Math.hpp>
60template <
class TypeTag>
72 static const int vtkFormat = getPropValue<TypeTag, Properties::VtkOutputFormat>();
75 enum { enableEnergy = getPropValue<TypeTag, Properties::EnableEnergy>() };
76 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
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");
113 if (!Parameters::Get<Parameters::EnableVtkOutput>()) {
120 if (rockInternalEnergyOutput_())
122 if (totalThermalConductivityOutput_())
124 if (fluidInternalEnergiesOutput_())
126 if (fluidEnthalpiesOutput_())
136 if (!Parameters::Get<Parameters::EnableVtkOutput>()) {
143 for (
unsigned dofIdx = 0; dofIdx < elemCtx.numPrimaryDof(0); ++dofIdx) {
144 const auto& intQuants = elemCtx.intensiveQuantities(dofIdx, 0);
145 unsigned globalDofIdx = elemCtx.globalSpaceIndex(dofIdx, 0);
147 if (rockInternalEnergyOutput_())
148 rockInternalEnergy_[globalDofIdx] =
149 scalarValue(intQuants.rockInternalEnergy());
151 if (totalThermalConductivityOutput_())
152 totalThermalConductivity_[globalDofIdx] =
153 scalarValue(intQuants.totalThermalConductivity());
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));
161 if (fluidEnthalpiesOutput_())
162 fluidEnthalpies_[phaseIdx][globalDofIdx] =
163 scalarValue(intQuants.fluidState().enthalpy(phaseIdx));
181 if (rockInternalEnergyOutput_())
182 this->
commitScalarBuffer_(baseWriter,
"volumetric internal energy rock", rockInternalEnergy_);
184 if (totalThermalConductivityOutput_())
185 this->
commitScalarBuffer_(baseWriter,
"total thermal conductivity", totalThermalConductivity_);
187 if (fluidInternalEnergiesOutput_())
190 if (fluidEnthalpiesOutput_())
195 static bool rockInternalEnergyOutput_()
197 static bool val = Parameters::Get<Parameters::VtkWriteRockInternalEnergy>();
201 static bool totalThermalConductivityOutput_()
203 static bool val = Parameters::Get<Parameters::VtkWriteTotalThermalConductivity>();
207 static bool fluidInternalEnergiesOutput_()
209 static bool val = Parameters::Get<Parameters::VtkWriteFluidInternalEnergies>();
213 static bool fluidEnthalpiesOutput_()
215 static bool val = Parameters::Get<Parameters::VtkWriteFluidEnthalpies>();
219 ScalarBuffer rockInternalEnergy_;
220 ScalarBuffer totalThermalConductivity_;
221 PhaseBuffer fluidInternalEnergies_;
222 PhaseBuffer fluidEnthalpies_;
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