Dune Core Modules (2.9.0)

subsamplingvtkwriter.hh
Go to the documentation of this file.
1 // SPDX-FileCopyrightText: Copyright (C) DUNE Project contributors, see file LICENSE.md in module root
2 // SPDX-License-Identifier: LicenseRef-GPL-2.0-only-with-DUNE-exception
3 // -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
4 // vi: set et ts=4 sw=2 sts=2:
5 
6 #ifndef DUNE_SUBSAMPLINGVTKWRITER_HH
7 #define DUNE_SUBSAMPLINGVTKWRITER_HH
8 
9 #include <ostream>
10 #include <memory>
11 
12 #include <dune/common/indent.hh>
13 #include <dune/geometry/type.hh>
16 #include <dune/grid/io/file/vtk/vtuwriter.hh>
17 
24 namespace Dune
25 {
37  template< class GridView >
39  : public VTKWriter<GridView>
40  {
41  typedef VTKWriter<GridView> Base;
42  constexpr static int dim = GridView::dimension;
43  constexpr static int dimw = GridView::dimensionworld;
44  typedef typename GridView::Grid::ctype ctype;
45  typedef typename GridView::template Codim< 0 >::Entity Entity;
47  typedef typename Refinement::IndexVector IndexVector;
48  typedef typename Refinement::ElementIterator SubElementIterator;
49  typedef typename Refinement::VertexIterator SubVertexIterator;
50 
51  typedef typename Base::CellIterator CellIterator;
52  typedef typename Base::FunctionIterator FunctionIterator;
53  using Base::cellBegin;
54  using Base::cellEnd;
55  using Base::celldata;
56  using Base::ncells;
57  using Base::ncorners;
58  using Base::nvertices;
59  using Base::outputtype;
60  using Base::vertexBegin;
61  using Base::vertexEnd;
62  using Base::vertexdata;
63 
64  public:
80  explicit SubsamplingVTKWriter (const GridView &gridView,
81  Dune::RefinementIntervals intervals_, bool coerceToSimplex_ = false,
82  VTK::Precision coordPrecision = VTK::Precision::float32)
83  : Base(gridView, VTK::nonconforming, coordPrecision)
84  , intervals(intervals_), coerceToSimplex(coerceToSimplex_)
85  {
86  if(intervals_.intervals() < 1) {
87  DUNE_THROW(Dune::IOError,"SubsamplingVTKWriter: Refinement intervals must be larger than zero! (One interval means no subsampling)");
88  }
89  }
90 
91  private:
92  GeometryType subsampledGeometryType(GeometryType geometryType)
93  {
94  return (geometryType.isCube() && !coerceToSimplex ? geometryType : GeometryTypes::simplex(dim));
95  }
96 
97  template<typename SubIterator>
98  struct IteratorSelector
99  {};
100 
101  SubElementIterator refinementBegin(const Refinement& refinement, Dune::RefinementIntervals intervals, IteratorSelector<SubElementIterator>)
102  {
103  return refinement.eBegin(intervals);
104  }
105 
106  SubVertexIterator refinementBegin(const Refinement& refinement, Dune::RefinementIntervals intervals, IteratorSelector<SubVertexIterator>)
107  {
108  return refinement.vBegin(intervals);
109  }
110 
111  SubElementIterator refinementEnd(const Refinement& refinement, Dune::RefinementIntervals intervals, IteratorSelector<SubElementIterator>)
112  {
113  return refinement.eEnd(intervals);
114  }
115 
116  SubVertexIterator refinementEnd(const Refinement& refinement, Dune::RefinementIntervals intervals, IteratorSelector<SubVertexIterator>)
117  {
118  return refinement.vEnd(intervals);
119  }
120 
121  template<typename Data, typename Iterator, typename SubIterator>
122  void writeData(VTK::VTUWriter& writer, const Data& data, const Iterator begin, const Iterator end, int nentries, IteratorSelector<SubIterator> sis)
123  {
124  for (auto it = data.begin(),
125  iend = data.end();
126  it != iend;
127  ++it)
128  {
129  const auto& f = *it;
130  VTK::FieldInfo fieldInfo = f.fieldInfo();
131  std::size_t writecomps = fieldInfo.size();
132  switch (fieldInfo.type())
133  {
135  break;
137  // vtk file format: a vector data always should have 3 comps (with
138  // 3rd comp = 0 in 2D case)
139  if (writecomps > 3)
140  DUNE_THROW(IOError,"Cannot write VTK vectors with more than 3 components (components was " << writecomps << ")");
141  writecomps = 3;
142  break;
144  DUNE_THROW(NotImplemented,"VTK output for tensors not implemented yet");
145  }
146  std::shared_ptr<VTK::DataArrayWriter> p
147  (writer.makeArrayWriter(f.name(), writecomps, nentries, fieldInfo.precision()));
148  if(!p->writeIsNoop())
149  for (Iterator eit = begin; eit!=end; ++eit)
150  {
151  const Entity & e = *eit;
152  f.bind(e);
153  Refinement &refinement =
154  buildRefinement<dim, ctype>(eit->type(),
155  subsampledGeometryType(eit->type()));
156  for(SubIterator sit = refinementBegin(refinement,intervals,sis),
157  send = refinementEnd(refinement,intervals,sis);
158  sit != send;
159  ++sit)
160  {
161  f.write(sit.coords(),*p);
162  // expand 2D-Vectors to 3D for VTK format
163  for(unsigned j = f.fieldInfo().size(); j < writecomps; j++)
164  p->write(0.0);
165  }
166  f.unbind();
167  }
168  }
169  }
170 
171 
172  protected:
174  virtual void countEntities(int &nvertices_, int &ncells_, int &ncorners_);
175 
177  virtual void writeCellData(VTK::VTUWriter& writer);
178 
180  virtual void writeVertexData(VTK::VTUWriter& writer);
181 
183  virtual void writeGridPoints(VTK::VTUWriter& writer);
184 
186  virtual void writeGridCells(VTK::VTUWriter& writer);
187 
188  public:
189  using Base::addVertexData;
190  using Base::addCellData;
191 
192  private:
193  // hide addVertexData -- adding raw data directly without a VTKFunction
194  // currently does not make sense for subsampled meshes, as the higher order
195  // information is missing. See FS#676.
196  template<class V>
197  void addVertexData (const V& v, const std::string &name, int ncomps=1);
198  template<class V>
199  void addCellData (const V& v, const std::string &name, int ncomps=1);
200 
201  Dune::RefinementIntervals intervals;
202  bool coerceToSimplex;
203  };
204 
206  template <class GridView>
207  void SubsamplingVTKWriter<GridView>::countEntities(int &nvertices_, int &ncells_, int &ncorners_)
208  {
209  nvertices_ = 0;
210  ncells_ = 0;
211  ncorners_ = 0;
212  for (CellIterator it=this->cellBegin(); it!=cellEnd(); ++it)
213  {
214  Refinement &refinement = buildRefinement<dim, ctype>(it->type(), subsampledGeometryType(it->type()));
215 
216  ncells_ += refinement.nElements(intervals);
217  nvertices_ += refinement.nVertices(intervals);
218  ncorners_ += refinement.nElements(intervals) * refinement.eBegin(intervals).vertexIndices().size();
219  }
220  }
221 
222 
224  template <class GridView>
226  {
227  if(celldata.size() == 0)
228  return;
229 
230  // Find the names of the first scalar and vector data fields.
231  // These will be marked as the default fields (the ones that ParaView shows
232  // when the file has just been opened).
233  std::string defaultScalarField, defaultVectorField;
234  std::tie(defaultScalarField, defaultVectorField) = this->getDataNames(celldata);
235 
236  writer.beginCellData(defaultScalarField, defaultVectorField);
237  writeData(writer,celldata,cellBegin(),cellEnd(),ncells,IteratorSelector<SubElementIterator>());
238  writer.endCellData();
239  }
240 
242  template <class GridView>
244  {
245  if(vertexdata.size() == 0)
246  return;
247 
248  // Find the names of the first scalar and vector data fields.
249  // These will be marked as the default fields (the ones that ParaView shows
250  // when the file has just been opened).
251  std::string defaultScalarField, defaultVectorField;
252  std::tie(defaultScalarField, defaultVectorField) = this->getDataNames(vertexdata);
253 
254  writer.beginPointData(defaultScalarField, defaultVectorField);
255  writeData(writer,vertexdata,cellBegin(),cellEnd(),nvertices,IteratorSelector<SubVertexIterator>());
256  writer.endPointData();
257  }
258 
260  template <class GridView>
262  {
263  writer.beginPoints();
264 
265  std::shared_ptr<VTK::DataArrayWriter> p
266  (writer.makeArrayWriter("Coordinates", 3, nvertices, this->coordPrecision()));
267  if(!p->writeIsNoop())
268  for (CellIterator i=cellBegin(); i!=cellEnd(); ++i)
269  {
270  Refinement &refinement =
271  buildRefinement<dim, ctype>(i->type(),
272  subsampledGeometryType(i->type()));
273  for(SubVertexIterator sit = refinement.vBegin(intervals),
274  send = refinement.vEnd(intervals);
275  sit != send; ++sit)
276  {
277  FieldVector<ctype, dimw> coords = i->geometry().global(sit.coords());
278  for (int j=0; j<std::min(int(dimw),3); j++)
279  p->write(coords[j]);
280  for (int j=std::min(int(dimw),3); j<3; j++)
281  p->write(0.0);
282  }
283  }
284  // free the VTK::DataArrayWriter before touching the stream
285  p.reset();
286 
287  writer.endPoints();
288  }
289 
291  template <class GridView>
293  {
294  writer.beginCells();
295 
296  // connectivity
297  {
298  std::shared_ptr<VTK::DataArrayWriter> p1
299  (writer.makeArrayWriter("connectivity", 1, ncorners, VTK::Precision::int32));
300  // The offset within the index numbering
301  if(!p1->writeIsNoop()) {
302  int offset = 0;
303  for (CellIterator i=cellBegin(); i!=cellEnd(); ++i)
304  {
305  GeometryType coercedToType = subsampledGeometryType(i->type());
306  Refinement &refinement =
307  buildRefinement<dim, ctype>(i->type(), coercedToType);
308  for(SubElementIterator sit = refinement.eBegin(intervals),
309  send = refinement.eEnd(intervals);
310  sit != send; ++sit)
311  {
312  IndexVector indices = sit.vertexIndices();
313  for(unsigned int ii = 0; ii < indices.size(); ++ii)
314  p1->write(offset+indices[VTK::renumber(coercedToType, ii)]);
315  }
316  offset += refinement.nVertices(intervals);
317  }
318  }
319  }
320 
321  // offsets
322  {
323  std::shared_ptr<VTK::DataArrayWriter> p2
324  (writer.makeArrayWriter("offsets", 1, ncells, VTK::Precision::int32));
325  if(!p2->writeIsNoop()) {
326  // The offset into the connectivity array
327  int offset = 0;
328  for (CellIterator i=cellBegin(); i!=cellEnd(); ++i)
329  {
330  Refinement &refinement =
331  buildRefinement<dim, ctype>(i->type(),
332  subsampledGeometryType(i->type()));
333  unsigned int verticesPerCell =
334  refinement.eBegin(intervals).vertexIndices().size();
335  for(int element = 0; element < refinement.nElements(intervals);
336  ++element)
337  {
338  offset += verticesPerCell;
339  p2->write(offset);
340  }
341  }
342  }
343  }
344 
345  // types
346  if (dim>1)
347  {
348  std::shared_ptr<VTK::DataArrayWriter> p3
349  (writer.makeArrayWriter("types", 1, ncells, VTK::Precision::uint8));
350  if(!p3->writeIsNoop())
351  for (CellIterator it=cellBegin(); it!=cellEnd(); ++it)
352  {
353  GeometryType coerceTo = subsampledGeometryType(it->type());
354  Refinement &refinement =
355  buildRefinement<dim, ctype>(it->type(), coerceTo);
356  int vtktype = VTK::geometryType(coerceTo);
357  for(int i = 0; i < refinement.nElements(intervals); ++i)
358  p3->write(vtktype);
359  }
360  }
361 
362  writer.endCells();
363  }
364 }
365 
366 #endif // DUNE_SUBSAMPLINGVTKWRITER_HH
vector space out of a tensor product of fields.
Definition: fvector.hh:95
Unique label for each type of entities that can occur in DUNE grids.
Definition: type.hh:125
Grid view abstract base class.
Definition: gridview.hh:66
Default exception class for I/O errors.
Definition: exceptions.hh:231
Holds the number of refined intervals per axis needed for virtual and static refinement.
Definition: base.cc:94
Writer for the output of subsampled grid functions in the vtk format.
Definition: subsamplingvtkwriter.hh:40
void addVertexData(const std::shared_ptr< const VTKFunction > &p)
Add a grid function that lives on the vertices of the grid to the visualization.
Definition: vtkwriter.hh:713
void addCellData(const std::shared_ptr< const VTKFunction > &p)
Add a grid function that lives on the cells of the grid to the visualization.
Definition: vtkwriter.hh:649
virtual void writeGridPoints(VTK::VTUWriter &writer)
write the positions of vertices
Definition: subsamplingvtkwriter.hh:261
virtual void writeVertexData(VTK::VTUWriter &writer)
write vertex data
Definition: subsamplingvtkwriter.hh:243
SubsamplingVTKWriter(const GridView &gridView, Dune::RefinementIntervals intervals_, bool coerceToSimplex_=false, VTK::Precision coordPrecision=VTK::Precision::float32)
Construct a SubsamplingVTKWriter working on a specific GridView.
Definition: subsamplingvtkwriter.hh:80
virtual void countEntities(int &nvertices_, int &ncells_, int &ncorners_)
count the vertices, cells and corners
Definition: subsamplingvtkwriter.hh:207
virtual void writeCellData(VTK::VTUWriter &writer)
write cell data
Definition: subsamplingvtkwriter.hh:225
virtual void writeGridCells(VTK::VTUWriter &writer)
write the connectivity array
Definition: subsamplingvtkwriter.hh:292
Iterator over the grids elements.
Definition: vtkwriter.hh:385
Writer for the ouput of grid functions in the vtk format.
Definition: vtkwriter.hh:95
void addVertexData(const std::shared_ptr< const VTKFunction > &p)
Add a grid function that lives on the vertices of the grid to the visualization.
Definition: vtkwriter.hh:713
void addCellData(const std::shared_ptr< const VTKFunction > &p)
Add a grid function that lives on the cells of the grid to the visualization.
Definition: vtkwriter.hh:649
VTK::Precision coordPrecision() const
get the precision with which coordinates are written out
Definition: vtkwriter.hh:782
@ tensor
tensor field (always 3x3)
@ vector
vector-valued field (always 3D, will be padded if necessary)
Dump a .vtu/.vtp files contents to a stream.
Definition: vtuwriter.hh:98
void endCellData()
finish CellData section
Definition: vtuwriter.hh:220
void beginCells()
start section for the grid cells/PolyData lines
Definition: vtuwriter.hh:274
void endPointData()
finish PointData section
Definition: vtuwriter.hh:182
void beginCellData(const std::string &scalars="", const std::string &vectors="")
start CellData section
Definition: vtuwriter.hh:205
void beginPointData(const std::string &scalars="", const std::string &vectors="")
start PointData section
Definition: vtuwriter.hh:167
DataArrayWriter * makeArrayWriter(const std::string &name, unsigned ncomps, unsigned nitems, Precision prec)
acquire a DataArrayWriter
Definition: vtuwriter.hh:380
void endPoints()
finish section for the point coordinates
Definition: vtuwriter.hh:249
void endCells()
start section for the grid cells/PolyData lines
Definition: vtuwriter.hh:285
void beginPoints()
start section for the point coordinates
Definition: vtuwriter.hh:238
VirtualRefinement base class.
Definition: virtualrefinement.hh:294
VertexIterator vBegin(Dune::RefinementIntervals tag) const
Get an VertexIterator.
Definition: virtualrefinement.cc:38
ElementIterator eBegin(Dune::RefinementIntervals tag) const
Get an ElementIterator.
Definition: virtualrefinement.cc:54
virtual int nElements(Dune::RefinementIntervals tag) const =0
Get the number of Elements.
ElementIterator eEnd(Dune::RefinementIntervals tag) const
Get an ElementIterator.
Definition: virtualrefinement.cc:62
Codim< 0 >::SubEntityIterator ElementIterator
The ElementIterator of the VirtualRefinement.
Definition: virtualrefinement.hh:301
VertexIterator vEnd(Dune::RefinementIntervals tag) const
Get an VertexIterator.
Definition: virtualrefinement.cc:46
virtual int nVertices(Dune::RefinementIntervals tag) const =0
Get the number of Vertices.
std::vector< int > IndexVector
The IndexVector of the VirtualRefinement.
Definition: virtualrefinement.hh:314
Codim< dimension >::SubEntityIterator VertexIterator
The VertexIterator of the VirtualRefinement.
Definition: virtualrefinement.hh:297
Precision
which precision to use when writing out data to vtk files
Definition: common.hh:271
#define DUNE_THROW(E, m)
Definition: exceptions.hh:218
constexpr static int dimension
The dimension of the grid.
Definition: gridview.hh:148
constexpr static int dimensionworld
The dimension of the world the grid lives in.
Definition: gridview.hh:151
constexpr GeometryType simplex(unsigned int dim)
Returns a GeometryType representing a simplex of dimension dim.
Definition: type.hh:463
auto min(ADLTag< 0 >, const V &v1, const V &v2)
implements binary Simd::min()
Definition: defaults.hh:89
Utility class for handling nested indentation in output.
Dune namespace.
Definition: alignedallocator.hh:13
@ nonconforming
Definition: declaration.hh:26
Static tag representing a codimension.
Definition: dimension.hh:24
A unique label for each type of element that can occur in a grid.
This file contains the virtual wrapper around refinement.
Provides file i/o for the visualization toolkit.
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.80.0 (Apr 29, 22:29, 2024)