vtkblackoilmicpmodule.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_MICP_MODULE_HPP
28#define OPM_VTK_BLACK_OIL_MICP_MODULE_HPP
29
30#include <dune/common/fvector.hh>
31
32#include <opm/material/densead/Math.hpp>
33
35
37
41
44
45namespace Opm {
51template <class TypeTag>
53{
55
61
62 static const int vtkFormat = getPropValue<TypeTag, Properties::VtkOutputFormat>();
64
65 enum { enableMICP = getPropValue<TypeTag, Properties::EnableMICP>() };
66
68
69public:
70 VtkBlackOilMICPModule(const Simulator& simulator)
71 : ParentType(simulator)
72 {
73 if constexpr (enableMICP) {
74 params_.read();
75 }
76 }
77
82 static void registerParameters()
83 {
84 if constexpr (enableMICP) {
86 }
87 }
88
94 {
95 if constexpr (enableMICP) {
96 if (!Parameters::Get<Parameters::EnableVtkOutput>()) {
97 return;
98 }
99
100 if (params_.microbialConcentrationOutput_) {
101 this->resizeScalarBuffer_(microbialConcentration_);
102 }
103 if (params_.oxygenConcentrationOutput_) {
104 this->resizeScalarBuffer_(oxygenConcentration_);
105 }
106 if (params_.ureaConcentrationOutput_) {
107 this->resizeScalarBuffer_(ureaConcentration_);
108 }
109 if (params_.biofilmConcentrationOutput_) {
110 this->resizeScalarBuffer_(biofilmConcentration_);
111 }
112 if (params_.calciteConcentrationOutput_) {
113 this->resizeScalarBuffer_(calciteConcentration_);
114 }
115 }
116 }
117
122 void processElement(const ElementContext& elemCtx)
123 {
124 if constexpr (enableMICP) {
125 if (!Parameters::Get<Parameters::EnableVtkOutput>()) {
126 return;
127 }
128
129 for (unsigned dofIdx = 0; dofIdx < elemCtx.numPrimaryDof(/*timeIdx=*/0); ++dofIdx) {
130 const auto& intQuants = elemCtx.intensiveQuantities(dofIdx, /*timeIdx=*/0);
131 unsigned globalDofIdx = elemCtx.globalSpaceIndex(dofIdx, /*timeIdx=*/0);
132
133 if (params_.microbialConcentrationOutput_) {
134 microbialConcentration_[globalDofIdx] =
135 scalarValue(intQuants.microbialConcentration());
136 }
137
138 if (params_.oxygenConcentrationOutput_) {
139 oxygenConcentration_[globalDofIdx] =
140 scalarValue(intQuants.oxygenConcentration());
141 }
142
143 if (params_.ureaConcentrationOutput_) {
144 ureaConcentration_[globalDofIdx] =
145 10 * scalarValue(intQuants.ureaConcentration());//Multypliging by scaling factor 10 (see WellInterface_impl.hpp)
146 }
147
148 if (params_.biofilmConcentrationOutput_) {
149 biofilmConcentration_[globalDofIdx] =
150 scalarValue(intQuants.biofilmConcentration());
151 }
152
153 if (params_.calciteConcentrationOutput_) {
154 calciteConcentration_[globalDofIdx] =
155 scalarValue(intQuants.calciteConcentration());
156 }
157 }
158 }
159 }
160
165 {
166 if constexpr (enableMICP) {
167 VtkMultiWriter* vtkWriter = dynamic_cast<VtkMultiWriter*>(&baseWriter);
168 if (!vtkWriter) {
169 return;
170 }
171
172 if (params_.microbialConcentrationOutput_) {
173 this->commitScalarBuffer_(baseWriter, "microbial concentration", microbialConcentration_);
174 }
175
176 if (params_.oxygenConcentrationOutput_) {
177 this->commitScalarBuffer_(baseWriter, "oxygen concentration", oxygenConcentration_);
178 }
179
180 if (params_.ureaConcentrationOutput_) {
181 this->commitScalarBuffer_(baseWriter, "urea concentration", ureaConcentration_);
182 }
183
184 if (params_.biofilmConcentrationOutput_) {
185 this->commitScalarBuffer_(baseWriter, "biofilm fraction", biofilmConcentration_);
186 }
187
188 if (params_.calciteConcentrationOutput_) {
189 this->commitScalarBuffer_(baseWriter, "calcite fraction", calciteConcentration_);
190 }
191 }
192 }
193
194private:
195 VtkBlackoilMICPParams params_{};
196 ScalarBuffer microbialConcentration_{};
197 ScalarBuffer oxygenConcentration_{};
198 ScalarBuffer ureaConcentration_{};
199 ScalarBuffer biofilmConcentration_{};
200 ScalarBuffer calciteConcentration_{};
201};
202
203} // namespace Opm
204
205#endif // OPM_VTK_BLACKOIL_MICP_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 MICP model's related quantities.
Definition: vtkblackoilmicpmodule.hpp:53
void processElement(const ElementContext &elemCtx)
Modify the internal buffers according to the intensive quantities relevant for an element.
Definition: vtkblackoilmicpmodule.hpp:122
void commitBuffers(BaseOutputWriter &baseWriter)
Add all buffers to the VTK output writer.
Definition: vtkblackoilmicpmodule.hpp:164
VtkBlackOilMICPModule(const Simulator &simulator)
Definition: vtkblackoilmicpmodule.hpp:70
void allocBuffers()
Allocate memory for the scalar fields we would like to write to the VTK file.
Definition: vtkblackoilmicpmodule.hpp:93
static void registerParameters()
Register all run-time parameters for the multi-phase VTK output module.
Definition: vtkblackoilmicpmodule.hpp:82
Simplifies writing multi-file VTK datasets.
Definition: vtkmultiwriter.hh:66
Declare the properties used by the infrastructure code of the finite volume discretizations.
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.
Struct holding the parameters for VtkBlackoilMICPModule.
Definition: vtkblackoilmicpparams.hpp:47
bool oxygenConcentrationOutput_
Definition: vtkblackoilmicpparams.hpp:55
static void registerParameters()
Registers the parameters in parameter system.
void read()
Reads the parameter values from the parameter system.
bool ureaConcentrationOutput_
Definition: vtkblackoilmicpparams.hpp:56
bool calciteConcentrationOutput_
Definition: vtkblackoilmicpparams.hpp:58
bool microbialConcentrationOutput_
Definition: vtkblackoilmicpparams.hpp:54
bool biofilmConcentrationOutput_
Definition: vtkblackoilmicpparams.hpp:57