vtkblackoilpolymermodule.hpp
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_BLACK_OIL_POLYMER_MODULE_HPP
28#define OPM_VTK_BLACK_OIL_POLYMER_MODULE_HPP
29
30#include <dune/common/fvector.hh>
31
32#include <opm/material/densead/Math.hpp>
33
35
37
41
44
45namespace Opm::Properties::TTag {
46
47// create new type tag for the VTK multi-phase output
49
50} // namespace Opm::Properties::TTag
51
52namespace Opm {
53
59template <class TypeTag>
61{
63
69
70 static const int vtkFormat = getPropValue<TypeTag, Properties::VtkOutputFormat>();
72
73 enum { enablePolymer = getPropValue<TypeTag, Properties::EnablePolymer>() };
74
76
77public:
78 VtkBlackOilPolymerModule(const Simulator& simulator)
79 : ParentType(simulator)
80 {
81 if constexpr (enablePolymer) {
82 params_.read();
83 }
84 }
85
90 static void registerParameters()
91 {
92 if constexpr (enablePolymer) {
94 }
95 }
96
102 {
103 if constexpr (enablePolymer) {
104 if (!Parameters::Get<Parameters::EnableVtkOutput>()) {
105 return;
106 }
107
108 if (params_.polymerConcentrationOutput_) {
109 this->resizeScalarBuffer_(polymerConcentration_);
110 }
111 if (params_.polymerDeadPoreVolumeOutput_) {
112 this->resizeScalarBuffer_(polymerDeadPoreVolume_);
113 }
114 if (params_.polymerRockDensityOutput_) {
115 this->resizeScalarBuffer_(polymerRockDensity_);
116 }
117 if (params_.polymerAdsorptionOutput_) {
118 this->resizeScalarBuffer_(polymerAdsorption_);
119 }
121 this->resizeScalarBuffer_(polymerViscosityCorrection_);
122 }
124 this->resizeScalarBuffer_(waterViscosityCorrection_);
125 }
126 }
127 }
128
133 void processElement(const ElementContext& elemCtx)
134 {
135 if constexpr (enablePolymer) {
136 if (!Parameters::Get<Parameters::EnableVtkOutput>()) {
137 return;
138 }
139
140 for (unsigned dofIdx = 0; dofIdx < elemCtx.numPrimaryDof(/*timeIdx=*/0); ++dofIdx) {
141 const auto& intQuants = elemCtx.intensiveQuantities(dofIdx, /*timeIdx=*/0);
142 unsigned globalDofIdx = elemCtx.globalSpaceIndex(dofIdx, /*timeIdx=*/0);
143
144 if (params_.polymerConcentrationOutput_) {
145 polymerConcentration_[globalDofIdx] =
146 scalarValue(intQuants.polymerConcentration());
147 }
148
149 if (params_.polymerDeadPoreVolumeOutput_) {
150 polymerDeadPoreVolume_[globalDofIdx] =
151 scalarValue(intQuants.polymerDeadPoreVolume());
152 }
153
154 if (params_.polymerRockDensityOutput_) {
155 polymerRockDensity_[globalDofIdx] =
156 scalarValue(intQuants.polymerRockDensity());
157 }
158
159 if (params_.polymerAdsorptionOutput_) {
160 polymerAdsorption_[globalDofIdx] =
161 scalarValue(intQuants.polymerAdsorption());
162 }
163
165 polymerViscosityCorrection_[globalDofIdx] =
166 scalarValue(intQuants.polymerViscosityCorrection());
167 }
168
170 waterViscosityCorrection_[globalDofIdx] =
171 scalarValue(intQuants.waterViscosityCorrection());
172 }
173 }
174 }
175 }
176
181 {
182 if constexpr (enablePolymer) {
183 VtkMultiWriter* vtkWriter = dynamic_cast<VtkMultiWriter*>(&baseWriter);
184 if (!vtkWriter) {
185 return;
186 }
187
188 if (params_.polymerConcentrationOutput_) {
189 this->commitScalarBuffer_(baseWriter, "polymer concentration", polymerConcentration_);
190 }
191
192 if (params_.polymerDeadPoreVolumeOutput_) {
193 this->commitScalarBuffer_(baseWriter, "dead pore volume fraction", polymerDeadPoreVolume_);
194 }
195
196 if (params_.polymerRockDensityOutput_) {
197 this->commitScalarBuffer_(baseWriter, "polymer rock density", polymerRockDensity_);
198 }
199
200 if (params_.polymerAdsorptionOutput_) {
201 this->commitScalarBuffer_(baseWriter, "polymer adsorption", polymerAdsorption_);
202 }
203
205 this->commitScalarBuffer_(baseWriter, "polymer viscosity correction", polymerViscosityCorrection_);
206 }
207
209 this->commitScalarBuffer_(baseWriter, "water viscosity correction", waterViscosityCorrection_);
210 }
211 }
212 }
213
214private:
215 VtkBlackoilPolymerParams params_{};
216 ScalarBuffer polymerConcentration_{};
217 ScalarBuffer polymerDeadPoreVolume_{};
218 ScalarBuffer polymerRockDensity_{};
219 ScalarBuffer polymerAdsorption_{};
220 ScalarBuffer polymerViscosityCorrection_{};
221 ScalarBuffer waterViscosityCorrection_{};
222};
223
224} // namespace Opm
225
226#endif // OPM_VTK_BLACK_OIL_POLYMER_MODULE_HPP
Declares the properties required by the black oil model.
The base class for writer modules.
Definition: baseoutputmodule.hh:67
BaseOutputWriter::ScalarBuffer ScalarBuffer
Definition: baseoutputmodule.hh:85
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 the black oil model's polymer related quantities.
Definition: vtkblackoilpolymermodule.hpp:61
void allocBuffers()
Allocate memory for the scalar fields we would like to write to the VTK file.
Definition: vtkblackoilpolymermodule.hpp:101
static void registerParameters()
Register all run-time parameters for the multi-phase VTK output module.
Definition: vtkblackoilpolymermodule.hpp:90
void processElement(const ElementContext &elemCtx)
Modify the internal buffers according to the intensive quantities relevant for an element.
Definition: vtkblackoilpolymermodule.hpp:133
VtkBlackOilPolymerModule(const Simulator &simulator)
Definition: vtkblackoilpolymermodule.hpp:78
void commitBuffers(BaseOutputWriter &baseWriter)
Add all buffers to the VTK output writer.
Definition: vtkblackoilpolymermodule.hpp:180
Simplifies writing multi-file VTK datasets.
Definition: vtkmultiwriter.hh:66
Declare the properties used by the infrastructure code of the finite volume discretizations.
The generic type tag for problems using the immiscible multi-phase model.
Definition: blackoilmodel.hh:74
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
This file provides the infrastructure to retrieve run-time parameters.
The Opm property system, traits with inheritance.
Definition: vtkblackoilpolymermodule.hpp:48
Struct holding the parameters for VtkBlackoilPolymerModule.
Definition: vtkblackoilpolymerparams.hpp:48
bool polymerRockDensityOutput_
Definition: vtkblackoilpolymerparams.hpp:57
bool polymerViscosityCorrectionOutput_
Definition: vtkblackoilpolymerparams.hpp:59
bool polymerDeadPoreVolumeOutput_
Definition: vtkblackoilpolymerparams.hpp:56
void read()
Reads the parameter values from the parameter system.
static void registerParameters()
Registers the parameters in parameter system.
bool waterViscosityCorrectionOutput_
Definition: vtkblackoilpolymerparams.hpp:60
bool polymerAdsorptionOutput_
Definition: vtkblackoilpolymerparams.hpp:58
bool polymerConcentrationOutput_
Definition: vtkblackoilpolymerparams.hpp:55