4 #ifndef DUNE_VTKWRITER_HH 5 #define DUNE_VTKWRITER_HH 17 #include <dune/common/deprecated.hh> 18 #include <dune/common/typetraits.hh> 19 #include <dune/common/exceptions.hh> 20 #include <dune/common/std/memory.hh> 21 #include <dune/common/indent.hh> 22 #include <dune/common/iteratorfacades.hh> 23 #include <dune/common/path.hh> 24 #include <dune/common/shared_ptr.hh> 25 #include <dune/geometry/referenceelements.hh> 53 template<
typename F,
typename =
int>
55 :
public std::false_type
59 struct _has_local_context<T,typename std::enable_if<(sizeof(static_cast<T*>(nullptr)->localContext()) > 0),int>::type>
60 :
public std::true_type
65 namespace VTKWriteTypeTraits {
73 template <
class Gr
idView>
75 template <
class Gr
idView>
86 template<
class Gr
idView >
101 typedef typename GridView::template Codim< 0 >::Entity Cell;
102 typedef typename GridView::template Codim< n >::Entity Vertex;
110 typedef typename GridView::template Codim< 0 >
111 ::template Partition< VTK_Partition >::Iterator
113 typedef typename GridView::template Codim< n >
114 ::template Partition< VTK_Partition >::Iterator
117 typedef typename GridCellIterator::Reference EntityReference;
119 typedef typename GridView::template Codim< 0 >
120 ::Entity::Geometry::LocalCoordinate Coordinate;
127 switch( VTK_Partition )
132 default: DUNE_THROW(NotImplemented,
"Add check for this partition type");
160 virtual void bind(
const Entity& e) = 0;
163 virtual void unbind() = 0;
169 virtual void write(
const Coordinate& pos, Writer& w, std::size_t count)
const = 0;
183 template<
typename F_>
185 : _f(
std::forward<F_>(f))
188 virtual void bind(
const Entity& e)
198 virtual void write(
const Coordinate& pos, Writer& w, std::size_t count)
const 202 do_write(w,r,count,is_indexable<decltype(r)>());
208 void do_write(Writer& w,
const R& r, std::size_t count, std::true_type)
const 210 for (std::size_t i = 0; i < count; ++i)
215 void do_write(Writer& w,
const R& r, std::size_t count, std::false_type)
const 233 virtual void bind(
const Entity& e)
243 virtual void write(
const Coordinate& pos, Writer& w, std::size_t count)
const 245 for (std::size_t i = 0; i < count; ++i)
246 w.
write(_f->evaluate(i,*_entity,pos));
252 const Entity* _entity;
261 , _fieldInfo(fieldInfo)
269 typename
std::decay<decltype(localFunction(
std::forward<F>(f)))>::type
270 > >(localFunction(
std::forward<F>(f))))
271 , _fieldInfo(fieldInfo)
278 vtkFunctionPtr->name(),
279 vtkFunctionPtr->ncomps() > 1 ? VTK::FieldInfo::Type::vector : VTK::FieldInfo::Type::scalar,
280 vtkFunctionPtr->ncomps()
287 return fieldInfo().name();
297 void bind(
const Entity& e)
const 309 void write(
const Coordinate& pos, Writer& w)
const 311 _f->write(pos,w,fieldInfo().size());
314 std::shared_ptr<FunctionWrapperBase>
_f;
336 return ReferenceElements<DT,n>::general((*this)->type()).position(0,0);
342 return gridView_.template begin< 0, VTK_Partition >();
347 return gridView_.template end< 0, VTK_Partition >();
366 public ForwardIteratorFacade<VertexIterator, const Entity, EntityReference, int>
368 GridCellIterator git;
369 GridCellIterator gend;
374 const VertexMapper & vertexmapper;
375 std::vector<bool> visited;
388 const int numCorners = git->subEntities(n);
389 if( cornerIndexDune == numCorners )
391 offset += numCorners;
395 while( (git != gend) && skipEntity( git->partitionType() ) )
401 const GridCellIterator & end,
403 const VertexMapper & vm) :
404 git(x), gend(end), datamode(dm), cornerIndexDune(0),
405 vertexmapper(vm), visited(vm.size(), false),
409 visited[vertexmapper.
subIndex(*git,cornerIndexDune,n)] =
true;
416 while(visited[vertexmapper.
subIndex(*git,cornerIndexDune,n)])
419 if (git == gend)
return;
421 visited[vertexmapper.
subIndex(*git,cornerIndexDune,n)] =
true;
430 return git == cit.git
431 && cornerIndexDune == cit.cornerIndexDune
432 && datamode == cit.datamode;
441 return cornerIndexDune;
446 return ReferenceElements<DT,n>::general(git->type())
447 .position(cornerIndexDune,n);
453 return VertexIterator( gridView_.template begin< 0, VTK_Partition >(),
454 gridView_.template end< 0, VTK_Partition >(),
455 datamode, *vertexmapper );
460 return VertexIterator( gridView_.template end< 0, VTK_Partition >(),
461 gridView_.template end< 0, VTK_Partition >(),
462 datamode, *vertexmapper );
481 public ForwardIteratorFacade<CornerIterator, const Entity, EntityReference, int>
483 GridCellIterator git;
484 GridCellIterator gend;
489 const VertexMapper & vertexmapper;
493 const std::vector<int> & number;
502 const GridCellIterator & end,
504 const VertexMapper & vm,
505 const std::vector<int> & num) :
506 git(x), gend(end), datamode(dm), cornerIndexVTK(0),
508 number(num), offset(0) {}
514 const int numCorners = git->subEntities(n);
515 if( cornerIndexVTK == numCorners )
517 offset += numCorners;
521 while( (git != gend) && skipEntity( git->partitionType() ) )
527 return git == cit.git
528 && cornerIndexVTK == cit.cornerIndexVTK
529 && datamode == cit.datamode;
551 DUNE_THROW(IOError,
"VTKWriter: unsupported DataMode" << datamode);
558 return CornerIterator( gridView_.template begin< 0, VTK_Partition >(),
559 gridView_.template end< 0, VTK_Partition >(),
560 datamode, *vertexmapper, number );
565 return CornerIterator( gridView_.template end< 0, VTK_Partition >(),
566 gridView_.template end< 0, VTK_Partition >(),
567 datamode, *vertexmapper, number );
580 : gridView_( gridView ),
604 void addCellData (VTKFunction* p) DUNE_DEPRECATED_MSG(
"Don't pass raw pointers, use the version with shared_ptr")
625 void addCellData (
const V& v,
const std::string &name,
int ncomps = 1)
628 for (
int c=0; c<ncomps; ++c) {
629 std::stringstream compName;
632 compName <<
"[" << c <<
"]";
633 VTKFunction* p =
new Function(gridView_, v, compName.str(), ncomps, c);
634 addCellData(VTKFunctionPtr(p));
643 void addVertexData (VTKFunction* p) DUNE_DEPRECATED_MSG(
"Don't pass raw pointers, use the version with shared_ptr")
683 for (
int c=0; c<ncomps; ++c) {
684 std::stringstream compName;
687 compName <<
"[" << c <<
"]";
688 VTKFunction* p =
new Function(gridView_, v, compName.str(), ncomps, c);
689 addVertexData(VTKFunctionPtr(p));
717 std::string
write (
const std::string &name,
720 return write( name, type, gridView_.comm().rank(), gridView_.comm().size() );
749 std::string
pwrite (
const std::string & name,
const std::string & path,
const std::string & extendpath,
752 return pwrite( name, path, extendpath, type, gridView_.comm().rank(), gridView_.comm().size() );
769 const std::string& path,
770 int commRank,
int commSize)
const 772 std::ostringstream s;
773 if(path.size() > 0) {
775 if(path[path.size()-1] !=
'/')
778 s <<
's' << std::setw(4) << std::setfill(
'0') << commSize <<
'-';
779 s <<
'p' << std::setw(4) << std::setfill(
'0') << commRank <<
'-';
800 const std::string& path,
803 std::ostringstream s;
804 if(path.size() > 0) {
806 if(path[path.size()-1] !=
'/')
809 s <<
's' << std::setw(4) << std::setfill(
'0') << commSize <<
'-';
832 const std::string& path)
const 834 static const std::string extension =
837 return concatPaths(path, name+extension);
855 std::string
write (
const std::string &name,
863 return pwrite(name,
"",
"", type, commRank, commSize);
869 std::string pieceName = getSerialPieceName(name,
"");
873 file.exceptions(std::ios_base::badbit | std::ios_base::failbit |
874 std::ios_base::eofbit);
877 file.open( pieceName.c_str(), std::ios::binary );
880 std::cerr <<
"Filename: " << pieceName <<
" could not be opened" << std::endl;
883 if (! file.is_open())
884 DUNE_THROW(IOError,
"Could not write to piece file " << pieceName);
885 writeDataFile( file );
915 std::string
pwrite(
const std::string& name,
const std::string& path,
916 const std::string& extendpath,
925 file.exceptions(std::ios_base::badbit | std::ios_base::failbit |
926 std::ios_base::eofbit);
927 std::string piecepath = concatPaths(path, extendpath);
928 std::string relpiecepath = relativePath(path, piecepath);
931 std::string fullname = getParallelPieceName(name, piecepath, commRank,
935 file.open(fullname.c_str(),std::ios::binary);
938 std::cerr <<
"Filename: " << fullname <<
" could not be opened" << std::endl;
941 if (! file.is_open())
942 DUNE_THROW(IOError,
"Could not write to piecefile file " << fullname);
945 gridView_.comm().barrier();
948 fullname = getParallelHeaderName(name, path, commSize);
951 file.open(fullname.c_str());
952 if (! file.is_open())
953 DUNE_THROW(IOError,
"Could not write to parallel file " << fullname);
954 writeParallelHeader(file,name,relpiecepath, commSize );
957 gridView_.comm().barrier();
980 void writeParallelHeader(std::ostream& s,
const std::string& piecename,
981 const std::string& piecepath,
const int commSize)
992 std::string scalars, vectors;
993 std::tie(scalars,vectors) = getDataNames(vertexdata);
996 for (
auto it = vertexdata.begin(),
997 end = vertexdata.end();
1001 unsigned writecomps = it->fieldInfo().size();
1002 if(writecomps == 2) writecomps = 3;
1003 writer.
addArray<
float>(it->name(), writecomps);
1009 std::string scalars, vectors;
1010 std::tie(scalars,vectors) = getDataNames(celldata);
1013 for (
auto it = celldata.begin(),
1014 end = celldata.end();
1018 unsigned writecomps = it->fieldInfo().size();
1019 if(writecomps == 2) writecomps = 3;
1020 writer.
addArray<
float>(it->name(), writecomps);
1026 writer.
addArray<
float>(
"Coordinates", 3);
1030 for(
int i = 0; i < commSize; ++i )
1032 const std::string& fullname = getParallelPieceName(piecename,
1042 void writeDataFile (std::ostream& s)
1050 vertexmapper =
new VertexMapper( gridView_ );
1053 number.resize(vertexmapper->size());
1054 for (std::vector<int>::size_type i=0; i<number.size(); i++) number[i] = -1;
1056 countEntities(nvertices, ncells, ncorners);
1059 writeAllData(writer);
1064 writeAllData(writer);
1067 delete vertexmapper; number.clear();
1072 writeVertexData(writer);
1075 writeCellData(writer);
1078 writeGridPoints(writer);
1081 writeGridCells(writer);
1095 DUNE_THROW(IOError,
"VTKWriter: unsupported OutputType" << outputtype);
1103 return "UnstructuredGrid";
1117 const int subEntities = it->subEntities(n);
1118 for (
int i=0; i<subEntities; ++i)
1123 int alpha = vertexmapper->subIndex(*it,i,n);
1124 if (number[alpha]<0)
1125 number[alpha] = nvertices++;
1135 template<
typename T>
1138 std::string scalars =
"";
1139 for (
auto it = data.begin(),
1145 scalars = it->name();
1149 std::string vectors =
"";
1150 for (
auto it = data.begin(),
1156 vectors = it->name();
1159 return std::make_tuple(scalars,vectors);
1162 template<
typename Data,
typename Iterator>
1165 for (
auto it = data.begin(),
1170 const auto& f = *it;
1172 std::size_t writecomps = fieldInfo.
size();
1173 switch (fieldInfo.
type())
1181 DUNE_THROW(IOError,
"Cannot write VTK vectors with more than 3 components (components was " << writecomps <<
")");
1185 DUNE_THROW(NotImplemented,
"VTK output for tensors not implemented yet");
1187 shared_ptr<VTK::DataArrayWriter<float> > p
1189 if(!p->writeIsNoop())
1190 for (Iterator eit = begin; eit!=end; ++eit)
1192 const Entity & e = *eit;
1194 f.write(eit.position(),*p);
1198 for (std::size_t j=fieldInfo.
size(); j < writecomps; ++j)
1207 if(celldata.size() == 0)
1210 std::string scalars, vectors;
1211 std::tie(scalars,vectors) = getDataNames(celldata);
1214 writeData(writer,celldata,cellBegin(),cellEnd(),ncells);
1221 if(vertexdata.size() == 0)
1224 std::string scalars, vectors;
1225 std::tie(scalars,vectors) = getDataNames(vertexdata);
1228 writeData(writer,vertexdata,vertexBegin(),vertexEnd(),nvertices);
1237 shared_ptr<VTK::DataArrayWriter<float> > p
1239 if(!p->writeIsNoop()) {
1244 for (
int j=0; j<
std::min(dimw,3); j++)
1245 p->write((*vit).geometry().corner(vit.localindex())[j]);
1246 for (
int j=std::min(dimw,3); j<3; j++)
1263 shared_ptr<VTK::DataArrayWriter<int> > p1
1265 if(!p1->writeIsNoop())
1272 shared_ptr<VTK::DataArrayWriter<int> > p2
1274 if(!p2->writeIsNoop()) {
1278 offset += it->subEntities(n);
1287 shared_ptr<VTK::DataArrayWriter<unsigned char> > p3
1289 if(!p3->writeIsNoop())
1313 VertexMapper* vertexmapper;
1316 std::vector<int> number;
void beginPointData(const std::string &scalars="", const std::string &vectors="")
start PointData section
Definition: pvtuwriter.hh:118
std::string write(const std::string &name, VTK::OutputType type=VTK::ascii)
write output (interface might change later)
Definition: vtkwriter.hh:717
CellIterator cellBegin() const
Definition: vtkwriter.hh:340
std::string getFormatString() const
Definition: vtkwriter.hh:1085
virtual void writeGridCells(VTK::VTUWriter &writer)
write the connectivity array
Definition: vtkwriter.hh:1257
void endAppended()
finish the appended data section
Definition: vtuwriter.hh:357
void writeData(VTK::VTUWriter &writer, const Data &data, const Iterator begin, const Iterator end, int nentries)
Definition: vtkwriter.hh:1163
virtual ~FunctionWrapperBase()
Definition: vtkwriter.hh:171
void increment()
Definition: vtkwriter.hh:509
void endCellData()
finish CellData section
Definition: pvtuwriter.hh:153
Base class for polymorphic container of underlying data set.
Definition: vtkwriter.hh:156
CellIterator(const GridCellIterator &x)
construct a CellIterator from the gridview's Iterator.
Definition: vtkwriter.hh:331
Type erasure implementation for legacy VTKFunctions.
Definition: vtkwriter.hh:225
void endMain()
finish the main PolyData/UnstructuredGrid section
Definition: pvtuwriter.hh:193
Dune::VTKFunction< GridView > VTKFunction
Definition: vtkwriter.hh:139
virtual ~VTKWriter()
destructor
Definition: vtkwriter.hh:701
Implementation class for a multiple codim and multiple geometry type mapper.
Definition: mcmgmapper.hh:102
Traits::Grid Grid
type of the grid
Definition: common/gridview.hh:77
const FieldVector< DT, n > position() const
Definition: vtkwriter.hh:334
void addCellData(const V &v, const std::string &name, int ncomps=1)
Add a grid function (represented by container) that lives on the cells of the grid to the visualizati...
Definition: vtkwriter.hh:625
VTK::OutputType outputtype
Definition: vtkwriter.hh:1319
Output to the file is in ascii.
Definition: common.hh:42
Functions for VTK output.
PartitionType
Attributes used in the generic overlap model.
Definition: gridenums.hh:28
EntityReference dereference() const
Definition: vtkwriter.hh:531
void beginCells()
start section for the grid cells/PolyData lines
Definition: vtuwriter.hh:272
VTKFunctionWrapper(const VTKFunctionPtr &f)
Definition: vtkwriter.hh:228
void unbind() const
Unbind the data set from the currently bound entity.
Definition: vtkwriter.hh:303
VTKLocalFunction(F &&f, VTK::FieldInfo fieldInfo, typename std::enable_if< detail::_has_local_context< F >::value, int >::type dummy=0)
Construct a VTKLocalFunction for a dune-functions style LocalFunction.
Definition: vtkwriter.hh:258
void endMain()
finish the main PolyData/UnstructuredGrid section
Definition: vtuwriter.hh:318
CellIterator cellEnd() const
Definition: vtkwriter.hh:345
Definition: vtkwriter.hh:67
virtual void bind(const Entity &e)
Bind data set to grid entity - must be called before evaluating (i.e. calling write()) ...
Definition: vtkwriter.hh:233
void beginPoints()
start section for the point coordinates
Definition: pvtuwriter.hh:164
void bind(const Entity &e) const
Bind the data set to grid entity e.
Definition: vtkwriter.hh:297
void beginCellData(const std::string &scalars="", const std::string &vectors="")
start CellData section
Definition: vtuwriter.hh:203
void beginPoints()
start section for the point coordinates
Definition: vtuwriter.hh:236
for .vtu files (UnstructuredGrid)
Definition: common.hh:294
VertexIterator vertexBegin() const
Definition: vtkwriter.hh:451
void beginMain(unsigned ncells, unsigned npoints)
start the main PolyData/UnstructuredGrid section
Definition: vtuwriter.hh:308
Descriptor struct for VTK fields.
Definition: common.hh:306
Ouput is to the file is appended raw binary.
Definition: common.hh:46
Writer for the ouput of grid functions in the vtk format.Writes arbitrary grid functions (living on c...
Definition: vtksequencewriter.hh:23
all entities
Definition: gridenums.hh:139
Definition: vtkwriter.hh:54
std::string getParallelPieceName(const std::string &name, const std::string &path, int commRank, int commSize) const
return name of a parallel piece file
Definition: vtkwriter.hh:768
void beginPointData(const std::string &scalars="", const std::string &vectors="")
start PointData section
Definition: vtuwriter.hh:165
shared_ptr< const VTKFunction > VTKFunctionPtr
Definition: vtkwriter.hh:140
virtual void bind(const Entity &e)
Bind data set to grid entity - must be called before evaluating (i.e. calling write()) ...
Definition: vtkwriter.hh:188
bool equals(const VertexIterator &cit) const
Definition: vtkwriter.hh:428
bool beginAppended()
start the appended data section
Definition: vtuwriter.hh:343
Take a vector and interpret it as point data for the VTKWriter.
Definition: function.hh:186
std::tuple< std::string, std::string > getDataNames(const T &data) const
Definition: vtkwriter.hh:1136
DataMode
Whether to produce conforming or non-conforming output.
Definition: common.hh:64
void increment()
Definition: vtkwriter.hh:411
void addVertexData(VTKFunction *p)
Add a grid function that lives on the vertices of the grid to the visualization.
Definition: vtkwriter.hh:643
std::string getSerialPieceName(const std::string &name, const std::string &path) const
return name of a serial piece file
Definition: vtkwriter.hh:831
OutputType
How the bulk data should be stored in the file.
Definition: common.hh:40
FunctionWrapper(F_ &&f)
Definition: vtkwriter.hh:184
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:198
VertexIterator vertexEnd() const
Definition: vtkwriter.hh:458
virtual void writeCellData(VTK::VTUWriter &writer)
write cell data
Definition: vtkwriter.hh:1205
int ncorners
Definition: vtkwriter.hh:1311
A base class for grid functions with any return type and dimension.
Definition: function.hh:38
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:309
all interior entities
Definition: gridenums.hh:29
int localindex() const
index of vertex within the entity, in Dune-numbering
Definition: vtkwriter.hh:439
void addVertexData(F &&f, VTK::FieldInfo vtkFieldInfo)
Definition: vtkwriter.hh:658
Iterate over the elements' corners.
Definition: vtkwriter.hh:480
vector-valued field (always 3D, will be padded if necessary)
void endCellData()
finish CellData section
Definition: vtuwriter.hh:218
Type type() const
The type of the data field.
Definition: common.hh:336
Type erasure wrapper for VTK data sets.
Definition: vtkwriter.hh:148
VertexIterator(const GridCellIterator &x, const GridCellIterator &end, const VTK::DataMode &dm, const VertexMapper &vm)
Definition: vtkwriter.hh:400
VTKLocalFunction(const VTKFunctionPtr &vtkFunctionPtr)
Construct a VTKLocalFunction for a legacy VTKFunction.
Definition: vtkwriter.hh:275
virtual void countEntities(int &nvertices, int &ncells, int &ncorners)
count the vertices, cells and corners
Definition: vtkwriter.hh:1107
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:915
interior and border entities
Definition: gridenums.hh:136
EntityReference dereference() const
Definition: vtkwriter.hh:434
void addArray(const std::string &name, unsigned ncomps)
Add an array to the output file.
Definition: pvtuwriter.hh:205
CornerIterator(const GridCellIterator &x, const GridCellIterator &end, const VTK::DataMode &dm, const VertexMapper &vm, const std::vector< int > &num)
Definition: vtkwriter.hh:501
Dump a .vtu/.vtp files contents to a stream.
Definition: vtuwriter.hh:96
void endPointData()
finish PointData section
Definition: vtuwriter.hh:180
Grid::ctype ctype
type used for coordinates in grid
Definition: common/gridview.hh:127
int nvertices
Definition: vtkwriter.hh:1310
CornerIterator cornerBegin() const
Definition: vtkwriter.hh:556
VTK::FieldInfo _fieldInfo
Definition: vtkwriter.hh:315
std::string name() const
Returns the name of the data set.
Definition: vtkwriter.hh:285
for .vtp files (PolyData)
Definition: common.hh:292
int ncells
Definition: vtkwriter.hh:1309
void addVertexData(const VTKFunctionPtr &p)
Add a grid function that lives on the vertices of the grid to the visualization.
Definition: vtkwriter.hh:652
Ouput is to the file is appended base64 binary.
Definition: common.hh:48
bool equals(const CornerIterator &cit) const
Definition: vtkwriter.hh:525
tensor field (always 3x3)
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:243
int id() const
Process-local consecutive zero-starting vertex id.
Definition: vtkwriter.hh:540
Grid view abstract base class.
Definition: common/gridview.hh:58
void beginCellData(const std::string &scalars="", const std::string &vectors="")
start CellData section
Definition: pvtuwriter.hh:144
void basicIncrement()
Definition: vtkwriter.hh:383
std::size_t size() const
The number of components in the data field.
Definition: common.hh:342
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:855
The dimension of the world the grid lives in.
Definition: common/gridview.hh:134
Output to the file is inline base64 binary.
Definition: common.hh:44
PartitionIteratorType
Parameter to be used for the parallel level- and leaf iterators.
Definition: gridenums.hh:134
int renumber(const Dune::GeometryType &t, int i)
renumber VTK <-> Dune
Definition: common.hh:224
VTKLocalFunction(F &&f, VTK::FieldInfo fieldInfo, typename std::enable_if< not detail::_has_local_context< F >::value, int >::type dummy=0)
Construct a VTKLocalFunction for a dune-functions style Function.
Definition: vtkwriter.hh:266
base class for data array writers
Definition: dataarraywriter.hh:53
std::decay< F >::type Function
Definition: vtkwriter.hh:181
Output conforming data.
Definition: common.hh:70
void endCells()
start section for the grid cells/PolyData lines
Definition: vtuwriter.hh:283
Include standard header files.
Definition: agrid.hh:59
VTKWriter(const GridView &gridView, VTK::DataMode dm=VTK::conforming)
Construct a VTKWriter working on a specific GridView.
Definition: vtkwriter.hh:578
DataArrayWriter< T > * makeArrayWriter(const std::string &name, unsigned ncomps, unsigned nitems)
aquire a DataArrayWriter
Definition: vtuwriter.hh:379
Mapper for multiple codim and multiple geometry types.
void addVertexData(const V &v, const std::string &name, int ncomps=1)
Add a grid function (represented by container) that lives on the vertices of the grid to the visualiz...
Definition: vtkwriter.hh:680
Take a vector and interpret it as cell data for the VTKWriter.
Definition: function.hh:88
CornerIterator cornerEnd() const
Definition: vtkwriter.hh:563
Iterator over the grids elements.
Definition: vtkwriter.hh:327
Iterate over the grid's vertices.
Definition: vtkwriter.hh:365
void endPoints()
finish section for the point coordinates
Definition: pvtuwriter.hh:169
virtual void unbind()
Unbind data set from current grid entity - mostly here for performance and symmetry reasons...
Definition: vtkwriter.hh:193
Dump a .vtu/.vtp files contents to a stream.
Definition: pvtuwriter.hh:60
FileType
which type of VTK file to write
Definition: common.hh:290
The dimension of the grid.
Definition: common/gridview.hh:130
virtual void writeGridPoints(VTK::VTUWriter &writer)
write the positions of vertices
Definition: vtkwriter.hh:1233
void addCellData(VTKFunction *p)
Add a grid function that lives on the cells of the grid to the visualization.
Definition: vtkwriter.hh:604
Writer for the ouput of grid functions in the vtk format.Writes arbitrary grid functions (living on c...
Definition: vtkwriter.hh:87
GridView gridView_
Definition: vtkwriter.hh:1306
std::string getParallelHeaderName(const std::string &name, const std::string &path, int commSize) const
return name of a parallel header file
Definition: vtkwriter.hh:799
int min(const DofVectorPointer< int > &dofVector)
Definition: dofvector.hh:346
Common stuff for the VTKWriter.
void clear()
clear list of registered functions
Definition: vtkwriter.hh:694
virtual void writeVertexData(VTK::VTUWriter &writer)
write vertex data
Definition: vtkwriter.hh:1219
Type erasure implementation for functions conforming to the dune-functions LocalFunction interface...
Definition: vtkwriter.hh:178
void beginMain(unsigned ghostLevel=0)
start the main PPolyData/PUnstructuredGrid section
Definition: pvtuwriter.hh:187
void endPointData()
finish PointData section
Definition: pvtuwriter.hh:127
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:749
virtual void unbind()
Unbind data set from current grid entity - mostly here for performance and symmetry reasons...
Definition: vtkwriter.hh:238
std::shared_ptr< FunctionWrapperBase > _f
Definition: vtkwriter.hh:314
std::list< VTKLocalFunction >::const_iterator FunctionIterator
Definition: vtkwriter.hh:319
Data array writers for the VTKWriter.
const FieldVector< DT, n > & position() const
position of vertex inside the entity
Definition: vtkwriter.hh:444
std::string getTypeString() const
Definition: vtkwriter.hh:1098
void addCellData(const VTKFunctionPtr &p)
Add a grid function that lives on the cells of the grid to the visualization.
Definition: vtkwriter.hh:588
GeometryType geometryType(const Dune::GeometryType &t)
mapping from GeometryType to VTKGeometryType
Definition: common.hh:195
void addCellData(F &&f, VTK::FieldInfo vtkFieldInfo)
Definition: vtkwriter.hh:594
Output non-conforming data.
Definition: common.hh:78
const VTK::FieldInfo & fieldInfo() const
Returns the VTK::FieldInfo for the data set.
Definition: vtkwriter.hh:291
Index subIndex(const typename GV::template Codim< 0 >::Entity &e, int i, unsigned int codim) const
Map subentity of codim 0 entity to array index.
Definition: mcmgmapper.hh:190
Traits::IndexSet IndexSet
type of the index set
Definition: common/gridview.hh:80
std::list< VTKLocalFunction > celldata
Definition: vtkwriter.hh:1302
virtual void write(T data)=0
write one data element
VTK::DataArrayWriter< float > Writer
Definition: vtkwriter.hh:153
Base class to write pvd-files which contains a list of all collected vtk-files.
Definition: vtksequencewriterbase.hh:32
void addPiece(const std::string &filename)
Add a serial piece to the output file.
Definition: pvtuwriter.hh:214
std::list< VTKLocalFunction > vertexdata
Definition: vtkwriter.hh:1303
void endPoints()
finish section for the point coordinates
Definition: vtuwriter.hh:247