12#ifndef DUMUX_STAGGERED_VTK_OUTPUT_MODULE_HH
13#define DUMUX_STAGGERED_VTK_OUTPUT_MODULE_HH
15#include <dune/common/fvector.hh>
24template<
class Scalar,
class GlobalPosition>
25class PointCloudVtkWriter;
35template<
class Gr
idVariables,
class SolutionVector>
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;
47 using Element =
typename GridView::template Codim<0>::Entity;
48 using GlobalPosition =
typename Element::Geometry::GlobalCoordinate;
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; };
53 struct FaceFieldScalarDataInfo
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;
60 struct FaceFieldVectorDataInfo
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;
72 const std::string&
name,
74 Dune::VTK::DataMode dm = Dune::VTK::conforming,
78 , sequenceWriter_(faceWriter_,
name +
"-face",
"",
"",
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.");
88 writeFaceVars_ = getParamFromGroup<bool>(
paramGroup,
"Vtk.WriteFaceData",
false);
89 coordinatesInitialized_ =
false;
102 if (v.size() == this->gridGeometry().gridView().size(1))
103 faceFieldScalarDataInfo_.emplace_back(v,
name);
105 DUNE_THROW(Dune::RangeError,
"Size mismatch of added field!");
113 if (v.size() == this->gridGeometry().gridView().size(1))
114 faceFieldVectorDataInfo_.emplace_back(v,
name);
116 DUNE_THROW(Dune::RangeError,
"Size mismatch of added field!");
124 faceVarScalarDataInfo_.push_back(FaceVarScalarDataInfo{f,
name});
130 void addFaceVariable(std::function<GlobalPosition(
const SubControlVolumeFace& scvf,
const FaceVariables&)>&& f,
const std::string&
name)
132 faceVarVectorDataInfo_.push_back(FaceVarVectorDataInfo{f,
name});
138 void write(
double time, Dune::VTK::OutputType type = Dune::VTK::ascii)
142 getFaceDataAndWrite_(time);
149 void updateCoordinates_()
151 coordinates_.resize(this->
gridGeometry().numFaceDofs());
152 for(
auto&& facet : facets(this->
gridGeometry().gridView()))
154 const int dofIdxGlobal = this->
gridGeometry().gridView().indexSet().index(facet);
155 coordinates_[dofIdxGlobal] = facet.geometry().center();
157 coordinatesInitialized_ =
true;
162 void getFaceDataAndWrite_(
const Scalar time)
164 const auto numPoints = this->
gridGeometry().numFaceDofs();
167 std::vector<bool> dofVisited(numPoints,
false);
170 if(!coordinatesInitialized_)
171 updateCoordinates_();
174 std::vector<std::vector<Scalar>> faceVarScalarData;
175 std::vector<std::vector<GlobalPosition>> faceVarVectorData;
177 if(!faceVarScalarDataInfo_.empty())
178 faceVarScalarData.resize(faceVarScalarDataInfo_.size(), std::vector<Scalar>(numPoints));
180 if(!faceVarVectorDataInfo_.empty())
181 faceVarVectorData.resize(faceVarVectorDataInfo_.size(), std::vector<GlobalPosition>(numPoints));
185 for (
const auto& element : elements(this->
gridGeometry().gridView(), Dune::Partitions::interior))
187 if (!faceVarScalarDataInfo_.empty() || !faceVarVectorDataInfo_.empty())
189 fvGeometry.bind(element);
190 elemFaceVars.bindElement(element, fvGeometry, this->
sol());
192 for (
auto&& scvf : scvfs(fvGeometry))
194 const auto dofIdxGlobal = scvf.dofIndex();
195 if(dofVisited[dofIdxGlobal])
198 dofVisited[dofIdxGlobal] =
true;
200 const auto& faceVars = elemFaceVars[scvf];
203 for (std::size_t i = 0; i < faceVarScalarDataInfo_.size(); ++i)
204 faceVarScalarData[i][dofIdxGlobal] = faceVarScalarDataInfo_[i].get(faceVars);
207 for (std::size_t i = 0; i < faceVarVectorDataInfo_.size(); ++i)
208 faceVarVectorData[i][dofIdxGlobal] = faceVarVectorDataInfo_[i].get(scvf, faceVars);
214 if(!faceVarScalarDataInfo_.empty())
215 for (std::size_t i = 0; i < faceVarScalarDataInfo_.size(); ++i)
216 faceWriter_->addPointData(faceVarScalarData[i], faceVarScalarDataInfo_[i].name);
218 if(!faceVarVectorDataInfo_.empty())
219 for (std::size_t i = 0; i < faceVarVectorDataInfo_.size(); ++i)
220 faceWriter_->addPointData(faceVarVectorData[i], faceVarVectorDataInfo_[i].name);
223 for(
const auto& field: faceFieldScalarDataInfo_)
224 faceWriter_->addPointData(field.data, field.name);
226 for(
const auto& field: faceFieldVectorDataInfo_)
227 faceWriter_->addPointData(field.data, field.name);
230 sequenceWriter_.write(time);
233 coordinates_.clear();
234 coordinates_.shrink_to_fit();
235 coordinatesInitialized_ =
false;
239 std::shared_ptr<PointCloudVtkWriter<Scalar, GlobalPosition>> faceWriter_;
241 VTKSequenceWriter<PointCloudVtkWriter<Scalar, GlobalPosition>> sequenceWriter_;
245 std::vector<GlobalPosition> coordinates_;
246 bool coordinatesInitialized_;
248 std::vector<FaceVarScalarDataInfo> faceVarScalarDataInfo_;
249 std::vector<FaceVarVectorDataInfo> faceVarVectorDataInfo_;
251 std::vector<FaceFieldScalarDataInfo> faceFieldScalarDataInfo_;
252 std::vector<FaceFieldVectorDataInfo> faceFieldVectorDataInfo_;
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 ¶mGroup="", 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.
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...