27#ifndef OPM_VTK_ENERGY_MODULE_HPP
28#define OPM_VTK_ENERGY_MODULE_HPP
30#include <opm/material/common/MathToolbox.hpp>
56template <
class TypeTag>
67 using ScalarBuffer =
typename ParentType::ScalarBuffer;
68 using PhaseBuffer =
typename ParentType::PhaseBuffer;
73 using Toolbox =
typename Opm::MathToolbox<Evaluation>;
97 if (params_.enthalpyOutput_) {
100 if (params_.internalEnergyOutput_) {
104 if (params_.solidInternalEnergyOutput_) {
107 if (params_.thermalConductivityOutput_) {
118 if (!Parameters::Get<Parameters::EnableVtkOutput>()) {
122 for (
unsigned i = 0; i < elemCtx.numPrimaryDof(0); ++i) {
123 unsigned I = elemCtx.globalSpaceIndex(i, 0);
124 const auto& intQuants = elemCtx.intensiveQuantities(i, 0);
125 const auto& fs = intQuants.fluidState();
127 if (params_.solidInternalEnergyOutput_) {
128 solidInternalEnergy_[I] = Toolbox::value(intQuants.solidInternalEnergy());
130 if (params_.thermalConductivityOutput_) {
131 thermalConductivity_[I] = Toolbox::value(intQuants.thermalConductivity());
134 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
135 if (params_.enthalpyOutput_) {
136 enthalpy_[phaseIdx][I] = Toolbox::value(fs.enthalpy(phaseIdx));
138 if (params_.internalEnergyOutput_) {
139 internalEnergy_[phaseIdx][I] = Toolbox::value(fs.internalEnergy(phaseIdx));
155 if (params_.solidInternalEnergyOutput_) {
158 if (params_.thermalConductivityOutput_) {
162 if (params_.enthalpyOutput_) {
165 if (params_.internalEnergyOutput_) {
172 PhaseBuffer enthalpy_{};
173 PhaseBuffer internalEnergy_{};
175 ScalarBuffer thermalConductivity_{};
176 ScalarBuffer solidInternalEnergy_{};
The base class for writer modules.
The base class for writer modules.
Definition baseoutputmodule.hh:67
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
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.hpp:58
void processElement(const ElementContext &elemCtx)
Modify the internal buffers according to the intensive quanties relevant for an element.
Definition vtkenergymodule.hpp:116
void allocBuffers()
Allocate memory for the scalar fields we would like to write to the VTK file.
Definition vtkenergymodule.hpp:95
void commitBuffers(BaseOutputWriter &baseWriter)
Add all buffers to the VTK output writer.
Definition vtkenergymodule.hpp:148
static void registerParameters()
Register all run-time parameters for the Vtk output module.
Definition vtkenergymodule.hpp:86
Simplifies writing multi-file VTK datasets.
Definition vtkmultiwriter.hh:66
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 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
constexpr auto getPropValue()
get the value data member of a property
Definition propertysystem.hh:242
This file provides the infrastructure to retrieve run-time parameters.
The Opm property system, traits with inheritance.
Struct holding the parameters for VtkEnergyModule.
Definition vtkenergyparams.hpp:46
static void registerParameters()
Registers the parameters in parameter system.
Definition vtkenergyparams.cpp:31
void read()
Reads the parameter values from the parameter system.
Definition vtkenergyparams.cpp:47
VTK output module for quantities which make sense for models which assume thermal equilibrium.
Simplifies writing multi-file VTK datasets.