27#ifndef OPM_VTK_ENERGY_MODULE_HPP
28#define OPM_VTK_ENERGY_MODULE_HPP
30#include <opm/material/common/MathToolbox.hpp>
56template <
class TypeTag>
57class VtkEnergyModule :
public BaseOutputModule<TypeTag>
59 using ParentType = BaseOutputModule<TypeTag>;
68 using ScalarBuffer =
typename ParentType::ScalarBuffer;
69 using PhaseBuffer =
typename ParentType::PhaseBuffer;
74 using Toolbox =
typename Opm::MathToolbox<Evaluation>;
78 explicit VtkEnergyModule(
const Simulator& simulator)
79 : ParentType(simulator)
98 if (params_.enthalpyOutput_) {
101 if (params_.internalEnergyOutput_) {
105 if (params_.solidInternalEnergyOutput_) {
108 if (params_.thermalConductivityOutput_) {
123 for (
unsigned i = 0; i < elemCtx.numPrimaryDof(0); ++i) {
124 const unsigned I = elemCtx.globalSpaceIndex(i, 0);
125 const auto& intQuants = elemCtx.intensiveQuantities(i, 0);
126 const auto& fs = intQuants.fluidState();
128 if (params_.solidInternalEnergyOutput_) {
129 solidInternalEnergy_[I] = Toolbox::value(intQuants.solidInternalEnergy());
131 if (params_.thermalConductivityOutput_) {
132 thermalConductivity_[I] = Toolbox::value(intQuants.thermalConductivity());
135 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
136 if (params_.enthalpyOutput_) {
137 enthalpy_[phaseIdx][I] = Toolbox::value(fs.enthalpy(phaseIdx));
139 if (params_.internalEnergyOutput_) {
140 internalEnergy_[phaseIdx][I] = Toolbox::value(fs.internalEnergy(phaseIdx));
151 if (!
dynamic_cast<VtkMultiWriter*
>(&baseWriter)) {
155 if (params_.solidInternalEnergyOutput_) {
157 solidInternalEnergy_, BufferType::Dof);
159 if (params_.thermalConductivityOutput_) {
161 thermalConductivity_, BufferType::Dof);
164 if (params_.enthalpyOutput_) {
166 enthalpy_, BufferType::Dof);
168 if (params_.internalEnergyOutput_) {
170 internalEnergy_, BufferType::Dof);
176 PhaseBuffer enthalpy_{};
177 PhaseBuffer internalEnergy_{};
179 ScalarBuffer thermalConductivity_{};
180 ScalarBuffer solidInternalEnergy_{};
The base class for writer modules.
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 allocBuffers() override
Allocate memory for the scalar fields we would like to write to the VTK file.
Definition vtkenergymodule.hpp:96
static void registerParameters()
Register all run-time parameters for the Vtk output module.
Definition vtkenergymodule.hpp:87
void processElement(const ElementContext &elemCtx) override
Modify the internal buffers according to the intensive quanties relevant for an element.
Definition vtkenergymodule.hpp:117
void commitBuffers(BaseOutputWriter &baseWriter) override
Add all buffers to the VTK output writer.
Definition vtkenergymodule.hpp:149
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 VtkEnergyModule.
Definition vtkenergyparams.hpp:46
static void registerParameters()
Registers the parameters in parameter system.
Definition vtkenergyparams.cpp:31
VTK output module for quantities which make sense for models which assume thermal equilibrium.
Simplifies writing multi-file VTK datasets.