24#ifndef DUMUX_STAGGERED_VTK_OUTPUT_MODULE_HH
25#define DUMUX_STAGGERED_VTK_OUTPUT_MODULE_HH
27#include <dune/common/fvector.hh>
36template<
class Scalar,
class GlobalPosition>
37class PointCloudVtkWriter;
47template<
class Gr
idVariables,
class SolutionVector>
52 using GridGeometry =
typename GridVariables::GridGeometry;
53 using GridView =
typename GridGeometry::GridView;
54 using Scalar =
typename GridVariables::Scalar;
55 using FaceVariables =
typename GridVariables::GridFaceVariables::FaceVariables;
56 using FVElementGeometry =
typename GridGeometry::LocalView;
57 using SubControlVolumeFace =
typename FVElementGeometry::SubControlVolumeFace;
59 using Element =
typename GridView::template Codim<0>::Entity;
60 using GlobalPosition =
typename Element::Geometry::GlobalCoordinate;
62 struct FaceVarScalarDataInfo { std::function<Scalar(
const FaceVariables&)> get; std::string
name; };
63 struct FaceVarVectorDataInfo { std::function<GlobalPosition(
const SubControlVolumeFace& scvf,
const FaceVariables&)> get; std::string
name; };
65 struct FaceFieldScalarDataInfo
67 FaceFieldScalarDataInfo(
const std::vector<Scalar>& f,
const std::string& n) : data(f),
name(n) {}
68 const std::vector<Scalar>& data;
69 const std::string
name;
72 struct FaceFieldVectorDataInfo
74 FaceFieldVectorDataInfo(
const std::vector<GlobalPosition>& f,
const std::string& n) : data(f),
name(n) {}
75 const std::vector<GlobalPosition>& data;
76 const std::string
name;
84 const std::string&
name,
86 Dune::VTK::DataMode dm = Dune::VTK::conforming,
90 , sequenceWriter_(faceWriter_,
name +
"-face",
"",
"",
95 static_assert(std::is_same<Sol, SolutionVector>::value,
"Make sure that sol has the same type as SolutionVector."
96 "Use StaggeredVtkOutputModule<GridVariables, decltype(sol)> when calling the constructor.");
100 writeFaceVars_ = getParamFromGroup<bool>(
paramGroup,
"Vtk.WriteFaceData",
false);
101 coordinatesInitialized_ =
false;
114 if (v.size() == this->gridGeometry().gridView().size(1))
115 faceFieldScalarDataInfo_.emplace_back(v,
name);
117 DUNE_THROW(Dune::RangeError,
"Size mismatch of added field!");
125 if (v.size() == this->gridGeometry().gridView().size(1))
126 faceFieldVectorDataInfo_.emplace_back(v,
name);
128 DUNE_THROW(Dune::RangeError,
"Size mismatch of added field!");
136 faceVarScalarDataInfo_.push_back(FaceVarScalarDataInfo{f,
name});
142 void addFaceVariable(std::function<GlobalPosition(
const SubControlVolumeFace& scvf,
const FaceVariables&)>&& f,
const std::string&
name)
144 faceVarVectorDataInfo_.push_back(FaceVarVectorDataInfo{f,
name});
150 void write(
double time, Dune::VTK::OutputType type = Dune::VTK::ascii)
154 getFaceDataAndWrite_(time);
161 void updateCoordinates_()
163 coordinates_.resize(this->
gridGeometry().numFaceDofs());
164 for(
auto&& facet : facets(this->
gridGeometry().gridView()))
166 const int dofIdxGlobal = this->
gridGeometry().gridView().indexSet().index(facet);
167 coordinates_[dofIdxGlobal] = facet.geometry().center();
169 coordinatesInitialized_ =
true;
174 void getFaceDataAndWrite_(
const Scalar time)
176 const auto numPoints = this->
gridGeometry().numFaceDofs();
179 std::vector<bool> dofVisited(numPoints,
false);
182 if(!coordinatesInitialized_)
183 updateCoordinates_();
186 std::vector<std::vector<Scalar>> faceVarScalarData;
187 std::vector<std::vector<GlobalPosition>> faceVarVectorData;
189 if(!faceVarScalarDataInfo_.empty())
190 faceVarScalarData.resize(faceVarScalarDataInfo_.size(), std::vector<Scalar>(numPoints));
192 if(!faceVarVectorDataInfo_.empty())
193 faceVarVectorData.resize(faceVarVectorDataInfo_.size(), std::vector<GlobalPosition>(numPoints));
195 for (
const auto& element : elements(this->
gridGeometry().gridView(), Dune::Partitions::interior))
200 if (!faceVarScalarDataInfo_.empty() || !faceVarVectorDataInfo_.empty())
202 fvGeometry.bind(element);
203 elemFaceVars.bindElement(element, fvGeometry, this->
sol());
205 for (
auto&& scvf : scvfs(fvGeometry))
207 const auto dofIdxGlobal = scvf.dofIndex();
208 if(dofVisited[dofIdxGlobal])
211 dofVisited[dofIdxGlobal] =
true;
213 const auto& faceVars = elemFaceVars[scvf];
216 for (std::size_t i = 0; i < faceVarScalarDataInfo_.size(); ++i)
217 faceVarScalarData[i][dofIdxGlobal] = faceVarScalarDataInfo_[i].get(faceVars);
220 for (std::size_t i = 0; i < faceVarVectorDataInfo_.size(); ++i)
221 faceVarVectorData[i][dofIdxGlobal] = faceVarVectorDataInfo_[i].get(scvf, faceVars);
227 if(!faceVarScalarDataInfo_.empty())
228 for (std::size_t i = 0; i < faceVarScalarDataInfo_.size(); ++i)
229 faceWriter_->addPointData(faceVarScalarData[i], faceVarScalarDataInfo_[i].name);
231 if(!faceVarVectorDataInfo_.empty())
232 for (std::size_t i = 0; i < faceVarVectorDataInfo_.size(); ++i)
233 faceWriter_->addPointData(faceVarVectorData[i], faceVarVectorDataInfo_[i].name);
236 for(
const auto& field: faceFieldScalarDataInfo_)
237 faceWriter_->addPointData(field.data, field.name);
239 for(
const auto& field: faceFieldVectorDataInfo_)
240 faceWriter_->addPointData(field.data, field.name);
243 sequenceWriter_.write(time);
246 coordinates_.clear();
247 coordinates_.shrink_to_fit();
248 coordinatesInitialized_ =
false;
252 std::shared_ptr<PointCloudVtkWriter<Scalar, GlobalPosition>> faceWriter_;
254 VTKSequenceWriter<PointCloudVtkWriter<Scalar, GlobalPosition>> sequenceWriter_;
258 std::vector<GlobalPosition> coordinates_;
259 bool coordinatesInitialized_;
261 std::vector<FaceVarScalarDataInfo> faceVarScalarDataInfo_;
262 std::vector<FaceVarVectorDataInfo> faceVarVectorDataInfo_;
264 std::vector<FaceFieldScalarDataInfo> faceFieldScalarDataInfo_;
265 std::vector<FaceFieldVectorDataInfo> faceFieldVectorDataInfo_;
Base class to write pvd-files which contains a list of all collected vtk-files. This is a modified ve...
A VTK writer specialized for staggered grid implementations with dofs on the faces.
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:38
Velocity output for staggered free-flow models.
Definition: discretization/staggered/freeflow/velocityoutput.hh:38
A VTK output module to simplify writing dumux simulation data to VTK format Specialization for stagge...
Definition: staggeredvtkoutputmodule.hh:50
void addFaceField(const std::vector< Scalar > &v, const std::string &name)
Definition: staggeredvtkoutputmodule.hh:112
void write(double time, Dune::VTK::OutputType type=Dune::VTK::ascii)
Definition: staggeredvtkoutputmodule.hh:150
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:82
void addFaceField(const std::vector< GlobalPosition > &v, const std::string &name)
Definition: staggeredvtkoutputmodule.hh:123
void addFaceVariable(std::function< GlobalPosition(const SubControlVolumeFace &scvf, const FaceVariables &)> &&f, const std::string &name)
Definition: staggeredvtkoutputmodule.hh:142
void addFaceVariable(std::function< Scalar(const FaceVariables &)> &&f, const std::string &name)
Definition: staggeredvtkoutputmodule.hh:134
A VTK output module to simplify writing dumux simulation data to VTK format.
Definition: pointcloudvtkwriter.hh:51
const std::string & paramGroup() const
the parameter group for getting parameter from the parameter tree
Definition: io/vtkoutputmodule.hh:94
void write(double time, Dune::VTK::OutputType type=Dune::VTK::ascii)
Definition: io/vtkoutputmodule.hh:171
const std::string & name() const
Definition: io/vtkoutputmodule.hh:193
bool verbose() const
Definition: io/vtkoutputmodule.hh:192
A VTK output module to simplify writing dumux simulation data to VTK format.
Definition: io/vtkoutputmodule.hh:305
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
A VTK output module to simplify writing dumux simulation data to VTK format.