27#ifndef OPM_VTK_BLACK_OIL_ENERGY_MODULE_HPP
28#define OPM_VTK_BLACK_OIL_ENERGY_MODULE_HPP
30#include <dune/common/fvector.hh>
32#include <opm/material/densead/Math.hpp>
52template <
class TypeTag>
53class VtkBlackOilEnergyModule :
public BaseOutputModule<TypeTag>
55 using ParentType = BaseOutputModule<TypeTag>;
71 using PhaseBuffer =
typename ParentType::PhaseBuffer;
72 using ScalarBuffer =
typename ParentType::ScalarBuffer;
75 explicit VtkBlackOilEnergyModule(
const Simulator& simulator)
76 : ParentType(simulator)
78 if constexpr (enableEnergy) {
89 if constexpr (enableEnergy) {
100 if constexpr (enableEnergy) {
105 if (params_.rockInternalEnergyOutput_) {
108 if (params_.totalThermalConductivityOutput_) {
111 if (params_.fluidInternalEnergiesOutput_) {
114 if (params_.fluidEnthalpiesOutput_) {
126 if constexpr (enableEnergy) {
131 for (
unsigned dofIdx = 0; dofIdx < elemCtx.numPrimaryDof(0); ++dofIdx) {
132 const auto& intQuants = elemCtx.intensiveQuantities(dofIdx, 0);
133 const unsigned globalDofIdx = elemCtx.globalSpaceIndex(dofIdx, 0);
135 if (params_.rockInternalEnergyOutput_) {
136 rockInternalEnergy_[globalDofIdx] =
137 scalarValue(intQuants.rockInternalEnergy());
140 if (params_.totalThermalConductivityOutput_) {
141 totalThermalConductivity_[globalDofIdx] =
142 scalarValue(intQuants.totalThermalConductivity());
145 for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
146 if (FluidSystem::phaseIsActive(phaseIdx)) {
147 if (params_.fluidInternalEnergiesOutput_) {
148 fluidInternalEnergies_[phaseIdx][globalDofIdx] =
149 scalarValue(intQuants.fluidState().internalEnergy(phaseIdx));
152 if (params_.fluidEnthalpiesOutput_) {
153 fluidEnthalpies_[phaseIdx][globalDofIdx] =
154 scalarValue(intQuants.fluidState().enthalpy(phaseIdx));
167 if constexpr (enableEnergy) {
168 if (!
dynamic_cast<VtkMultiWriter*
>(&baseWriter)) {
172 if (params_.rockInternalEnergyOutput_) {
174 rockInternalEnergy_, BufferType::Dof);
177 if (params_.totalThermalConductivityOutput_) {
179 totalThermalConductivity_, BufferType::Dof);
182 if (params_.fluidInternalEnergiesOutput_) {
184 fluidInternalEnergies_, BufferType::Dof);
187 if (params_.fluidEnthalpiesOutput_) {
189 fluidEnthalpies_, BufferType::Dof);
197 ScalarBuffer rockInternalEnergy_{};
198 ScalarBuffer totalThermalConductivity_{};
199 PhaseBuffer fluidInternalEnergies_{};
200 PhaseBuffer fluidEnthalpies_{};
The base class for writer modules.
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:124
void allocBuffers() override
Allocate memory for the scalar fields we would like to write to the VTK file.
Definition vtkblackoilenergymodule.hpp:98
void commitBuffers(BaseOutputWriter &baseWriter) override
Add all buffers to the VTK output writer.
Definition vtkblackoilenergymodule.hpp:165
static void registerParameters()
Register all run-time parameters for the multi-phase VTK output module.
Definition vtkblackoilenergymodule.hpp:87
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:43
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.