vtkptflashmodule.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3/*
4 This file is part of the Open Porous Media project (OPM).
5
6 OPM is free software: you can redistribute it and/or modify
7 it under the terms of the GNU General Public License as published by
8 the Free Software Foundation, either version 2 of the License, or
9 (at your option) any later version.
10
11 OPM is distributed in the hope that it will be useful,
12 but WITHOUT ANY WARRANTY; without even the implied warranty of
13 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14 GNU General Public License for more details.
15
16 You should have received a copy of the GNU General Public License
17 along with OPM. If not, see <http://www.gnu.org/licenses/>.
18
19 Consult the COPYING file in the top-level source directory of this
20 module for the precise wording of the license and the list of
21 copyright holders.
22*/
27#ifndef OPM_VTK_PTFLASH_MODULE_HH
28#define OPM_VTK_PTFLASH_MODULE_HH
29
30#include "vtkmultiwriter.hh"
31#include "baseoutputmodule.hh"
32
35
36#include <opm/material/common/MathToolbox.hpp>
37
38namespace Opm::Properties {
39
40namespace TTag {
41
42// create new type tag for the VTK PTFlash output
43struct VtkPTFlash {};
44
45} // namespace TTag
46
47// create the property tags needed for the composition module
48template<class TypeTag, class MyTypeTag>
50template<class TypeTag, class MyTypeTag>
52
53// set default values for what quantities to output
54template<class TypeTag>
55struct VtkWriteLiquidMoleFractions<TypeTag, TTag::VtkPTFlash> { static constexpr bool value = false; };
56template<class TypeTag>
57struct VtkWriteEquilibriumConstants<TypeTag, TTag::VtkPTFlash> { static constexpr bool value = false; };
58
59} // namespace Opm::Properties
60
61namespace Opm {
62
71template <class TypeTag>
72class VtkPTFlashModule: public BaseOutputModule<TypeTag>
73{
75
80
82
83 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
84 enum { numComponents = getPropValue<TypeTag, Properties::NumComponents>() };
85
86 static const int vtkFormat = getPropValue<TypeTag, Properties::VtkOutputFormat>();
88
91
92public:
93 explicit VtkPTFlashModule(const Simulator& simulator)
94 : ParentType(simulator)
95 { }
96
100 static void registerParameters()
101 {
102 Parameters::registerParam<TypeTag, Properties::VtkWriteLiquidMoleFractions>
103 ("Include liquid mole fractions (L) in the VTK output files");
104 Parameters::registerParam<TypeTag, Properties::VtkWriteEquilibriumConstants>
105 ("Include equilibrium constants (K) in the VTK output files");
106 }
107
113 {
114 if (LOutput_())
115 this->resizeScalarBuffer_(L_);
116 if (equilConstOutput_())
117 this->resizeComponentBuffer_(K_);
118 }
119
124 void processElement(const ElementContext& elemCtx)
125 {
126 using Toolbox = MathToolbox<Evaluation>;
127
128 if (!Parameters::get<TypeTag, Properties::EnableVtkOutput>())
129 return;
130
131 for (unsigned i = 0; i < elemCtx.numPrimaryDof(/*timeIdx=*/0); ++i) {
132 unsigned I = elemCtx.globalSpaceIndex(i, /*timeIdx=*/0);
133 const auto& intQuants = elemCtx.intensiveQuantities(i, /*timeIdx=*/0);
134 const auto& fs = intQuants.fluidState();
135
136 if (LOutput_())
137 L_[I] = Toolbox::value(fs.L());
138
139 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
140 if (equilConstOutput_())
141 K_[compIdx][I] = Toolbox::value(fs.K(compIdx));
142 }
143 }
144 }
145
150 {
151 auto *vtkWriter = dynamic_cast<VtkMultiWriter*>(&baseWriter);
152 if (!vtkWriter) {
153 return;
154 }
155
156 if (equilConstOutput_())
157 this->commitComponentBuffer_(baseWriter, "K^%s", K_);
158 if (LOutput_())
159 this->commitScalarBuffer_(baseWriter, "L", L_);
160 }
161
162private:
163 static bool LOutput_()
164 {
165 static bool val = Parameters::get<TypeTag, Properties::VtkWriteLiquidMoleFractions>();
166 return val;
167 }
168
169 static bool equilConstOutput_()
170 {
171 static bool val = Parameters::get<TypeTag, Properties::VtkWriteEquilibriumConstants>();
172 return val;
173 }
174
175 ComponentBuffer K_;
176 ScalarBuffer L_;
177};
178
179} // namespace Opm
180
181#endif
The base class for writer modules.
Definition: baseoutputmodule.hh:73
BaseOutputWriter::ScalarBuffer ScalarBuffer
Definition: baseoutputmodule.hh:91
void commitComponentBuffer_(BaseOutputWriter &baseWriter, const char *pattern, ComponentBuffer &buffer, BufferType bufferType=DofBuffer)
Add a component-specific buffer to the result file.
Definition: baseoutputmodule.hh:442
std::array< ScalarBuffer, numComponents > ComponentBuffer
Definition: baseoutputmodule.hh:97
void resizeComponentBuffer_(ComponentBuffer &buffer, BufferType bufferType=DofBuffer)
Allocate the space for a buffer storing a component specific quantity.
Definition: baseoutputmodule.hh:269
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:316
void resizeScalarBuffer_(ScalarBuffer &buffer, BufferType bufferType=DofBuffer)
Allocate the space for a buffer storing a scalar quantity.
Definition: baseoutputmodule.hh:162
The base class for all output writers.
Definition: baseoutputwriter.hh:44
Simplifies writing multi-file VTK datasets.
Definition: vtkmultiwriter.hh:66
VTK output module for the PT Flash calculation This module deals with the following quantities: K,...
Definition: vtkptflashmodule.hh:73
static void registerParameters()
Register all run-time parameters for the Vtk output module.
Definition: vtkptflashmodule.hh:100
void commitBuffers(BaseOutputWriter &baseWriter)
Add all buffers to the VTK output writer.
Definition: vtkptflashmodule.hh:149
void allocBuffers()
Allocate memory for the scalar fields we would like to write to the VTK file.
Definition: vtkptflashmodule.hh:112
void processElement(const ElementContext &elemCtx)
Modify the internal buffers according to the intensive quantities relevant for an element.
Definition: vtkptflashmodule.hh:124
VtkPTFlashModule(const Simulator &simulator)
Definition: vtkptflashmodule.hh:93
Definition: blackoilmodel.hh:72
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:242
This file provides the infrastructure to retrieve run-time parameters.
The Opm property system, traits with inheritance.
Definition: vtkptflashmodule.hh:43
a tag to mark properties as undefined
Definition: propertysystem.hh:40
Definition: vtkptflashmodule.hh:51
Definition: vtkptflashmodule.hh:49