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 <opm/material/densead/Math.hpp>
31
32#include "vtkmultiwriter.hh"
33#include "baseoutputmodule.hh"
34
38
39#include <dune/common/fvector.hh>
40
41#include <cstdio>
42
43namespace Opm::Properties {
44
45namespace TTag {
46
47// create new type tag for the VTK multi-phase output
49
50} // namespace TTag
51
52// create the property tags needed for the MICP output module
53template<class TypeTag, class MyTypeTag>
55template<class TypeTag, class MyTypeTag>
57template<class TypeTag, class MyTypeTag>
59template<class TypeTag, class MyTypeTag>
61template<class TypeTag, class MyTypeTag>
63
64// set default values for what quantities to output
65template<class TypeTag>
66struct VtkWriteMicrobialConcentration<TypeTag, TTag::VtkBlackOilMICP> { static constexpr bool value = true; };
67template<class TypeTag>
68struct VtkWriteOxygenConcentration<TypeTag, TTag::VtkBlackOilMICP> { static constexpr bool value = true; };
69template<class TypeTag>
70struct VtkWriteUreaConcentration<TypeTag, TTag::VtkBlackOilMICP> { static constexpr bool value = true; };
71template<class TypeTag>
72struct VtkWriteBiofilmConcentration<TypeTag, TTag::VtkBlackOilMICP> { static constexpr bool value = true; };
73template<class TypeTag>
74struct VtkWriteCalciteConcentration<TypeTag, TTag::VtkBlackOilMICP> { static constexpr bool value = true; };
75
76} // namespace Opm::Properties
77
78namespace Opm {
84template <class TypeTag>
86{
88
94
95 static const int vtkFormat = getPropValue<TypeTag, Properties::VtkOutputFormat>();
97
98 enum { enableMICP = getPropValue<TypeTag, Properties::EnableMICP>() };
99
100 using ScalarBuffer = typename ParentType::ScalarBuffer;
101
102public:
103 VtkBlackOilMICPModule(const Simulator& simulator)
104 : ParentType(simulator)
105 { }
106
111 static void registerParameters()
112 {
113 if (!enableMICP)
114 return;
115
116 Parameters::registerParam<TypeTag, Properties::VtkWriteMicrobialConcentration>
117 ("Include the concentration of the microbial component in the water phase "
118 "in the VTK output files");
119 Parameters::registerParam<TypeTag, Properties::VtkWriteOxygenConcentration>
120 ("Include the concentration of the oxygen component in the water phase "
121 "in the VTK output files");
122 Parameters::registerParam<TypeTag, Properties::VtkWriteUreaConcentration>
123 ("Include the concentration of the urea component in the water phase "
124 "in the VTK output files");
125 Parameters::registerParam<TypeTag, Properties::VtkWriteBiofilmConcentration>
126 ("Include the biofilm volume fraction in the VTK output files");
127 Parameters::registerParam<TypeTag, Properties::VtkWriteCalciteConcentration>
128 ("Include the calcite volume fraction in the VTK output files");
129 }
130
136 {
137 if (!Parameters::get<TypeTag, Properties::EnableVtkOutput>())
138 return;
139
140 if (!enableMICP)
141 return;
142
143 if (microbialConcentrationOutput_())
144 this->resizeScalarBuffer_(microbialConcentration_);
145 if (oxygenConcentrationOutput_())
146 this->resizeScalarBuffer_(oxygenConcentration_);
147 if (ureaConcentrationOutput_())
148 this->resizeScalarBuffer_(ureaConcentration_);
149 if (biofilmConcentrationOutput_())
150 this->resizeScalarBuffer_(biofilmConcentration_);
151 if (calciteConcentrationOutput_())
152 this->resizeScalarBuffer_(calciteConcentration_);
153 }
154
159 void processElement(const ElementContext& elemCtx)
160 {
161 if (!Parameters::get<TypeTag, Properties::EnableVtkOutput>())
162 return;
163
164 if (!enableMICP)
165 return;
166
167 for (unsigned dofIdx = 0; dofIdx < elemCtx.numPrimaryDof(/*timeIdx=*/0); ++dofIdx) {
168 const auto& intQuants = elemCtx.intensiveQuantities(dofIdx, /*timeIdx=*/0);
169 unsigned globalDofIdx = elemCtx.globalSpaceIndex(dofIdx, /*timeIdx=*/0);
170
171 if (microbialConcentrationOutput_())
172 microbialConcentration_[globalDofIdx] =
173 scalarValue(intQuants.microbialConcentration());
174
175 if (oxygenConcentrationOutput_())
176 oxygenConcentration_[globalDofIdx] =
177 scalarValue(intQuants.oxygenConcentration());
178
179 if (ureaConcentrationOutput_())
180 ureaConcentration_[globalDofIdx] =
181 10 * scalarValue(intQuants.ureaConcentration());//Multypliging by scaling factor 10 (see WellInterface_impl.hpp)
182
183 if (biofilmConcentrationOutput_())
184 biofilmConcentration_[globalDofIdx] =
185 scalarValue(intQuants.biofilmConcentration());
186
187 if (calciteConcentrationOutput_())
188 calciteConcentration_[globalDofIdx] =
189 scalarValue(intQuants.calciteConcentration());
190
191 }
192 }
193
198 {
199 VtkMultiWriter *vtkWriter = dynamic_cast<VtkMultiWriter*>(&baseWriter);
200 if (!vtkWriter)
201 return;
202
203 if (!enableMICP)
204 return;
205
206 if (microbialConcentrationOutput_())
207 this->commitScalarBuffer_(baseWriter, "microbial concentration", microbialConcentration_);
208
209 if (oxygenConcentrationOutput_())
210 this->commitScalarBuffer_(baseWriter, "oxygen concentration", oxygenConcentration_);
211
212 if (ureaConcentrationOutput_())
213 this->commitScalarBuffer_(baseWriter, "urea concentration", ureaConcentration_);
214
215 if (biofilmConcentrationOutput_())
216 this->commitScalarBuffer_(baseWriter, "biofilm fraction", biofilmConcentration_);
217
218 if (calciteConcentrationOutput_())
219 this->commitScalarBuffer_(baseWriter, "calcite fraction", calciteConcentration_);
220
221 }
222
223private:
224 static bool microbialConcentrationOutput_()
225 {
226 static bool val = Parameters::get<TypeTag, Properties::VtkWriteMicrobialConcentration>();
227 return val;
228 }
229
230 static bool oxygenConcentrationOutput_()
231 {
232 static bool val = Parameters::get<TypeTag, Properties::VtkWriteOxygenConcentration>();
233 return val;
234 }
235
236 static bool ureaConcentrationOutput_()
237 {
238 static bool val = Parameters::get<TypeTag, Properties::VtkWriteUreaConcentration>();
239 return val;
240 }
241
242 static bool biofilmConcentrationOutput_()
243 {
244 static bool val = Parameters::get<TypeTag, Properties::VtkWriteBiofilmConcentration>();
245 return val;
246 }
247
248 static bool calciteConcentrationOutput_()
249 {
250 static bool val = Parameters::get<TypeTag, Properties::VtkWriteCalciteConcentration>();
251 return val;
252 }
253
254 ScalarBuffer microbialConcentration_;
255 ScalarBuffer oxygenConcentration_;
256 ScalarBuffer ureaConcentration_;
257 ScalarBuffer biofilmConcentration_;
258 ScalarBuffer calciteConcentration_;
259};
260} // namespace Opm
261
262#endif
Declares the properties required by the black oil model.
The base class for writer modules.
Definition: baseoutputmodule.hh:73
BaseOutputWriter::ScalarBuffer ScalarBuffer
Definition: baseoutputmodule.hh:91
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
VTK output module for the MICP model's related quantities.
Definition: vtkblackoilmicpmodule.hh:86
void processElement(const ElementContext &elemCtx)
Modify the internal buffers according to the intensive quantities relevant for an element.
Definition: vtkblackoilmicpmodule.hh:159
void commitBuffers(BaseOutputWriter &baseWriter)
Add all buffers to the VTK output writer.
Definition: vtkblackoilmicpmodule.hh:197
VtkBlackOilMICPModule(const Simulator &simulator)
Definition: vtkblackoilmicpmodule.hh:103
void allocBuffers()
Allocate memory for the scalar fields we would like to write to the VTK file.
Definition: vtkblackoilmicpmodule.hh:135
static void registerParameters()
Register all run-time parameters for the multi-phase VTK output module.
Definition: vtkblackoilmicpmodule.hh:111
Simplifies writing multi-file VTK datasets.
Definition: vtkmultiwriter.hh:66
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: vtkblackoilmicpmodule.hh:48
a tag to mark properties as undefined
Definition: propertysystem.hh:40
Definition: vtkblackoilmicpmodule.hh:60
Definition: vtkblackoilmicpmodule.hh:62
Definition: vtkblackoilmicpmodule.hh:54
Definition: vtkblackoilmicpmodule.hh:56
Definition: vtkblackoilmicpmodule.hh:58