version 3.10-dev
staggeredvtkoutputmodule.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// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
12#ifndef DUMUX_STAGGERED_VTK_OUTPUT_MODULE_HH
13#define DUMUX_STAGGERED_VTK_OUTPUT_MODULE_HH
14
15#include <dune/common/fvector.hh>
16
21
22namespace Dumux {
23
24template<class Scalar, class GlobalPosition>
25class PointCloudVtkWriter;
26
35template<class GridVariables, class SolutionVector>
37: public VtkOutputModule<GridVariables, SolutionVector>
38{
40 using GridGeometry = typename GridVariables::GridGeometry;
41 using GridView = typename GridGeometry::GridView;
42 using Scalar = typename GridVariables::Scalar;
43 using FaceVariables = typename GridVariables::GridFaceVariables::FaceVariables;
44 using FVElementGeometry = typename GridGeometry::LocalView;
45 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
46
47 using Element = typename GridView::template Codim<0>::Entity;
48 using GlobalPosition = typename Element::Geometry::GlobalCoordinate;
49
50 struct FaceVarScalarDataInfo { std::function<Scalar(const FaceVariables&)> get; std::string name; };
51 struct FaceVarVectorDataInfo { std::function<GlobalPosition(const SubControlVolumeFace& scvf, const FaceVariables&)> get; std::string name; };
52
53 struct FaceFieldScalarDataInfo
54 {
55 FaceFieldScalarDataInfo(const std::vector<Scalar>& f, const std::string& n) : data(f), name(n) {}
56 const std::vector<Scalar>& data;
57 const std::string name;
58 };
59
60 struct FaceFieldVectorDataInfo
61 {
62 FaceFieldVectorDataInfo(const std::vector<GlobalPosition>& f, const std::string& n) : data(f), name(n) {}
63 const std::vector<GlobalPosition>& data;
64 const std::string name;
65 };
66
67public:
68
69 template<class Sol>
71 const Sol& sol,
72 const std::string& name,
73 const std::string& paramGroup = "",
74 Dune::VTK::DataMode dm = Dune::VTK::conforming,
75 bool verbose = true)
77 , faceWriter_(std::make_shared<PointCloudVtkWriter<Scalar, GlobalPosition>>(coordinates_))
78 , sequenceWriter_(faceWriter_, name + "-face", "","",
79 gridVariables.curGridVolVars().problem().gridGeometry().gridView().comm().rank(),
80 gridVariables.curGridVolVars().problem().gridGeometry().gridView().comm().size() )
81
82 {
83 static_assert(std::is_same<Sol, SolutionVector>::value, "Make sure that sol has the same type as SolutionVector."
84 "Use StaggeredVtkOutputModule<GridVariables, decltype(sol)> when calling the constructor.");
85
86 // enable velocity output per default
88 writeFaceVars_ = getParamFromGroup<bool>(paramGroup, "Vtk.WriteFaceData", false);
89 coordinatesInitialized_ = false;
90 }
91
96
100 void addFaceField(const std::vector<Scalar>& v, const std::string& name)
101 {
102 if (v.size() == this->gridGeometry().gridView().size(1))
103 faceFieldScalarDataInfo_.emplace_back(v, name);
104 else
105 DUNE_THROW(Dune::RangeError, "Size mismatch of added field!");
106 }
107
111 void addFaceField(const std::vector<GlobalPosition>& v, const std::string& name)
112 {
113 if (v.size() == this->gridGeometry().gridView().size(1))
114 faceFieldVectorDataInfo_.emplace_back(v, name);
115 else
116 DUNE_THROW(Dune::RangeError, "Size mismatch of added field!");
117 }
118
122 void addFaceVariable(std::function<Scalar(const FaceVariables&)>&& f, const std::string& name)
123 {
124 faceVarScalarDataInfo_.push_back(FaceVarScalarDataInfo{f, name});
125 }
126
130 void addFaceVariable(std::function<GlobalPosition(const SubControlVolumeFace& scvf, const FaceVariables&)>&& f, const std::string& name)
131 {
132 faceVarVectorDataInfo_.push_back(FaceVarVectorDataInfo{f, name});
133 }
134
138 void write(double time, Dune::VTK::OutputType type = Dune::VTK::ascii)
139 {
140 ParentType::write(time, type);
141 if(writeFaceVars_)
142 getFaceDataAndWrite_(time);
143 }
144
145
146private:
147
149 void updateCoordinates_()
150 {
151 coordinates_.resize(this->gridGeometry().numFaceDofs());
152 for(auto&& facet : facets(this->gridGeometry().gridView()))
153 {
154 const int dofIdxGlobal = this->gridGeometry().gridView().indexSet().index(facet);
155 coordinates_[dofIdxGlobal] = facet.geometry().center();
156 }
157 coordinatesInitialized_ = true;
158 }
159
162 void getFaceDataAndWrite_(const Scalar time)
163 {
164 const auto numPoints = this->gridGeometry().numFaceDofs();
165
166 // make sure not to iterate over the same dofs twice
167 std::vector<bool> dofVisited(numPoints, false);
168
169 // get fields for all primary coordinates and variables
170 if(!coordinatesInitialized_)
171 updateCoordinates_();
172
173 // prepare some containers to store the relevant data
174 std::vector<std::vector<Scalar>> faceVarScalarData;
175 std::vector<std::vector<GlobalPosition>> faceVarVectorData;
176
177 if(!faceVarScalarDataInfo_.empty())
178 faceVarScalarData.resize(faceVarScalarDataInfo_.size(), std::vector<Scalar>(numPoints));
179
180 if(!faceVarVectorDataInfo_.empty())
181 faceVarVectorData.resize(faceVarVectorDataInfo_.size(), std::vector<GlobalPosition>(numPoints));
182
183 auto fvGeometry = localView(this->gridGeometry());
184 auto elemFaceVars = localView(this->gridVariables().curGridFaceVars());
185 for (const auto& element : elements(this->gridGeometry().gridView(), Dune::Partitions::interior))
186 {
187 if (!faceVarScalarDataInfo_.empty() || !faceVarVectorDataInfo_.empty())
188 {
189 fvGeometry.bind(element);
190 elemFaceVars.bindElement(element, fvGeometry, this->sol());
191
192 for (auto&& scvf : scvfs(fvGeometry))
193 {
194 const auto dofIdxGlobal = scvf.dofIndex();
195 if(dofVisited[dofIdxGlobal])
196 continue;
197
198 dofVisited[dofIdxGlobal] = true;
199
200 const auto& faceVars = elemFaceVars[scvf];
201
202 // get the scalar-valued data
203 for (std::size_t i = 0; i < faceVarScalarDataInfo_.size(); ++i)
204 faceVarScalarData[i][dofIdxGlobal] = faceVarScalarDataInfo_[i].get(faceVars);
205
206 // get the vector-valued data
207 for (std::size_t i = 0; i < faceVarVectorDataInfo_.size(); ++i)
208 faceVarVectorData[i][dofIdxGlobal] = faceVarVectorDataInfo_[i].get(scvf, faceVars);
209 }
210 }
211 }
212
213 // transfer the data to the point writer
214 if(!faceVarScalarDataInfo_.empty())
215 for (std::size_t i = 0; i < faceVarScalarDataInfo_.size(); ++i)
216 faceWriter_->addPointData(faceVarScalarData[i], faceVarScalarDataInfo_[i].name);
217
218 if(!faceVarVectorDataInfo_.empty())
219 for (std::size_t i = 0; i < faceVarVectorDataInfo_.size(); ++i)
220 faceWriter_->addPointData(faceVarVectorData[i], faceVarVectorDataInfo_[i].name);
221
222 // account for the custom fields
223 for(const auto& field: faceFieldScalarDataInfo_)
224 faceWriter_->addPointData(field.data, field.name);
225
226 for(const auto& field: faceFieldVectorDataInfo_)
227 faceWriter_->addPointData(field.data, field.name);
228
229 // write for the current time step
230 sequenceWriter_.write(time);
231
232 // clear coordinates to save some memory
233 coordinates_.clear();
234 coordinates_.shrink_to_fit();
235 coordinatesInitialized_ = false;
236 }
237
238
239 std::shared_ptr<PointCloudVtkWriter<Scalar, GlobalPosition>> faceWriter_;
240
241 VTKSequenceWriter<PointCloudVtkWriter<Scalar, GlobalPosition>> sequenceWriter_;
242
243 bool writeFaceVars_;
244
245 std::vector<GlobalPosition> coordinates_;
246 bool coordinatesInitialized_;
247
248 std::vector<FaceVarScalarDataInfo> faceVarScalarDataInfo_;
249 std::vector<FaceVarVectorDataInfo> faceVarVectorDataInfo_;
250
251 std::vector<FaceFieldScalarDataInfo> faceFieldScalarDataInfo_;
252 std::vector<FaceFieldVectorDataInfo> faceFieldVectorDataInfo_;
253
254};
255
256} // end namespace Dumux
257
258#endif
A VTK output module to simplify writing dumux simulation data to VTK format.
Definition: pointcloudvtkwriter.hh:39
Velocity output for staggered free-flow models.
Definition: discretization/staggered/freeflow/velocityoutput.hh:26
A VTK output module to simplify writing dumux simulation data to VTK format Specialization for stagge...
Definition: staggeredvtkoutputmodule.hh:38
void addFaceField(const std::vector< Scalar > &v, const std::string &name)
Definition: staggeredvtkoutputmodule.hh:100
void write(double time, Dune::VTK::OutputType type=Dune::VTK::ascii)
Definition: staggeredvtkoutputmodule.hh:138
StaggeredVtkOutputModule(const GridVariables &gridVariables, const Sol &sol, const std::string &name, const std::string &paramGroup="", Dune::VTK::DataMode dm=Dune::VTK::conforming, bool verbose=true)
Definition: staggeredvtkoutputmodule.hh:70
void addFaceField(const std::vector< GlobalPosition > &v, const std::string &name)
Definition: staggeredvtkoutputmodule.hh:111
void addFaceVariable(std::function< GlobalPosition(const SubControlVolumeFace &scvf, const FaceVariables &)> &&f, const std::string &name)
Definition: staggeredvtkoutputmodule.hh:130
void addFaceVariable(std::function< Scalar(const FaceVariables &)> &&f, const std::string &name)
Definition: staggeredvtkoutputmodule.hh:122
const std::string & paramGroup() const
the parameter group for getting parameter from the parameter tree
Definition: io/vtkoutputmodule.hh:79
void write(double time, Dune::VTK::OutputType type=Dune::VTK::ascii)
Definition: io/vtkoutputmodule.hh:165
const std::string & name() const
Definition: io/vtkoutputmodule.hh:187
bool verbose() const
Definition: io/vtkoutputmodule.hh:186
A VTK output module to simplify writing dumux simulation data to VTK format.
Definition: io/vtkoutputmodule.hh:300
const auto & problem() const
Definition: io/vtkoutputmodule.hh:383
const GridVariables & gridVariables() const
Definition: io/vtkoutputmodule.hh:384
void addVelocityOutput(std::shared_ptr< VelocityOutputType > velocityOutput)
Add a velocity output policy.
Definition: io/vtkoutputmodule.hh:358
const SolutionVector & sol() const
Definition: io/vtkoutputmodule.hh:386
const GridGeometry & gridGeometry() const
Definition: io/vtkoutputmodule.hh:385
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:26
A VTK output module to simplify writing dumux simulation data to VTK format.
Definition: adapt.hh:17
A VTK writer specialized for staggered grid implementations with dofs on the faces.
Base class to write pvd-files which contains a list of all collected vtk-files. This is a modified ve...