DUNE-FEM (unstable)
vtkwriter.hh
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright © DUNE Project contributors, see file LICENSE.md in module root
62 struct IsBindable< F, E, std::void_t< decltype( std::declval< F & >().bind( std::declval< const E & >() ) ),
313 // That is, a function that you can create a LocalFunction for, and evaluate that in element coordinates
314 template<typename F, std::enable_if_t<not Impl::IsBindable<F, Entity>::value && Impl::HasLocalFunction<F>::value, int> = 0>
324 template<typename F, std::enable_if_t<not Impl::IsBindable<F, Entity>::value && not Impl::HasLocalFunction<F>::value, int> = 0>
335 (vtkFunctionPtr->ncomps() == 2 || vtkFunctionPtr->ncomps() == 3) ? VTK::FieldInfo::Type::vector : VTK::FieldInfo::Type::scalar,
835 std::string pwrite ( const std::string & name, const std::string & path, const std::string & extendpath,
838 return pwrite( name, path, extendpath, type, gridView_.comm().rank(), gridView_.comm().size() );
1295 void writeData(VTK::VTUWriter& writer, const Data& data, const Iterator begin, const Iterator end, int nentries)
1313 DUNE_THROW(IOError,"Cannot write VTK vectors with more than 3 components (components was " << writecomps << ")");
Base class for stl conformant forward iterators.
Definition: iteratorfacades.hh:142
IndexType index(const typename Traits::template Codim< cc >::Entity &e) const
Map entity to index. The result of calling this method with an entity that is not in the index set is...
Definition: indexidset.hh:113
Implementation class for a multiple codim and multiple geometry type mapper.
Definition: mcmgmapper.hh:129
size_type size() const
Return total number of entities in the entity set managed by the mapper.
Definition: mcmgmapper.hh:204
Index subIndex(const typename GV::template Codim< 0 >::Entity &e, int i, unsigned int codim) const
Map subentity of codim 0 entity to starting index in array for dof block.
Definition: mcmgmapper.hh:185
Take a vector and interpret it as cell data for the VTKWriter.
Definition: function.hh:97
Take a vector and interpret it as point data for the VTKWriter.
Definition: function.hh:205
A base class for grid functions with any return type and dimension.
Definition: function.hh:42
Base class to write pvd-files which contains a list of all collected vtk-files.
Definition: vtksequencewriterbase.hh:34
Writer for the output of grid functions in the vtk format.
Definition: vtksequencewriter.hh:29
CellIterator(const GridCellIterator &x)
construct a CellIterator from the gridview's Iterator.
Definition: vtkwriter.hh:388
const FieldVector< DT, n > position() const
Definition: vtkwriter.hh:391
int id() const
Process-local consecutive zero-starting vertex id.
Definition: vtkwriter.hh:597
Type erasure wrapper for VTK data sets.
Definition: vtkwriter.hh:156
void unbind() const
Unbind the data set from the currently bound entity.
Definition: vtkwriter.hh:360
VTKLocalFunction(F &&f, VTK::FieldInfo fieldInfo)
Construct a VTKLocalFunction for a dune-functions style LocalFunction.
Definition: vtkwriter.hh:307
std::string name() const
Returns the name of the data set.
Definition: vtkwriter.hh:342
const VTK::FieldInfo & fieldInfo() const
Returns the VTK::FieldInfo for the data set.
Definition: vtkwriter.hh:348
void bind(const Entity &e) const
Bind the data set to grid entity e.
Definition: vtkwriter.hh:354
VTKLocalFunction(const std::shared_ptr< const VTKFunction > &vtkFunctionPtr)
Construct a VTKLocalFunction for a legacy VTKFunction.
Definition: vtkwriter.hh:331
void write(const Coordinate &pos, Writer &w) const
Write the value of the data set at local coordinate pos to the writer w.
Definition: vtkwriter.hh:366
FieldVector< DT, n > position() const
position of vertex inside the entity
Definition: vtkwriter.hh:501
int localindex() const
index of vertex within the entity, in Dune-numbering
Definition: vtkwriter.hh:496
Writer for the output of grid functions in the vtk format.
Definition: vtkwriter.hh:95
void addCellData(const Container &v, const std::string &name, int ncomps=1, VTK::Precision prec=VTK::Precision::float32)
Add a grid function (represented by container) that lives on the cells of the grid to the visualizati...
Definition: vtkwriter.hh:695
std::string write(const std::string &name, VTK::OutputType type=VTK::ascii)
write output (interface might change later)
Definition: vtkwriter.hh:803
std::string getParallelHeaderName(const std::string &name, const std::string &path, int commSize) const
return name of a parallel header file
Definition: vtkwriter.hh:920
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
std::string getSerialPieceName(const std::string &name, const std::string &path) const
return name of a serial piece file
Definition: vtkwriter.hh:940
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
void addVertexData(F &&f, VTK::FieldInfo vtkFieldInfo)
Add a function by sampling it on the grid vertices.
Definition: vtkwriter.hh:738
virtual void writeCellData(VTK::VTUWriter &writer)
write cell data
Definition: vtkwriter.hh:1337
virtual void countEntities(int &nvertices_, int &ncells_, int &ncorners_)
count the vertices, cells and corners
Definition: vtkwriter.hh:1239
std::string getParallelPieceName(const std::string &name, const std::string &path, int commRank, int commSize) const
return name of a parallel piece file (or header name)
Definition: vtkwriter.hh:855
virtual void writeGridCells(VTK::VTUWriter &writer)
write the connectivity array
Definition: vtkwriter.hh:1389
virtual void writeCellFaces(VTK::VTUWriter &writer)
write the connectivity array
Definition: vtkwriter.hh:1459
std::string write(const std::string &name, VTK::OutputType type, const int commRank, const int commSize)
write output (interface might change later)
Definition: vtkwriter.hh:965
VTK::Precision coordPrecision() const
get the precision with which coordinates are written out
Definition: vtkwriter.hh:782
virtual void writeGridPoints(VTK::VTUWriter &writer)
write the positions of vertices
Definition: vtkwriter.hh:1365
virtual void writeVertexData(VTK::VTUWriter &writer)
write vertex data
Definition: vtkwriter.hh:1351
void addCellData(F &&f, VTK::FieldInfo vtkFieldInfo)
Add a function by sampling it on the element centers.
Definition: vtkwriter.hh:674
void addVertexData(const Container &v, const std::string &name, int ncomps=1, VTK::Precision prec=VTK::Precision::float32)
Add a grid function (represented by container) that lives on the vertices of the grid to the visualiz...
Definition: vtkwriter.hh:760
VTKWriter(const GridView &gridView, VTK::DataMode dm=VTK::conforming, VTK::Precision coordPrecision=VTK::Precision::float32)
Construct a VTKWriter working on a specific GridView.
Definition: vtkwriter.hh:636
std::string pwrite(const std::string &name, const std::string &path, const std::string &extendpath, VTK::OutputType ot, const int commRank, const int commSize)
write output; interface might change later
Definition: vtkwriter.hh:1043
std::string pwrite(const std::string &name, const std::string &path, const std::string &extendpath, VTK::OutputType type=VTK::ascii)
write output (interface might change later)
Definition: vtkwriter.hh:835
void write(T data)
write one element of data
Definition: dataarraywriter.hh:69
@ vector
vector-valued field (always 3D, will be padded if necessary)
std::string name() const
The name of the data field.
Definition: common.hh:352
DataArrayWriter * makeArrayWriter(const std::string &name, unsigned ncomps, unsigned nitems, Precision prec)
acquire a DataArrayWriter
Definition: vtuwriter.hh:380
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
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
Data array writers for the VTKWriter.
A few common exception classes.
Common stuff for the VTKWriter.
Precision
which precision to use when writing out data to vtk files
Definition: common.hh:271
OutputType
How the bulk data should be stored in the file.
Definition: common.hh:43
DataMode
Whether to produce conforming or non-conforming output.
Definition: common.hh:67
Functions for VTK output.
typename Impl::voider< Types... >::type void_t
Is void for all valid input types. The workhorse for C++11 SFINAE-techniques.
Definition: typetraits.hh:40
const IndexSet & indexSet() const
obtain the index set
Definition: gridview.hh:177
const Communication & comm() const
obtain communication object
Definition: gridview.hh:273
static constexpr int dimension
The dimension of the grid.
Definition: gridview.hh:134
static constexpr int dimensionworld
The dimension of the world the grid lives in.
Definition: gridview.hh:137
PartitionIteratorType
Parameter to be used for the parallel level- and leaf iterators.
Definition: gridenums.hh:136
PartitionType
Attributes used in the generic overlap model.
Definition: gridenums.hh:30
@ InteriorBorder_Partition
interior and border entities
Definition: gridenums.hh:138
MCMGLayout mcmgVertexLayout()
layout for vertices (dim-0 entities)
Definition: mcmgmapper.hh:107
std::string relativePath(const std::string &newbase, const std::string &p)
compute a relative path between two paths
Definition: path.cc:149
std::string concatPaths(const std::string &base, const std::string &p)
concatenate two paths
Definition: path.cc:28
Utility class for handling nested indentation in output.
This file implements iterator facade classes for writing stl conformant iterators.
Mapper for multiple codim and multiple geometry types.
constexpr std::integral_constant< std::size_t, sizeof...(II)> size(std::integer_sequence< T, II... >)
Return the size of the sequence.
Definition: integersequence.hh:75
STL namespace.
Utilities for handling filesystem paths.
static const ReferenceElement & general(const GeometryType &type)
get general reference elements
Definition: referenceelements.hh:156
Type trait to determine whether an instance of T has an operator[](I), i.e. whether it can be indexed...
Definition: typetraits.hh:250
Base class for polymorphic container of underlying data set.
Definition: vtkwriter.hh:164
virtual void write(const Coordinate &pos, Writer &w, std::size_t count) const =0
Evaluate data set at local position pos inside the current entity and write result to w.
virtual void unbind()=0
Unbind data set from current grid entity - mostly here for performance and symmetry reasons.
virtual void bind(const Entity &e)=0
Bind data set to grid entity - must be called before evaluating (i.e. calling write())
Type erasure implementation for functions conforming to the dune-functions LocalFunction interface.
Definition: vtkwriter.hh:188
virtual void write(const Coordinate &pos, Writer &w, std::size_t count) const
Evaluate data set at local position pos inside the current entity and write result to w.
Definition: vtkwriter.hh:206
virtual void unbind()
Unbind data set from current grid entity - mostly here for performance and symmetry reasons.
Definition: vtkwriter.hh:201
virtual void bind(const Entity &e)
Bind data set to grid entity - must be called before evaluating (i.e. calling write())
Definition: vtkwriter.hh:196
Type erasure implementation for C++ functions, i.e., functions that can be evaluated in global coordi...
Definition: vtkwriter.hh:236
virtual void unbind()
Unbind data set from current grid entity - mostly here for performance and symmetry reasons.
Definition: vtkwriter.hh:250
virtual void write(const Coordinate &pos, Writer &w, std::size_t count) const
Evaluate data set at local position pos inside the current entity and write result to w.
Definition: vtkwriter.hh:255
virtual void bind(const Entity &e)
Bind data set to grid entity - must be called before evaluating (i.e. calling write())
Definition: vtkwriter.hh:245
Type erasure implementation for legacy VTKFunctions.
Definition: vtkwriter.hh:276
virtual void unbind()
Unbind data set from current grid entity - mostly here for performance and symmetry reasons.
Definition: vtkwriter.hh:287
virtual void write(const Coordinate &pos, Writer &w, std::size_t count) const
Evaluate data set at local position pos inside the current entity and write result to w.
Definition: vtkwriter.hh:292
virtual void bind(const Entity &e)
Bind data set to grid entity - must be called before evaluating (i.e. calling write())
Definition: vtkwriter.hh:282
Traits for type conversions and type information.
Definition of macros controlling symbol visibility at the ABI level.
#define DUNE_PRIVATE
Mark a symbol as being for internal use within the current DSO only.
Definition: visibility.hh:28
|
Legal Statements / Impressum |
Hosted by TU Dresden |
generated with Hugo v0.111.3
(Nov 12, 23:30, 2024)