27#ifndef EWOMS_VTK_ENERGY_MODULE_HH
28#define EWOMS_VTK_ENERGY_MODULE_HH
30#include <opm/material/common/MathToolbox.hpp>
64template <
class TypeTag>
78 static const int vtkFormat = getPropValue<TypeTag, Properties::VtkOutputFormat>();
79 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
81 using Toolbox =
typename Opm::MathToolbox<Evaluation>;
95 Parameters::Register<Parameters::VtkWriteSolidInternalEnergy>
96 (
"Include the volumetric internal energy of solid"
97 "matrix in the VTK output files");
98 Parameters::Register<Parameters::VtkWriteThermalConductivity>
99 (
"Include the total thermal conductivity of the"
100 "medium in the VTK output files");
101 Parameters::Register<Parameters::VtkWriteEnthalpies>
102 (
"Include the specific enthalpy of the phases in "
103 "the VTK output files");
104 Parameters::Register<Parameters::VtkWriteInternalEnergies>
105 (
"Include the specific internal energy of the "
106 "phases in the VTK output files");
115 if (enthalpyOutput_())
117 if (internalEnergyOutput_())
120 if (solidInternalEnergyOutput_())
122 if (thermalConductivityOutput_())
132 if (!Parameters::Get<Parameters::EnableVtkOutput>()) {
136 for (
unsigned i = 0; i < elemCtx.numPrimaryDof(0); ++i) {
137 unsigned I = elemCtx.globalSpaceIndex(i, 0);
138 const auto& intQuants = elemCtx.intensiveQuantities(i, 0);
139 const auto& fs = intQuants.fluidState();
141 if (solidInternalEnergyOutput_())
142 solidInternalEnergy_[I] = Toolbox::value(intQuants.solidInternalEnergy());
143 if (thermalConductivityOutput_())
144 thermalConductivity_[I] = Toolbox::value(intQuants.thermalConductivity());
146 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
147 if (enthalpyOutput_())
148 enthalpy_[phaseIdx][I] = Toolbox::value(fs.enthalpy(phaseIdx));
149 if (internalEnergyOutput_())
150 internalEnergy_[phaseIdx][I] = Toolbox::value(fs.internalEnergy(phaseIdx));
165 if (solidInternalEnergyOutput_())
167 if (thermalConductivityOutput_())
170 if (enthalpyOutput_())
172 if (internalEnergyOutput_())
177 static bool solidInternalEnergyOutput_()
179 static bool val = Parameters::Get<Parameters::VtkWriteSolidInternalEnergy>();
183 static bool thermalConductivityOutput_()
185 static bool val = Parameters::Get<Parameters::VtkWriteThermalConductivity>();
189 static bool enthalpyOutput_()
191 static bool val = Parameters::Get<Parameters::VtkWriteEnthalpies>();
195 static bool internalEnergyOutput_()
197 static bool val = Parameters::Get<Parameters::VtkWriteInternalEnergies>();
201 PhaseBuffer enthalpy_;
202 PhaseBuffer internalEnergy_;
204 ScalarBuffer thermalConductivity_;
205 ScalarBuffer solidInternalEnergy_;
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 quantities which make sense for models which assume thermal equilibrium.
Definition: vtkenergymodule.hh:66
void processElement(const ElementContext &elemCtx)
Modify the internal buffers according to the intensive quanties relevant for an element.
Definition: vtkenergymodule.hh:130
void allocBuffers()
Allocate memory for the scalar fields we would like to write to the VTK file.
Definition: vtkenergymodule.hh:113
void commitBuffers(BaseOutputWriter &baseWriter)
Add all buffers to the VTK output writer.
Definition: vtkenergymodule.hh:158
static void registerParameters()
Register all run-time parameters for the Vtk output module.
Definition: vtkenergymodule.hh:93
VtkEnergyModule(const Simulator &simulator)
Definition: vtkenergymodule.hh:85
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: vtkenergymodule.hh:46
static constexpr bool value
Definition: vtkenergymodule.hh:46
Definition: vtkenergymodule.hh:45
static constexpr bool value
Definition: vtkenergymodule.hh:45
Definition: vtkenergymodule.hh:43
static constexpr bool value
Definition: vtkenergymodule.hh:43
Definition: vtkenergymodule.hh:44
static constexpr bool value
Definition: vtkenergymodule.hh:44