vtkblackoilmicpmodule.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_BLACK_OIL_MICP_MODULE_HH
28#define EWOMS_VTK_BLACK_OIL_MICP_MODULE_HH
29
30#include <dune/common/fvector.hh>
31
32#include <opm/material/densead/Math.hpp>
33
35
37
40
43
44namespace Opm::Parameters {
45
46// set default values for what quantities to output
47struct VtkWriteMicrobialConcentration { static constexpr bool value = true; };
48struct VtkWriteOxygenConcentration { static constexpr bool value = true; };
49struct VtkWriteUreaConcentration { static constexpr bool value = true; };
50struct VtkWriteBiofilmConcentration { static constexpr bool value = true; };
51struct VtkWriteCalciteConcentration { static constexpr bool value = true; };
52
53} // namespace Opm::Parameters
54
55namespace Opm {
61template <class TypeTag>
63{
65
71
72 static const int vtkFormat = getPropValue<TypeTag, Properties::VtkOutputFormat>();
74
75 enum { enableMICP = getPropValue<TypeTag, Properties::EnableMICP>() };
76
78
79public:
80 VtkBlackOilMICPModule(const Simulator& simulator)
81 : ParentType(simulator)
82 { }
83
88 static void registerParameters()
89 {
90 if (!enableMICP)
91 return;
92
93 Parameters::Register<Parameters::VtkWriteMicrobialConcentration>
94 ("Include the concentration of the microbial component in the water phase "
95 "in the VTK output files");
96 Parameters::Register<Parameters::VtkWriteOxygenConcentration>
97 ("Include the concentration of the oxygen component in the water phase "
98 "in the VTK output files");
99 Parameters::Register<Parameters::VtkWriteUreaConcentration>
100 ("Include the concentration of the urea component in the water phase "
101 "in the VTK output files");
102 Parameters::Register<Parameters::VtkWriteBiofilmConcentration>
103 ("Include the biofilm volume fraction in the VTK output files");
104 Parameters::Register<Parameters::VtkWriteCalciteConcentration>
105 ("Include the calcite volume fraction in the VTK output files");
106 }
107
113 {
114 if (!Parameters::Get<Parameters::EnableVtkOutput>()) {
115 return;
116 }
117
118 if (!enableMICP)
119 return;
120
121 if (microbialConcentrationOutput_())
122 this->resizeScalarBuffer_(microbialConcentration_);
123 if (oxygenConcentrationOutput_())
124 this->resizeScalarBuffer_(oxygenConcentration_);
125 if (ureaConcentrationOutput_())
126 this->resizeScalarBuffer_(ureaConcentration_);
127 if (biofilmConcentrationOutput_())
128 this->resizeScalarBuffer_(biofilmConcentration_);
129 if (calciteConcentrationOutput_())
130 this->resizeScalarBuffer_(calciteConcentration_);
131 }
132
137 void processElement(const ElementContext& elemCtx)
138 {
139 if (!Parameters::Get<Parameters::EnableVtkOutput>()) {
140 return;
141 }
142
143 if (!enableMICP)
144 return;
145
146 for (unsigned dofIdx = 0; dofIdx < elemCtx.numPrimaryDof(/*timeIdx=*/0); ++dofIdx) {
147 const auto& intQuants = elemCtx.intensiveQuantities(dofIdx, /*timeIdx=*/0);
148 unsigned globalDofIdx = elemCtx.globalSpaceIndex(dofIdx, /*timeIdx=*/0);
149
150 if (microbialConcentrationOutput_())
151 microbialConcentration_[globalDofIdx] =
152 scalarValue(intQuants.microbialConcentration());
153
154 if (oxygenConcentrationOutput_())
155 oxygenConcentration_[globalDofIdx] =
156 scalarValue(intQuants.oxygenConcentration());
157
158 if (ureaConcentrationOutput_())
159 ureaConcentration_[globalDofIdx] =
160 10 * scalarValue(intQuants.ureaConcentration());//Multypliging by scaling factor 10 (see WellInterface_impl.hpp)
161
162 if (biofilmConcentrationOutput_())
163 biofilmConcentration_[globalDofIdx] =
164 scalarValue(intQuants.biofilmConcentration());
165
166 if (calciteConcentrationOutput_())
167 calciteConcentration_[globalDofIdx] =
168 scalarValue(intQuants.calciteConcentration());
169
170 }
171 }
172
177 {
178 VtkMultiWriter *vtkWriter = dynamic_cast<VtkMultiWriter*>(&baseWriter);
179 if (!vtkWriter)
180 return;
181
182 if (!enableMICP)
183 return;
184
185 if (microbialConcentrationOutput_())
186 this->commitScalarBuffer_(baseWriter, "microbial concentration", microbialConcentration_);
187
188 if (oxygenConcentrationOutput_())
189 this->commitScalarBuffer_(baseWriter, "oxygen concentration", oxygenConcentration_);
190
191 if (ureaConcentrationOutput_())
192 this->commitScalarBuffer_(baseWriter, "urea concentration", ureaConcentration_);
193
194 if (biofilmConcentrationOutput_())
195 this->commitScalarBuffer_(baseWriter, "biofilm fraction", biofilmConcentration_);
196
197 if (calciteConcentrationOutput_())
198 this->commitScalarBuffer_(baseWriter, "calcite fraction", calciteConcentration_);
199
200 }
201
202private:
203 static bool microbialConcentrationOutput_()
204 {
205 static bool val = Parameters::Get<Parameters::VtkWriteMicrobialConcentration>();
206 return val;
207 }
208
209 static bool oxygenConcentrationOutput_()
210 {
211 static bool val = Parameters::Get<Parameters::VtkWriteOxygenConcentration>();
212 return val;
213 }
214
215 static bool ureaConcentrationOutput_()
216 {
217 static bool val = Parameters::Get<Parameters::VtkWriteUreaConcentration>();
218 return val;
219 }
220
221 static bool biofilmConcentrationOutput_()
222 {
223 static bool val = Parameters::Get<Parameters::VtkWriteBiofilmConcentration>();
224 return val;
225 }
226
227 static bool calciteConcentrationOutput_()
228 {
229 static bool val = Parameters::Get<Parameters::VtkWriteCalciteConcentration>();
230 return val;
231 }
232
233 ScalarBuffer microbialConcentration_;
234 ScalarBuffer oxygenConcentration_;
235 ScalarBuffer ureaConcentration_;
236 ScalarBuffer biofilmConcentration_;
237 ScalarBuffer calciteConcentration_;
238};
239
240} // namespace Opm
241
242#endif
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 MICP model's related quantities.
Definition: vtkblackoilmicpmodule.hh:63
void processElement(const ElementContext &elemCtx)
Modify the internal buffers according to the intensive quantities relevant for an element.
Definition: vtkblackoilmicpmodule.hh:137
void commitBuffers(BaseOutputWriter &baseWriter)
Add all buffers to the VTK output writer.
Definition: vtkblackoilmicpmodule.hh:176
VtkBlackOilMICPModule(const Simulator &simulator)
Definition: vtkblackoilmicpmodule.hh:80
void allocBuffers()
Allocate memory for the scalar fields we would like to write to the VTK file.
Definition: vtkblackoilmicpmodule.hh:112
static void registerParameters()
Register all run-time parameters for the multi-phase VTK output module.
Definition: vtkblackoilmicpmodule.hh:88
Simplifies writing multi-file VTK datasets.
Definition: vtkmultiwriter.hh:66
Declare the properties used by the infrastructure code of the finite volume discretizations.
Definition: blackoilnewtonmethodparameters.hh:31
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: vtkblackoilmicpmodule.hh:50
static constexpr bool value
Definition: vtkblackoilmicpmodule.hh:50
Definition: vtkblackoilmicpmodule.hh:51
static constexpr bool value
Definition: vtkblackoilmicpmodule.hh:51
Definition: vtkblackoilmicpmodule.hh:47
static constexpr bool value
Definition: vtkblackoilmicpmodule.hh:47
Definition: vtkblackoilmicpmodule.hh:48
static constexpr bool value
Definition: vtkblackoilmicpmodule.hh:48
Definition: vtkblackoilmicpmodule.hh:49
static constexpr bool value
Definition: vtkblackoilmicpmodule.hh:49