VtkTracerModule.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_TRACER_MODULE_HPP
28#define OPM_VTK_TRACER_MODULE_HPP
29
30#include <dune/common/fvector.hh>
31
38
39#include <string>
40#include <vector>
41
42namespace Opm::Parameters {
43
44// set default values for what quantities to output
45struct VtkWriteTracerConcentration { static constexpr bool value = false; };
46
47} // namespace Opm::Parameters
48
49namespace Opm {
55 template <class TypeTag>
56 class VtkTracerModule : public BaseOutputModule<TypeTag>
57 {
59
64
67
68 static constexpr int vtkFormat = getPropValue<TypeTag, Properties::VtkOutputFormat>();
70
71
72 using ScalarBuffer = typename ParentType::ScalarBuffer;
73
74 public:
75 VtkTracerModule(const Simulator& simulator)
76 : ParentType(simulator)
77 { }
78
83 static void registerParameters()
84 {
85 Parameters::Register<Parameters::VtkWriteTracerConcentration>
86 ("Include the tracer concentration in the VTK output files");
87 }
88
94 {
95 if (eclTracerConcentrationOutput_()){
96 const auto& tracerModel = this->simulator_.problem().tracerModel();
97 eclFreeTracerConcentration_.resize(tracerModel.numTracers());
98 eclSolTracerConcentration_.resize(tracerModel.numTracers());
99 const auto& enableSolTracers = tracerModel.enableSolTracers();
100
101 for (std::size_t tracerIdx = 0; tracerIdx < eclFreeTracerConcentration_.size(); ++tracerIdx) {
102 this->resizeScalarBuffer_(eclFreeTracerConcentration_[tracerIdx]);
103 if (enableSolTracers[tracerIdx])
104 this->resizeScalarBuffer_(eclSolTracerConcentration_[tracerIdx]);
105 }
106 }
107
108 }
109
114 void processElement(const ElementContext& elemCtx)
115 {
116 if (!Parameters::Get<Parameters::EnableVtkOutput>()) {
117 return;
118 }
119
120 if (eclTracerConcentrationOutput_()) {
121 const auto& tracerModel = elemCtx.problem().tracerModel();
122 const auto& enableSolTracers = tracerModel.enableSolTracers();
123
124 for (std::size_t tracerIdx = 0; tracerIdx < eclFreeTracerConcentration_.size(); ++tracerIdx) {
125 // free tracer
126 for (unsigned dofIdx = 0; dofIdx < elemCtx.numPrimaryDof(/*timeIdx=*/0); ++dofIdx) {
127 unsigned globalDofIdx = elemCtx.globalSpaceIndex(dofIdx, /*timeIdx=*/0);
128 eclFreeTracerConcentration_[tracerIdx][globalDofIdx] = tracerModel.freeTracerConcentration(tracerIdx, globalDofIdx);
129 }
130 // solution tracer (only if it exist)
131 if (enableSolTracers[tracerIdx]) {
132 for (unsigned dofIdx = 0; dofIdx < elemCtx.numPrimaryDof(/*timeIdx=*/0); ++dofIdx) {
133 unsigned globalDofIdx = elemCtx.globalSpaceIndex(dofIdx, /*timeIdx=*/0);
134 eclSolTracerConcentration_[tracerIdx][globalDofIdx] = tracerModel.solTracerConcentration(tracerIdx, globalDofIdx);
135 }
136 }
137 }
138 }
139 }
140
145 {
146 VtkMultiWriter *vtkWriter = dynamic_cast<VtkMultiWriter*>(&baseWriter);
147 if (!vtkWriter)
148 return;
149
150 if (eclTracerConcentrationOutput_()){
151 const auto& tracerModel = this->simulator_.problem().tracerModel();
152 const auto& enableSolTracers = tracerModel.enableSolTracers();
153
154 for (std::size_t tracerIdx = 0; tracerIdx < eclFreeTracerConcentration_.size(); ++tracerIdx) {
155 const std::string tmp = "freeTracerConcentration_" + tracerModel.name(tracerIdx);
156 this->commitScalarBuffer_(baseWriter, tmp.c_str(), eclFreeTracerConcentration_[tracerIdx]);
157 if (enableSolTracers[tracerIdx]) {
158 const std::string tmp2 = "solTracerConcentration_" + tracerModel.name(tracerIdx);
159 this->commitScalarBuffer_(baseWriter, tmp2.c_str(), eclSolTracerConcentration_[tracerIdx]);
160 }
161 }
162 }
163 }
164
165 private:
166 static bool eclTracerConcentrationOutput_()
167 {
168 static bool val = Parameters::Get<Parameters::VtkWriteTracerConcentration>();
169 return val;
170 }
171
172 std::vector<ScalarBuffer> eclFreeTracerConcentration_;
173 std::vector<ScalarBuffer> eclSolTracerConcentration_;
174 };
175
176} // namespace Opm
177
178#endif // OPM_VTK_TRACER_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
const Simulator & simulator_
Definition: baseoutputmodule.hh:434
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
Simplifies writing multi-file VTK datasets.
Definition: vtkmultiwriter.hh:66
VTK output module for the tracer model's parameters.
Definition: VtkTracerModule.hpp:57
static void registerParameters()
Register all run-time parameters for the tracer VTK output module.
Definition: VtkTracerModule.hpp:83
void commitBuffers(BaseOutputWriter &baseWriter)
Add all buffers to the VTK output writer.
Definition: VtkTracerModule.hpp:144
void processElement(const ElementContext &elemCtx)
Modify the internal buffers according to the intensive quantities relevant for an element.
Definition: VtkTracerModule.hpp:114
VtkTracerModule(const Simulator &simulator)
Definition: VtkTracerModule.hpp:75
void allocBuffers()
Allocate memory for the scalar fields we would like to write to the VTK file.
Definition: VtkTracerModule.hpp:93
Declare the properties used by the infrastructure code of the finite volume discretizations.
Definition: blackoilnewtonmethodparams.hpp: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: VtkTracerModule.hpp:45
static constexpr bool value
Definition: VtkTracerModule.hpp:45