6#ifndef DUNE_VTKWRITER_HH
7#define DUNE_VTKWRITER_HH
27#include <dune/geometry/referenceelements.hh>
29#include <dune/grid/common/gridenums.hh>
33#include <dune/grid/io/file/vtk/pvtuwriter.hh>
34#include <dune/grid/io/file/vtk/streams.hh>
35#include <dune/grid/io/file/vtk/vtuwriter.hh>
56 template<
class F,
class E,
class =
void >
61 template<
class F,
class E >
62 struct IsBindable< F, E,
std::
void_t< decltype( std::declval< F & >().bind( std::declval< const E & >() ) ),
63 decltype( std::declval< F & >().unbind() ) > >
68 template<
class F,
class =
void >
69 struct HasLocalFunction
74 struct HasLocalFunction< F,
std::
void_t< decltype( localFunction( std::declval< F& >() ) ) > >
81 template <
class Gr
idView>
82 class VTKSequenceWriterBase;
83 template <
class Gr
idView>
84 class VTKSequenceWriter;
94 template<
class Gr
idView >
118 typedef typename GridView::template
Codim< 0 >
119 ::template Partition< VTK_Partition >::Iterator
121 typedef typename GridView::template
Codim< n >
122 ::template Partition< VTK_Partition >::Iterator
125 typedef typename GridCellIterator::Reference EntityReference;
127 typedef typename GridView::template
Codim< 0 >
128 ::Entity::Geometry::LocalCoordinate Coordinate;
135 switch( VTK_Partition )
167 virtual void bind(
const Entity& e) = 0;
176 virtual void write(
const Coordinate& pos,
Writer& w, std::size_t count)
const = 0;
189 using Function =
typename std::decay<F>::type;
191 template<
typename F_>
193 : _f(std::forward<F_>(f))
196 virtual void bind(
const Entity& e)
206 virtual void write(
const Coordinate& pos,
Writer& w, std::size_t count)
const
216 void do_write(
Writer& w,
const R& r, std::size_t count, std::true_type)
const
218 for (std::size_t i = 0; i < count; ++i)
223 void do_write(Writer& w,
const R& r, std::size_t count, std::false_type)
const
237 using Function =
typename std::decay<F>::type;
239 template<
typename F_>
241 : _f(std::forward<F_>(f))
245 virtual void bind(
const Entity& e)
255 virtual void write(
const Coordinate& pos,
Writer& w, std::size_t count)
const
257 auto globalPos = element_->geometry().global(pos);
258 auto r = _f(globalPos);
260 for (std::size_t i = 0; i < count; ++i)
270 const Entity* element_;
282 virtual void bind(
const Entity& e)
292 virtual void write(
const Coordinate& pos,
Writer& w, std::size_t count)
const
294 for (std::size_t i = 0; i < count; ++i)
295 w.write(_f->evaluate(i,*_entity,pos));
300 std::shared_ptr< const VTKFunction > _f;
301 const Entity* _entity;
306 template<typename F, std::enable_if_t<Impl::IsBindable<F, Entity>::value,
int> = 0>
314 template<typename F, std::enable_if_t<not Impl::IsBindable<F, Entity>::value && Impl::HasLocalFunction<F>::value,
int> = 0>
317 typename
std::decay<decltype(localFunction(
std::forward<F>(f)))>::type
318 > >(localFunction(
std::forward<F>(f))))
324 template<typename F, std::enable_if_t<not Impl::IsBindable<F, Entity>::value && not Impl::HasLocalFunction<F>::value,
int> = 0>
334 vtkFunctionPtr->
name(),
335 (vtkFunctionPtr->ncomps() == 2 || vtkFunctionPtr->ncomps() == 3) ? VTK::FieldInfo::Type::vector : VTK::FieldInfo::Type::scalar,
336 vtkFunctionPtr->ncomps(),
337 vtkFunctionPtr->precision()
354 void bind(
const Entity& e)
const
371 std::shared_ptr<FunctionWrapperBase> _f;
376 typedef typename std::list<VTKLocalFunction>::const_iterator FunctionIterator;
397 CellIterator cellBegin()
const
399 return gridView_.template begin< 0, VTK_Partition >();
402 CellIterator cellEnd()
const
404 return gridView_.template end< 0, VTK_Partition >();
425 GridCellIterator git;
426 GridCellIterator gend;
432 std::vector<bool> visited;
440 void basicIncrement ()
445 const int numCorners = git->subEntities(n);
446 if( cornerIndexDune == numCorners )
448 offset += numCorners;
452 while( (git != gend) && skipEntity( git->partitionType() ) )
458 const GridCellIterator & end,
461 git(x), gend(end), datamode(dm), cornerIndexDune(0),
462 vertexmapper(vm), visited(vm.
size(),
false),
465 if (datamode == VTK::conforming && git != gend)
466 visited[vertexmapper.
subIndex(*git,cornerIndexDune,n)] =
true;
472 case VTK::conforming :
473 while(visited[vertexmapper.
subIndex(*git,cornerIndexDune,n)])
476 if (git == gend)
return;
478 visited[vertexmapper.
subIndex(*git,cornerIndexDune,n)] =
true;
480 case VTK::nonconforming :
487 return git == cit.git
488 && cornerIndexDune == cit.cornerIndexDune
489 && datamode == cit.datamode;
491 EntityReference dereference()
const
498 return cornerIndexDune;
503 return referenceElement<DT,n>(git->type())
504 .position(cornerIndexDune,n);
508 VertexIterator vertexBegin ()
const
510 return VertexIterator( gridView_.template begin< 0, VTK_Partition >(),
511 gridView_.template end< 0, VTK_Partition >(),
512 datamode, *vertexmapper );
515 VertexIterator vertexEnd ()
const
517 return VertexIterator( gridView_.template end< 0, VTK_Partition >(),
518 gridView_.template end< 0, VTK_Partition >(),
519 datamode, *vertexmapper );
540 GridCellIterator git;
541 GridCellIterator gend;
550 const std::vector<int> & number;
559 const GridCellIterator & end,
562 const std::vector<int> & num) :
563 git(x), gend(end), datamode(dm), cornerIndexVTK(0),
565 number(num), offset(0) {}
571 const int numCorners = git->subEntities(n);
572 if( cornerIndexVTK == numCorners )
574 offset += numCorners;
578 while( (git != gend) && skipEntity( git->partitionType() ) )
584 return git == cit.git
585 && cornerIndexVTK == cit.cornerIndexVTK
586 && datamode == cit.datamode;
588 EntityReference dereference()
const
601 case VTK::conforming :
603 number[vertexmapper.
subIndex(*git,VTK::renumber(*git,cornerIndexVTK),
605 case VTK::nonconforming :
606 return offset + VTK::renumber(*git,cornerIndexVTK);
613 CornerIterator cornerBegin ()
const
615 return CornerIterator( gridView_.template begin< 0, VTK_Partition >(),
616 gridView_.template end< 0, VTK_Partition >(),
617 datamode, *vertexmapper, number );
620 CornerIterator cornerEnd ()
const
622 return CornerIterator( gridView_.template end< 0, VTK_Partition >(),
623 gridView_.template end< 0, VTK_Partition >(),
624 datamode, *vertexmapper, number );
639 : gridView_( gridView ),
642 polyhedralCellsPresent_( checkForPolyhedralCells() )
694 template<
class Container>
695 void addCellData (
const Container& v,
const std::string &name,
int ncomps = 1,
699 for (
int c=0; c<ncomps; ++c) {
700 std::stringstream compName;
703 compName <<
"[" << c <<
"]";
705 addCellData(std::shared_ptr< const VTKFunction >(p));
759 template<
class Container>
760 void addVertexData (
const Container& v,
const std::string &name,
int ncomps=1,
764 for (
int c=0; c<ncomps; ++c) {
765 std::stringstream compName;
768 compName <<
"[" << c <<
"]";
783 {
return coordPrec; }
803 std::string
write (
const std::string &name,
806 return write( name, type, gridView_.
comm().rank(), gridView_.
comm().size() );
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() );
856 const std::string& path,
857 int commRank,
int commSize)
const
859 std::ostringstream s;
864 if(path[path.size()-1] !=
'/')
868 std::string fileprefix;
871 auto pos = name.rfind(
'/');
872 if( pos != std::string::npos )
875 fileprefix = name.substr( pos+1 );
878 std::string newpath = name.substr(0, pos);
880 if(newpath[name.size()-1] !=
'/')
889 s <<
's' << std::setw(4) << std::setfill(
'0') << commSize <<
'-';
890 const bool writeHeader = commRank < 0;
893 s <<
'p' << std::setw(4) << std::setfill(
'0') << commRank <<
'-';
896 s << fileprefix <<
".";
921 const std::string& path,
941 const std::string& path)
const
943 static const std::string extension =
965 std::string
write (
const std::string &name,
974 std::string filename = name;
975 std::string path = std::string(
"");
979 auto pos = name.rfind(
'/');
980 if( pos != std::string::npos )
983 filename = name.substr( pos+1 );
987 path = name.substr(0, pos);
990 return pwrite(filename, path,
"", type, commRank, commSize);
1001 file.exceptions(std::ios_base::badbit | std::ios_base::failbit |
1002 std::ios_base::eofbit);
1005 file.open( pieceName.c_str(), std::ios::binary );
1008 std::cerr <<
"Filename: " << pieceName <<
" could not be opened" << std::endl;
1011 if (! file.is_open())
1013 writeDataFile( file );
1043 std::string
pwrite(
const std::string& name,
const std::string& path,
1044 const std::string& extendpath,
1046 const int commSize )
1053 file.exceptions(std::ios_base::badbit | std::ios_base::failbit |
1054 std::ios_base::eofbit);
1055 std::string piecepath =
concatPaths(path, extendpath);
1056 std::string relpiecepath =
relativePath(path, piecepath);
1063 file.open(fullname.c_str(),std::ios::binary);
1066 std::cerr <<
"Filename: " << fullname <<
" could not be opened" << std::endl;
1069 if (! file.is_open())
1071 writeDataFile(file);
1073 gridView_.
comm().barrier();
1079 file.open(fullname.c_str());
1080 if (! file.is_open())
1082 writeParallelHeader(file,name,relpiecepath, commSize );
1085 gridView_.
comm().barrier();
1108 void writeParallelHeader(std::ostream& s,
const std::string& piecename,
1109 const std::string& piecepath,
const int commSize)
1112 (n == 1) ? VTK::polyData : VTK::unstructuredGrid;
1120 std::string scalars, vectors;
1121 std::tie(scalars,vectors) = getDataNames(vertexdata);
1122 writer.beginPointData(scalars, vectors);
1124 for (
auto it = vertexdata.begin(),
1125 end = vertexdata.end();
1129 unsigned writecomps = it->fieldInfo().size();
1133 writer.addArray(it->name(), writecomps, it->fieldInfo().precision());
1135 writer.endPointData();
1139 std::string scalars, vectors;
1140 std::tie(scalars,vectors) = getDataNames(celldata);
1141 writer.beginCellData(scalars, vectors);
1143 for (
auto it = celldata.begin(),
1144 end = celldata.end();
1148 unsigned writecomps = it->fieldInfo().size();
1152 writer.addArray(it->name(), writecomps, it->fieldInfo().precision());
1154 writer.endCellData();
1157 writer.beginPoints();
1158 writer.addArray(
"Coordinates", 3, coordPrec);
1162 for(
int i = 0; i < commSize; ++i )
1167 writer.addPiece(fullname);
1174 void writeDataFile (std::ostream& s)
1176 VTK::FileType fileType =
1177 (n == 1) ? VTK::polyData : VTK::unstructuredGrid;
1179 VTK::VTUWriter writer(s, outputtype, fileType);
1183 if (datamode == VTK::conforming)
1185 number.resize(vertexmapper->
size());
1186 for (std::vector<int>::size_type i=0; i<number.size(); i++) number[i] = -1;
1190 writer.beginMain(ncells, nvertices);
1191 writeAllData(writer);
1195 if(writer.beginAppended())
1196 writeAllData(writer);
1197 writer.endAppended();
1199 delete vertexmapper; number.clear();
1202 void writeAllData(VTK::VTUWriter& writer) {
1217 std::string getFormatString()
const
1219 if (outputtype==VTK::ascii)
1221 if (outputtype==VTK::base64)
1223 if (outputtype==VTK::appendedraw)
1225 if (outputtype==VTK::appendedbase64)
1227 DUNE_THROW(IOError,
"VTKWriter: unsupported OutputType" << outputtype);
1230 std::string getTypeString()
const
1235 return "UnstructuredGrid";
1249 const int subEntities = it->subEntities(n);
1250 for (
int i=0; i<subEntities; ++i)
1253 if (datamode == VTK::conforming)
1255 int alpha = vertexmapper->
subIndex(*it,i,n);
1256 if (number[alpha]<0)
1257 number[alpha] = nvertices_++;
1267 template<
typename T>
1268 std::tuple<std::string,std::string> getDataNames(
const T& data)
const
1270 std::string scalars =
"";
1271 for (
auto it = data.begin(),
1277 scalars = it->name();
1281 std::string vectors =
"";
1282 for (
auto it = data.begin(),
1288 vectors = it->name();
1291 return std::make_tuple(scalars,vectors);
1294 template<
typename Data,
typename Iterator>
1295 void writeData(VTK::VTUWriter& writer,
const Data& data,
const Iterator begin,
const Iterator end,
int nentries)
1297 for (
auto it = data.begin(),
1302 const auto& f = *it;
1303 VTK::FieldInfo fieldInfo = f.fieldInfo();
1304 std::size_t writecomps = fieldInfo.size();
1305 switch (fieldInfo.type())
1313 DUNE_THROW(IOError,
"Cannot write VTK vectors with more than 3 components (components was " << writecomps <<
")");
1317 DUNE_THROW(NotImplemented,
"VTK output for tensors not implemented yet");
1319 std::shared_ptr<VTK::DataArrayWriter> p
1320 (writer.makeArrayWriter(f.name(), writecomps, nentries, fieldInfo.precision()));
1321 if(!p->writeIsNoop())
1322 for (Iterator eit = begin; eit!=end; ++eit)
1324 const Entity & e = *eit;
1326 f.write(eit.position(),*p);
1330 for (std::size_t j=fieldInfo.size(); j < writecomps; ++j)
1339 if(celldata.size() == 0)
1342 std::string scalars, vectors;
1343 std::tie(scalars,vectors) = getDataNames(celldata);
1346 writeData(writer,celldata,cellBegin(),cellEnd(),ncells);
1353 if(vertexdata.size() == 0)
1356 std::string scalars, vectors;
1357 std::tie(scalars,vectors) = getDataNames(vertexdata);
1360 writeData(writer,vertexdata,vertexBegin(),vertexEnd(),nvertices);
1369 std::shared_ptr<VTK::DataArrayWriter> p
1371 if(!p->writeIsNoop()) {
1376 for (
int j=0; j<
std::min(dimw,3); j++)
1377 p->write((*vit).geometry().corner(vit.localindex())[j]);
1378 for (
int j=
std::min(dimw,3); j<3; j++)
1395 std::shared_ptr<VTK::DataArrayWriter> p1
1396 (writer.
makeArrayWriter(
"connectivity", 1, ncorners, VTK::Precision::int32));
1397 if(!p1->writeIsNoop())
1404 std::shared_ptr<VTK::DataArrayWriter> p2
1405 (writer.
makeArrayWriter(
"offsets", 1, ncells, VTK::Precision::int32));
1406 if(!p2->writeIsNoop()) {
1410 offset += it->subEntities(n);
1420 std::shared_ptr<VTK::DataArrayWriter> p3
1423 if(!p3->writeIsNoop())
1427 int vtktype = VTK::geometryType(it->type());
1435 if( polyhedralCellsPresent_ )
1445 bool checkForPolyhedralCells()
const
1448 for(
const auto& geomType : gridView_.
indexSet().types( 0 ) )
1450 if( VTK::geometryType( geomType ) == VTK::polyhedron )
1461 if( ! faceVertices_ )
1463 faceVertices_.reset(
new std::pair< std::vector<int>, std::vector<int> > () );
1465 fillFaceVertices( cornerBegin(), cornerEnd(), gridView_.
indexSet(),
1466 faceVertices_->first, faceVertices_->second );
1469 std::vector< int >& faces = faceVertices_->first;
1470 std::vector< int >& faceOffsets = faceVertices_->second;
1471 assert(
int(faceOffsets.size()) == ncells );
1474 std::shared_ptr<VTK::DataArrayWriter> p4
1475 (writer.
makeArrayWriter(
"faces", 1, faces.size(), VTK::Precision::int32));
1476 if(!p4->writeIsNoop())
1478 for(
const auto& face : faces )
1484 std::shared_ptr<VTK::DataArrayWriter> p5
1485 (writer.
makeArrayWriter(
"faceoffsets", 1, ncells, VTK::Precision::int32));
1486 if(!p5->writeIsNoop())
1488 for(
const auto& offset : faceOffsets )
1489 p5->write( offset );
1492 faceVertices_.reset();
1497 template <
class CornerIterator,
class IndexSet,
class T>
1498 inline void fillFaceVertices( CornerIterator it,
1499 const CornerIterator end,
1501 std::vector<T>& faces,
1502 std::vector<T>& faceOffsets )
1504 if( n == 3 && it != end )
1508 faces.reserve( 15 * ncells );
1509 faceOffsets.clear();
1510 faceOffsets.reserve( ncells );
1515 int elIndex = indexSet.
index( element );
1516 std::vector< T > vertices;
1517 vertices.reserve( 30 );
1518 for( ; it != end; ++it )
1520 const Cell& cell = *it ;
1521 const int cellIndex = indexSet.
index( cell ) ;
1522 if( elIndex != cellIndex )
1524 fillFacesForElement( element, indexSet, vertices, offset, faces, faceOffsets );
1528 elIndex = cellIndex ;
1530 vertices.push_back( it.id() );
1534 fillFacesForElement( element, indexSet, vertices, offset, faces, faceOffsets );
1538 template <
class Entity,
class IndexSet,
class T>
1539 static void fillFacesForElement(
const Entity& element,
1540 const IndexSet& indexSet,
1541 const std::vector<T>& vertices,
1543 std::vector<T>& faces,
1544 std::vector<T>& faceOffsets )
1548 std::map< T, T > vxMap;
1551 const int nVertices = element.subEntities( dim );
1552 for(
int vx = 0; vx < nVertices; ++ vx )
1554 const int vxIdx = indexSet.subIndex( element, vx, dim );
1555 vxMap[ vxIdx ] = vertices[ vx ];
1559 const int nFaces = element.subEntities( 1 );
1561 faces.push_back( nFaces );
1564 for(
int fce = 0; fce < nFaces; ++ fce )
1567 const auto face = element.template subEntity< 1 > ( fce );
1570 const int nVxFace = face.subEntities( dim );
1571 faces.push_back( nVxFace );
1573 for(
int i=0; i<nVxFace; ++i )
1575 const T vxIndex = indexSet.subIndex( face, i, dim );
1576 assert( vxMap.find( vxIndex ) != vxMap.end() );
1577 faces.push_back( vxMap[ vxIndex ] );
1583 faceOffsets.push_back( offset );
1588 std::list<VTKLocalFunction> celldata;
1589 std::list<VTKLocalFunction> vertexdata;
1599 VertexMapper* vertexmapper;
1602 std::vector<int> number;
1603 VTK::DataMode datamode;
1604 VTK::Precision coordPrec;
1607 const bool polyhedralCellsPresent_;
1610 std::shared_ptr< std::pair< std::vector<int>, std::vector<int> > > faceVertices_;
1613 VTK::OutputType outputtype;
vector space out of a tensor product of fields.
Definition: fvector.hh:95
Base class for stl conformant forward iterators.
Definition: iteratorfacades.hh:141
Base class template for function classes.
Definition: function.hh:41
Grid view abstract base class.
Definition: gridview.hh:66
Default exception class for I/O errors.
Definition: exceptions.hh:231
Index Set Interface base class.
Definition: indexidset.hh:78
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
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
Default exception for dummy implementations.
Definition: exceptions.hh:263
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 ouput of grid functions in the vtk format.
Definition: vtksequencewriter.hh:29
Iterator over the grids elements.
Definition: vtkwriter.hh:385
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
Iterate over the elements' corners.
Definition: vtkwriter.hh:539
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
Iterate over the grid's vertices.
Definition: vtkwriter.hh:424
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 ouput 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
void clear()
clear list of registered functions
Definition: vtkwriter.hh:775
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
virtual ~VTKWriter()
destructor
Definition: vtkwriter.hh:786
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
base class for data array writers
Definition: dataarraywriter.hh:56
Descriptor struct for VTK fields.
Definition: common.hh:328
@ tensor
tensor field (always 3x3)
@ 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
Dump a .vtu/.vtp files contents to a stream.
Definition: pvtuwriter.hh:62
Dump a .vtu/.vtp files contents to a stream.
Definition: vtuwriter.hh:98
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.
Traits for type conversions and type information.
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
FileType
which type of VTK file to write
Definition: common.hh:252
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
#define DUNE_THROW(E, m)
Definition: exceptions.hh:218
const IndexSet & indexSet() const
obtain the index set
Definition: gridview.hh:191
Traits::Grid Grid
type of the grid
Definition: gridview.hh:83
Traits::IndexSet IndexSet
type of the index set
Definition: gridview.hh:86
const Communication & comm() const
obtain communication object
Definition: gridview.hh:280
static constexpr int dimension
The dimension of the grid.
Definition: gridview.hh:148
Grid::ctype ctype
type used for coordinates in grid
Definition: gridview.hh:145
static constexpr int dimensionworld
The dimension of the world the grid lives in.
Definition: gridview.hh:151
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:153
std::string concatPaths(const std::string &base, const std::string &p)
concatenate two paths
Definition: path.cc:32
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.
This file implements iterator facade classes for writing stl conformant iterators.
Mapper for multiple codim and multiple geometry types.
Dune namespace.
Definition: alignedallocator.hh:13
Utilities for handling filesystem paths.
Static tag representing a codimension.
Definition: dimension.hh:24
static const ReferenceElement & general(const GeometryType &type)
get general reference elements
Definition: referenceelements.hh:198
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
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