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
32#include <opm/models/blackoil/blackoilproperties.hh>
33#include <opm/models/io/baseoutputmodule.hh>
34#include <opm/models/io/vtkmultiwriter.hh>
35#include <opm/models/utils/parametersystem.hh>
36#include <opm/models/utils/propertysystem.hh>
37
38#include <string>
39#include <vector>
40
41namespace Opm::Properties {
42
43// create new type tag for the VTK tracer output
44namespace TTag {
45struct VtkTracer {};
46}
47
48// create the property tags needed for the tracer model
49template<class TypeTag, class MyTypeTag>
51 using type = UndefinedProperty;
52};
53
54// set default values for what quantities to output
55template<class TypeTag>
56struct VtkWriteTracerConcentration<TypeTag, TTag::VtkTracer> {
57 static constexpr bool value = false;
58};
59
60} // namespace Opm::Properties
61
62namespace Opm {
68 template <class TypeTag>
69 class VtkTracerModule : public BaseOutputModule<TypeTag>
70 {
71 using ParentType = BaseOutputModule<TypeTag>;
72
73 using Simulator = GetPropType<TypeTag, Properties::Simulator>;
74 using Scalar = GetPropType<TypeTag, Properties::Scalar>;
75 using Evaluation = GetPropType<TypeTag, Properties::Evaluation>;
76 using ElementContext = GetPropType<TypeTag, Properties::ElementContext>;
77
78 using GridView = GetPropType<TypeTag, Properties::GridView>;
79 using FluidSystem = GetPropType<TypeTag, Properties::FluidSystem>;
80
81 static constexpr int vtkFormat = getPropValue<TypeTag, Properties::VtkOutputFormat>();
82 using VtkMultiWriter = ::Opm::VtkMultiWriter<GridView, vtkFormat>;
83
84
85 using ScalarBuffer = typename ParentType::ScalarBuffer;
86
87 public:
88 VtkTracerModule(const Simulator& simulator)
89 : ParentType(simulator)
90 { }
91
96 static void registerParameters()
97 {
98 Parameters::registerParam<TypeTag, Properties::VtkWriteTracerConcentration>
99 ("Include the tracer concentration in the VTK output files");
100 }
101
107 {
108 if (eclTracerConcentrationOutput_()){
109 const auto& tracerModel = this->simulator_.problem().tracerModel();
110 eclTracerConcentration_.resize(tracerModel.numTracers());
111 for (std::size_t tracerIdx = 0; tracerIdx < eclTracerConcentration_.size(); ++tracerIdx) {
112
113 this->resizeScalarBuffer_(eclTracerConcentration_[tracerIdx]);
114 }
115 }
116
117 }
118
123 void processElement(const ElementContext& elemCtx)
124 {
125 if (!Parameters::get<TypeTag, Properties::EnableVtkOutput>())
126 return;
127
128 const auto& tracerModel = elemCtx.problem().tracerModel();
129
130 for (unsigned dofIdx = 0; dofIdx < elemCtx.numPrimaryDof(/*timeIdx=*/0); ++dofIdx) {
131 unsigned globalDofIdx = elemCtx.globalSpaceIndex(dofIdx, /*timeIdx=*/0);
132
133 if (eclTracerConcentrationOutput_()){
134 for (std::size_t tracerIdx = 0; tracerIdx < eclTracerConcentration_.size(); ++tracerIdx) {
135 eclTracerConcentration_[tracerIdx][globalDofIdx] = tracerModel.tracerConcentration(tracerIdx, globalDofIdx);
136 }
137 }
138 }
139 }
140
144 void commitBuffers(BaseOutputWriter& baseWriter)
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 for (std::size_t tracerIdx = 0; tracerIdx < eclTracerConcentration_.size(); ++tracerIdx) {
153 const std::string tmp = "tracerConcentration_" + tracerModel.name(tracerIdx);
154 this->commitScalarBuffer_(baseWriter,tmp.c_str(), eclTracerConcentration_[tracerIdx]);
155 }
156 }
157
158
159
160 }
161
162 private:
163 static bool eclTracerConcentrationOutput_()
164 {
165 static bool val = Parameters::get<TypeTag, Properties::VtkWriteTracerConcentration>();
166 return val;
167 }
168
169
170 std::vector<ScalarBuffer> eclTracerConcentration_;
171 };
172} // namespace Opm
173
174#endif // OPM_VTK_TRACER_MODULE_HPP
VTK output module for the tracer model's parameters.
Definition: VtkTracerModule.hpp:70
static void registerParameters()
Register all run-time parameters for the tracer VTK output module.
Definition: VtkTracerModule.hpp:96
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:123
VtkTracerModule(const Simulator &simulator)
Definition: VtkTracerModule.hpp:88
void allocBuffers()
Allocate memory for the scalar fields we would like to write to the VTK file.
Definition: VtkTracerModule.hpp:106
Definition: AluGridVanguard.hpp:57
Definition: BlackoilPhases.hpp:27
Definition: VtkTracerModule.hpp:45
Definition: VtkTracerModule.hpp:50
UndefinedProperty type
Definition: VtkTracerModule.hpp:51