vtkprimaryvarsmodule.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 EWOMS_VTK_PRIMARY_VARS_MODULE_HH
28#define EWOMS_VTK_PRIMARY_VARS_MODULE_HH
29
32
35
36namespace Opm::Properties {
37
38namespace TTag {
39
40// create new type tag for the VTK primary variables output
42
43} // namespace TTag
44
45// create the property tags needed for the primary variables module
46template<class TypeTag, class MyTypeTag>
48template<class TypeTag, class MyTypeTag>
50template<class TypeTag, class MyTypeTag>
52
53template<class TypeTag>
54struct VtkWritePrimaryVars<TypeTag, TTag::VtkPrimaryVars> { static constexpr bool value = false; };
55template<class TypeTag>
56struct VtkWriteProcessRank<TypeTag, TTag::VtkPrimaryVars> { static constexpr bool value = false; };
57template<class TypeTag>
58struct VtkWriteDofIndex<TypeTag, TTag::VtkPrimaryVars> { static constexpr bool value = false; };
59
60} // namespace Opm::Properties
61
62namespace Opm {
63
69template<class TypeTag>
71{
73
77
78 static const int vtkFormat = getPropValue<TypeTag, Properties::VtkOutputFormat>();
80
82 using EqBuffer = typename ParentType::EqBuffer;
83
84 enum { numEq = getPropValue<TypeTag, Properties::NumEq>() };
85
86public:
87 VtkPrimaryVarsModule(const Simulator& simulator)
88 : ParentType(simulator)
89 { }
90
94 static void registerParameters()
95 {
96 Parameters::registerParam<TypeTag, Properties::VtkWritePrimaryVars>
97 ("Include the primary variables into the VTK output files");
98 Parameters::registerParam<TypeTag, Properties::VtkWriteProcessRank>
99 ("Include the MPI process rank into the VTK output files");
100 Parameters::registerParam<TypeTag, Properties::VtkWriteDofIndex>
101 ("Include the index of the degrees of freedom into the VTK output files");
102 }
103
109 {
110 if (primaryVarsOutput_())
111 this->resizeEqBuffer_(primaryVars_);
112 if (processRankOutput_())
113 this->resizeScalarBuffer_(processRank_,
114 /*bufferType=*/ParentType::ElementBuffer);
115 if (dofIndexOutput_())
116 this->resizeScalarBuffer_(dofIndex_);
117 }
118
123 void processElement(const ElementContext& elemCtx)
124 {
125 if (!Parameters::get<TypeTag, Properties::EnableVtkOutput>())
126 return;
127
128 const auto& elementMapper = elemCtx.model().elementMapper();
129 unsigned elemIdx = static_cast<unsigned>(elementMapper.index(elemCtx.element()));
130 if (processRankOutput_() && !processRank_.empty())
131 processRank_[elemIdx] = static_cast<unsigned>(this->simulator_.gridView().comm().rank());
132
133 for (unsigned i = 0; i < elemCtx.numPrimaryDof(/*timeIdx=*/0); ++i) {
134 unsigned I = elemCtx.globalSpaceIndex(i, /*timeIdx=*/0);
135 const auto& priVars = elemCtx.primaryVars(i, /*timeIdx=*/0);
136
137 if (dofIndexOutput_())
138 dofIndex_[I] = I;
139
140 for (unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
141 if (primaryVarsOutput_() && !primaryVars_[eqIdx].empty())
142 primaryVars_[eqIdx][I] = priVars[eqIdx];
143 }
144 }
145 }
146
151 {
152 VtkMultiWriter *vtkWriter = dynamic_cast<VtkMultiWriter*>(&baseWriter);
153 if (!vtkWriter) {
154 return;
155 }
156
157 if (primaryVarsOutput_())
158 this->commitPriVarsBuffer_(baseWriter, "PV_%s", primaryVars_);
159 if (processRankOutput_())
160 this->commitScalarBuffer_(baseWriter,
161 "process rank",
162 processRank_,
163 /*bufferType=*/ParentType::ElementBuffer);
164 if (dofIndexOutput_())
165 this->commitScalarBuffer_(baseWriter, "DOF index", dofIndex_);
166 }
167
168private:
169 static bool primaryVarsOutput_()
170 {
171 static bool val = Parameters::get<TypeTag, Properties::VtkWritePrimaryVars>();
172 return val;
173 }
174 static bool processRankOutput_()
175 {
176 static bool val = Parameters::get<TypeTag, Properties::VtkWriteProcessRank>();
177 return val;
178 }
179 static bool dofIndexOutput_()
180 {
181 static bool val = Parameters::get<TypeTag, Properties::VtkWriteDofIndex>();
182 return val;
183 }
184
185 EqBuffer primaryVars_;
186 ScalarBuffer processRank_;
187 ScalarBuffer dofIndex_;
188};
189
190} // namespace Opm
191
192#endif
The base class for writer modules.
Definition: baseoutputmodule.hh:73
BaseOutputWriter::ScalarBuffer ScalarBuffer
Definition: baseoutputmodule.hh:91
const Simulator & simulator_
Definition: baseoutputmodule.hh:519
@ ElementBuffer
Buffer contains data associated with the grid's elements.
Definition: baseoutputmodule.hh:156
void commitPriVarsBuffer_(BaseOutputWriter &baseWriter, const char *pattern, EqBuffer &buffer, BufferType bufferType=DofBuffer)
Add a buffer with as many variables as PDEs to the result file.
Definition: baseoutputmodule.hh:370
void resizeEqBuffer_(EqBuffer &buffer, BufferType bufferType=DofBuffer)
Allocate the space for a buffer storing a equation specific quantity.
Definition: baseoutputmodule.hh:223
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
std::array< ScalarBuffer, numEq > EqBuffer
Definition: baseoutputmodule.hh:95
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 fluid composition.
Definition: vtkprimaryvarsmodule.hh:71
VtkPrimaryVarsModule(const Simulator &simulator)
Definition: vtkprimaryvarsmodule.hh:87
void processElement(const ElementContext &elemCtx)
Modify the internal buffers according to the intensive quantities relevant for an element.
Definition: vtkprimaryvarsmodule.hh:123
void commitBuffers(BaseOutputWriter &baseWriter)
Add all buffers to the VTK output writer.
Definition: vtkprimaryvarsmodule.hh:150
void allocBuffers()
Allocate memory for the scalar fields we would like to write to the VTK file.
Definition: vtkprimaryvarsmodule.hh:108
static void registerParameters()
Register all run-time parameters for the Vtk output module.
Definition: vtkprimaryvarsmodule.hh:94
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: vtkprimaryvarsmodule.hh:41
a tag to mark properties as undefined
Definition: propertysystem.hh:40
Definition: vtkprimaryvarsmodule.hh:51
Definition: vtkprimaryvarsmodule.hh:47
Definition: vtkprimaryvarsmodule.hh:49