4#ifndef DUNE_VTKWRITER_HH
5#define DUNE_VTKWRITER_HH
25#include <dune/geometry/referenceelements.hh>
27#include <dune/grid/common/gridenums.hh>
31#include <dune/grid/io/file/vtk/pvtuwriter.hh>
32#include <dune/grid/io/file/vtk/streams.hh>
33#include <dune/grid/io/file/vtk/vtuwriter.hh>
54 template<
class F,
class E,
class =
void >
59 template<
class F,
class E >
60 struct IsBindable< F, E,
std::
void_t< decltype( std::declval< F & >().bind( std::declval< const E & >() ) ),
61 decltype( std::declval< F & >().unbind() ) > >
66 template<
class F,
class =
void >
67 struct HasLocalFunction
72 struct HasLocalFunction< F,
std::
void_t< decltype( localFunction( std::declval< F& >() ) ) > >
79 template <
class Gr
idView>
80 class VTKSequenceWriterBase;
81 template <
class Gr
idView>
82 class VTKSequenceWriter;
92 template<
class Gr
idView >
116 typedef typename GridView::template
Codim< 0 >
117 ::template Partition< VTK_Partition >::Iterator
119 typedef typename GridView::template
Codim< n >
120 ::template Partition< VTK_Partition >::Iterator
123 typedef typename GridCellIterator::Reference EntityReference;
125 typedef typename GridView::template
Codim< 0 >
126 ::Entity::Geometry::LocalCoordinate Coordinate;
133 switch( VTK_Partition )
165 virtual void bind(
const Entity& e) = 0;
174 virtual void write(
const Coordinate& pos,
Writer& w, std::size_t count)
const = 0;
187 using Function =
typename std::decay<F>::type;
189 template<
typename F_>
191 : _f(std::forward<F_>(f))
194 virtual void bind(
const Entity& e)
204 virtual void write(
const Coordinate& pos,
Writer& w, std::size_t count)
const
214 void do_write(
Writer& w,
const R& r, std::size_t count, std::true_type)
const
216 for (std::size_t i = 0; i < count; ++i)
221 void do_write(Writer& w,
const R& r, std::size_t count, std::false_type)
const
235 using Function =
typename std::decay<F>::type;
237 template<
typename F_>
239 : _f(std::forward<F_>(f))
243 virtual void bind(
const Entity& e)
253 virtual void write(
const Coordinate& pos,
Writer& w, std::size_t count)
const
255 auto globalPos = element_->geometry().global(pos);
256 auto r = _f(globalPos);
258 for (std::size_t i = 0; i < count; ++i)
268 const Entity* element_;
280 virtual void bind(
const Entity& e)
290 virtual void write(
const Coordinate& pos,
Writer& w, std::size_t count)
const
292 for (std::size_t i = 0; i < count; ++i)
293 w.write(_f->evaluate(i,*_entity,pos));
298 std::shared_ptr< const VTKFunction > _f;
299 const Entity* _entity;
304 template<typename F, std::enable_if_t<Impl::IsBindable<F, Entity>::value,
int> = 0>
312 template<typename F, std::enable_if_t<not Impl::IsBindable<F, Entity>::value && Impl::HasLocalFunction<F>::value,
int> = 0>
322 template<typename F, std::enable_if_t<not Impl::IsBindable<F, Entity>::value && not Impl::HasLocalFunction<F>::value,
int> = 0>
332 vtkFunctionPtr->
name(),
333 (vtkFunctionPtr->ncomps() == 2 || vtkFunctionPtr->ncomps() == 3) ? VTK::FieldInfo::Type::vector : VTK::FieldInfo::Type::scalar,
334 vtkFunctionPtr->ncomps(),
335 vtkFunctionPtr->precision()
352 void bind(
const Entity& e)
const
369 std::shared_ptr<FunctionWrapperBase> _f;
374 typedef typename std::list<VTKLocalFunction>::const_iterator FunctionIterator;
395 CellIterator cellBegin()
const
397 return gridView_.template begin< 0, VTK_Partition >();
400 CellIterator cellEnd()
const
402 return gridView_.template end< 0, VTK_Partition >();
423 GridCellIterator git;
424 GridCellIterator gend;
430 std::vector<bool> visited;
438 void basicIncrement ()
443 const int numCorners = git->subEntities(n);
444 if( cornerIndexDune == numCorners )
446 offset += numCorners;
450 while( (git != gend) && skipEntity( git->partitionType() ) )
456 const GridCellIterator & end,
459 git(x), gend(end), datamode(dm), cornerIndexDune(0),
460 vertexmapper(vm), visited(vm.
size(),
false),
463 if (datamode == VTK::conforming && git != gend)
464 visited[vertexmapper.
subIndex(*git,cornerIndexDune,n)] =
true;
470 case VTK::conforming :
471 while(visited[vertexmapper.
subIndex(*git,cornerIndexDune,n)])
474 if (git == gend)
return;
476 visited[vertexmapper.
subIndex(*git,cornerIndexDune,n)] =
true;
478 case VTK::nonconforming :
485 return git == cit.git
486 && cornerIndexDune == cit.cornerIndexDune
487 && datamode == cit.datamode;
489 EntityReference dereference()
const
496 return cornerIndexDune;
501 return referenceElement<DT,n>(git->type())
502 .position(cornerIndexDune,n);
506 VertexIterator vertexBegin ()
const
508 return VertexIterator( gridView_.template begin< 0, VTK_Partition >(),
509 gridView_.template end< 0, VTK_Partition >(),
510 datamode, *vertexmapper );
513 VertexIterator vertexEnd ()
const
515 return VertexIterator( gridView_.template end< 0, VTK_Partition >(),
516 gridView_.template end< 0, VTK_Partition >(),
517 datamode, *vertexmapper );
538 GridCellIterator git;
539 GridCellIterator gend;
548 const std::vector<int> & number;
557 const GridCellIterator & end,
560 const std::vector<int> & num) :
561 git(x), gend(end), datamode(dm), cornerIndexVTK(0),
563 number(num), offset(0) {}
569 const int numCorners = git->subEntities(n);
570 if( cornerIndexVTK == numCorners )
572 offset += numCorners;
576 while( (git != gend) && skipEntity( git->partitionType() ) )
582 return git == cit.git
583 && cornerIndexVTK == cit.cornerIndexVTK
584 && datamode == cit.datamode;
586 EntityReference dereference()
const
599 case VTK::conforming :
601 number[vertexmapper.
subIndex(*git,VTK::renumber(*git,cornerIndexVTK),
603 case VTK::nonconforming :
604 return offset + VTK::renumber(*git,cornerIndexVTK);
611 CornerIterator cornerBegin ()
const
613 return CornerIterator( gridView_.template begin< 0, VTK_Partition >(),
614 gridView_.template end< 0, VTK_Partition >(),
615 datamode, *vertexmapper, number );
618 CornerIterator cornerEnd ()
const
620 return CornerIterator( gridView_.template end< 0, VTK_Partition >(),
621 gridView_.template end< 0, VTK_Partition >(),
622 datamode, *vertexmapper, number );
637 : gridView_( gridView ),
640 polyhedralCellsPresent_( checkForPolyhedralCells() )
692 template<
class Container>
693 void addCellData (
const Container& v,
const std::string &name,
int ncomps = 1,
697 for (
int c=0; c<ncomps; ++c) {
698 std::stringstream compName;
701 compName <<
"[" << c <<
"]";
703 addCellData(std::shared_ptr< const VTKFunction >(p));
757 template<
class Container>
758 void addVertexData (
const Container& v,
const std::string &name,
int ncomps=1,
762 for (
int c=0; c<ncomps; ++c) {
763 std::stringstream compName;
766 compName <<
"[" << c <<
"]";
781 {
return coordPrec; }
801 std::string
write (
const std::string &name,
804 return write( name, type, gridView_.
comm().rank(), gridView_.
comm().size() );
833 std::string
pwrite (
const std::string & name,
const std::string & path,
const std::string & extendpath,
836 return pwrite( name, path, extendpath, type, gridView_.
comm().rank(), gridView_.
comm().size() );
854 const std::string& path,
855 int commRank,
int commSize)
const
857 std::ostringstream s;
862 if(path[path.size()-1] !=
'/')
866 std::string fileprefix;
869 auto pos = name.rfind(
'/');
870 if( pos != std::string::npos )
873 fileprefix = name.substr( pos+1 );
876 std::string newpath = name.substr(0, pos);
878 if(newpath[name.size()-1] !=
'/')
887 s <<
's' << std::setw(4) << std::setfill(
'0') << commSize <<
'-';
888 const bool writeHeader = commRank < 0;
891 s <<
'p' << std::setw(4) << std::setfill(
'0') << commRank <<
'-';
894 s << fileprefix <<
".";
919 const std::string& path,
939 const std::string& path)
const
941 static const std::string extension =
963 std::string
write (
const std::string &name,
972 std::string filename = name;
973 std::string path = std::string(
"");
977 auto pos = name.rfind(
'/');
978 if( pos != std::string::npos )
981 filename = name.substr( pos+1 );
985 path = name.substr(0, pos);
988 return pwrite(filename, path,
"", type, commRank, commSize);
999 file.exceptions(std::ios_base::badbit | std::ios_base::failbit |
1000 std::ios_base::eofbit);
1003 file.open( pieceName.c_str(), std::ios::binary );
1006 std::cerr <<
"Filename: " << pieceName <<
" could not be opened" << std::endl;
1009 if (! file.is_open())
1011 writeDataFile( file );
1041 std::string
pwrite(
const std::string& name,
const std::string& path,
1042 const std::string& extendpath,
1044 const int commSize )
1051 file.exceptions(std::ios_base::badbit | std::ios_base::failbit |
1052 std::ios_base::eofbit);
1053 std::string piecepath =
concatPaths(path, extendpath);
1054 std::string relpiecepath =
relativePath(path, piecepath);
1061 file.open(fullname.c_str(),std::ios::binary);
1064 std::cerr <<
"Filename: " << fullname <<
" could not be opened" << std::endl;
1067 if (! file.is_open())
1069 writeDataFile(file);
1071 gridView_.
comm().barrier();
1077 file.open(fullname.c_str());
1078 if (! file.is_open())
1080 writeParallelHeader(file,name,relpiecepath, commSize );
1083 gridView_.
comm().barrier();
1106 void writeParallelHeader(std::ostream& s,
const std::string& piecename,
1107 const std::string& piecepath,
const int commSize)
1110 (n == 1) ? VTK::polyData : VTK::unstructuredGrid;
1118 std::string scalars, vectors;
1119 std::tie(scalars,vectors) = getDataNames(vertexdata);
1120 writer.beginPointData(scalars, vectors);
1122 for (
auto it = vertexdata.begin(),
1123 end = vertexdata.end();
1127 unsigned writecomps = it->fieldInfo().size();
1128 if(writecomps == 2) writecomps = 3;
1129 writer.addArray(it->name(), writecomps, it->fieldInfo().precision());
1131 writer.endPointData();
1135 std::string scalars, vectors;
1136 std::tie(scalars,vectors) = getDataNames(celldata);
1137 writer.beginCellData(scalars, vectors);
1139 for (
auto it = celldata.begin(),
1140 end = celldata.end();
1144 unsigned writecomps = it->fieldInfo().size();
1145 if(writecomps == 2) writecomps = 3;
1146 writer.addArray(it->name(), writecomps, it->fieldInfo().precision());
1148 writer.endCellData();
1151 writer.beginPoints();
1152 writer.addArray(
"Coordinates", 3, coordPrec);
1156 for(
int i = 0; i < commSize; ++i )
1161 writer.addPiece(fullname);
1168 void writeDataFile (std::ostream& s)
1170 VTK::FileType fileType =
1171 (n == 1) ? VTK::polyData : VTK::unstructuredGrid;
1173 VTK::VTUWriter writer(s, outputtype, fileType);
1177 if (datamode == VTK::conforming)
1179 number.resize(vertexmapper->
size());
1180 for (std::vector<int>::size_type i=0; i<number.size(); i++) number[i] = -1;
1184 writer.beginMain(ncells, nvertices);
1185 writeAllData(writer);
1189 if(writer.beginAppended())
1190 writeAllData(writer);
1191 writer.endAppended();
1193 delete vertexmapper; number.clear();
1196 void writeAllData(VTK::VTUWriter& writer) {
1211 std::string getFormatString()
const
1213 if (outputtype==VTK::ascii)
1215 if (outputtype==VTK::base64)
1217 if (outputtype==VTK::appendedraw)
1219 if (outputtype==VTK::appendedbase64)
1221 DUNE_THROW(IOError,
"VTKWriter: unsupported OutputType" << outputtype);
1224 std::string getTypeString()
const
1229 return "UnstructuredGrid";
1243 const int subEntities = it->subEntities(n);
1244 for (
int i=0; i<subEntities; ++i)
1247 if (datamode == VTK::conforming)
1249 int alpha = vertexmapper->
subIndex(*it,i,n);
1250 if (number[alpha]<0)
1251 number[alpha] = nvertices_++;
1261 template<
typename T>
1262 std::tuple<std::string,std::string> getDataNames(
const T& data)
const
1264 std::string scalars =
"";
1265 for (
auto it = data.begin(),
1271 scalars = it->name();
1275 std::string vectors =
"";
1276 for (
auto it = data.begin(),
1282 vectors = it->name();
1285 return std::make_tuple(scalars,vectors);
1288 template<
typename Data,
typename Iterator>
1289 void writeData(VTK::VTUWriter& writer,
const Data& data,
const Iterator begin,
const Iterator end,
int nentries)
1291 for (
auto it = data.begin(),
1296 const auto& f = *it;
1297 VTK::FieldInfo fieldInfo = f.fieldInfo();
1298 std::size_t writecomps = fieldInfo.size();
1299 switch (fieldInfo.type())
1307 DUNE_THROW(IOError,
"Cannot write VTK vectors with more than 3 components (components was " << writecomps <<
")");
1311 DUNE_THROW(NotImplemented,
"VTK output for tensors not implemented yet");
1313 std::shared_ptr<VTK::DataArrayWriter> p
1314 (writer.makeArrayWriter(f.name(), writecomps, nentries, fieldInfo.precision()));
1315 if(!p->writeIsNoop())
1316 for (Iterator eit = begin; eit!=end; ++eit)
1318 const Entity & e = *eit;
1320 f.write(eit.position(),*p);
1324 for (std::size_t j=fieldInfo.size(); j < writecomps; ++j)
1333 if(celldata.size() == 0)
1336 std::string scalars, vectors;
1337 std::tie(scalars,vectors) = getDataNames(celldata);
1340 writeData(writer,celldata,cellBegin(),cellEnd(),ncells);
1347 if(vertexdata.size() == 0)
1350 std::string scalars, vectors;
1351 std::tie(scalars,vectors) = getDataNames(vertexdata);
1354 writeData(writer,vertexdata,vertexBegin(),vertexEnd(),nvertices);
1363 std::shared_ptr<VTK::DataArrayWriter> p
1365 if(!p->writeIsNoop()) {
1370 for (
int j=0; j<
std::min(dimw,3); j++)
1371 p->write((*vit).geometry().corner(vit.localindex())[j]);
1372 for (
int j=
std::min(dimw,3); j<3; j++)
1389 std::shared_ptr<VTK::DataArrayWriter> p1
1390 (writer.
makeArrayWriter(
"connectivity", 1, ncorners, VTK::Precision::int32));
1391 if(!p1->writeIsNoop())
1398 std::shared_ptr<VTK::DataArrayWriter> p2
1399 (writer.
makeArrayWriter(
"offsets", 1, ncells, VTK::Precision::int32));
1400 if(!p2->writeIsNoop()) {
1404 offset += it->subEntities(n);
1414 std::shared_ptr<VTK::DataArrayWriter> p3
1417 if(!p3->writeIsNoop())
1421 int vtktype = VTK::geometryType(it->type());
1429 if( polyhedralCellsPresent_ )
1439 bool checkForPolyhedralCells()
const
1442 for(
const auto& geomType : gridView_.
indexSet().types( 0 ) )
1444 if( VTK::geometryType( geomType ) == VTK::polyhedron )
1455 if( ! faceVertices_ )
1457 faceVertices_.reset(
new std::pair< std::vector<int>, std::vector<int> > () );
1459 fillFaceVertices( cornerBegin(), cornerEnd(), gridView_.
indexSet(),
1460 faceVertices_->first, faceVertices_->second );
1463 std::vector< int >& faces = faceVertices_->first;
1464 std::vector< int >& faceOffsets = faceVertices_->second;
1465 assert(
int(faceOffsets.size()) == ncells );
1468 std::shared_ptr<VTK::DataArrayWriter> p4
1469 (writer.
makeArrayWriter(
"faces", 1, faces.size(), VTK::Precision::int32));
1470 if(!p4->writeIsNoop())
1472 for(
const auto& face : faces )
1478 std::shared_ptr<VTK::DataArrayWriter> p5
1479 (writer.
makeArrayWriter(
"faceoffsets", 1, ncells, VTK::Precision::int32));
1480 if(!p5->writeIsNoop())
1482 for(
const auto& offset : faceOffsets )
1483 p5->write( offset );
1486 faceVertices_.reset();
1491 template <
class CornerIterator,
class IndexSet,
class T>
1492 inline void fillFaceVertices( CornerIterator it,
1493 const CornerIterator end,
1495 std::vector<T>& faces,
1496 std::vector<T>& faceOffsets )
1498 if( n == 3 && it != end )
1502 faces.reserve( 15 * ncells );
1503 faceOffsets.clear();
1504 faceOffsets.reserve( ncells );
1509 int elIndex = indexSet.
index( element );
1510 std::vector< T > vertices;
1511 vertices.reserve( 30 );
1512 for( ; it != end; ++it )
1514 const Cell& cell = *it ;
1515 const int cellIndex = indexSet.
index( cell ) ;
1516 if( elIndex != cellIndex )
1518 fillFacesForElement( element, indexSet, vertices, offset, faces, faceOffsets );
1522 elIndex = cellIndex ;
1524 vertices.push_back( it.id() );
1528 fillFacesForElement( element, indexSet, vertices, offset, faces, faceOffsets );
1532 template <
class Entity,
class IndexSet,
class T>
1533 static void fillFacesForElement(
const Entity& element,
1534 const IndexSet& indexSet,
1535 const std::vector<T>& vertices,
1537 std::vector<T>& faces,
1538 std::vector<T>& faceOffsets )
1542 std::map< T, T > vxMap;
1545 const int nVertices = element.subEntities( dim );
1546 for(
int vx = 0; vx < nVertices; ++ vx )
1548 const int vxIdx = indexSet.subIndex( element, vx, dim );
1549 vxMap[ vxIdx ] = vertices[ vx ];
1553 const int nFaces = element.subEntities( 1 );
1555 faces.push_back( nFaces );
1558 for(
int fce = 0; fce < nFaces; ++ fce )
1561 const auto face = element.template subEntity< 1 > ( fce );
1564 const int nVxFace = face.subEntities( dim );
1565 faces.push_back( nVxFace );
1567 for(
int i=0; i<nVxFace; ++i )
1569 const T vxIndex = indexSet.subIndex( face, i, dim );
1570 assert( vxMap.find( vxIndex ) != vxMap.end() );
1571 faces.push_back( vxMap[ vxIndex ] );
1577 faceOffsets.push_back( offset );
1582 std::list<VTKLocalFunction> celldata;
1583 std::list<VTKLocalFunction> vertexdata;
1593 VertexMapper* vertexmapper;
1596 std::vector<int> number;
1597 VTK::DataMode datamode;
1598 VTK::Precision coordPrec;
1601 const bool polyhedralCellsPresent_;
1604 std::shared_ptr< std::pair< std::vector<int>, std::vector<int> > > faceVertices_;
1607 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:139
Base class template for function classes.
Definition: function.hh:39
Grid view abstract base class.
Definition: gridview.hh:63
Default exception class for I/O errors.
Definition: exceptions.hh:229
Index Set Interface base class.
Definition: indexidset.hh:76
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:111
size_type size() const
Return total number of entities in the entity set managed by the mapper.
Definition: mcmgmapper.hh:202
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:183
Default exception for dummy implementations.
Definition: exceptions.hh:261
Take a vector and interpret it as cell data for the VTKWriter.
Definition: function.hh:95
Take a vector and interpret it as point data for the VTKWriter.
Definition: function.hh:203
A base class for grid functions with any return type and dimension.
Definition: function.hh:40
Base class to write pvd-files which contains a list of all collected vtk-files.
Definition: vtksequencewriterbase.hh:32
Writer for the ouput of grid functions in the vtk format.
Definition: vtksequencewriter.hh:27
Iterator over the grids elements.
Definition: vtkwriter.hh:383
CellIterator(const GridCellIterator &x)
construct a CellIterator from the gridview's Iterator.
Definition: vtkwriter.hh:386
const FieldVector< DT, n > position() const
Definition: vtkwriter.hh:389
Iterate over the elements' corners.
Definition: vtkwriter.hh:537
int id() const
Process-local consecutive zero-starting vertex id.
Definition: vtkwriter.hh:595
Type erasure wrapper for VTK data sets.
Definition: vtkwriter.hh:154
void unbind() const
Unbind the data set from the currently bound entity.
Definition: vtkwriter.hh:358
VTKLocalFunction(F &&f, VTK::FieldInfo fieldInfo)
Construct a VTKLocalFunction for a dune-functions style LocalFunction.
Definition: vtkwriter.hh:305
std::string name() const
Returns the name of the data set.
Definition: vtkwriter.hh:340
const VTK::FieldInfo & fieldInfo() const
Returns the VTK::FieldInfo for the data set.
Definition: vtkwriter.hh:346
void bind(const Entity &e) const
Bind the data set to grid entity e.
Definition: vtkwriter.hh:352
VTKLocalFunction(const std::shared_ptr< const VTKFunction > &vtkFunctionPtr)
Construct a VTKLocalFunction for a legacy VTKFunction.
Definition: vtkwriter.hh:329
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:364
Iterate over the grid's vertices.
Definition: vtkwriter.hh:422
FieldVector< DT, n > position() const
position of vertex inside the entity
Definition: vtkwriter.hh:499
int localindex() const
index of vertex within the entity, in Dune-numbering
Definition: vtkwriter.hh:494
Writer for the ouput of grid functions in the vtk format.
Definition: vtkwriter.hh:93
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:693
void clear()
clear list of registered functions
Definition: vtkwriter.hh:773
std::string write(const std::string &name, VTK::OutputType type=VTK::ascii)
write output (interface might change later)
Definition: vtkwriter.hh:801
std::string getParallelHeaderName(const std::string &name, const std::string &path, int commSize) const
return name of a parallel header file
Definition: vtkwriter.hh:918
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:711
std::string getSerialPieceName(const std::string &name, const std::string &path) const
return name of a serial piece file
Definition: vtkwriter.hh:938
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:647
void addVertexData(F &&f, VTK::FieldInfo vtkFieldInfo)
Add a function by sampling it on the grid vertices.
Definition: vtkwriter.hh:736
virtual void writeCellData(VTK::VTUWriter &writer)
write cell data
Definition: vtkwriter.hh:1331
virtual void countEntities(int &nvertices_, int &ncells_, int &ncorners_)
count the vertices, cells and corners
Definition: vtkwriter.hh:1233
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:853
virtual void writeGridCells(VTK::VTUWriter &writer)
write the connectivity array
Definition: vtkwriter.hh:1383
virtual void writeCellFaces(VTK::VTUWriter &writer)
write the connectivity array
Definition: vtkwriter.hh:1453
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:963
VTK::Precision coordPrecision() const
get the precision with which coordinates are written out
Definition: vtkwriter.hh:780
virtual void writeGridPoints(VTK::VTUWriter &writer)
write the positions of vertices
Definition: vtkwriter.hh:1359
virtual void writeVertexData(VTK::VTUWriter &writer)
write vertex data
Definition: vtkwriter.hh:1345
void addCellData(F &&f, VTK::FieldInfo vtkFieldInfo)
Add a function by sampling it on the element centers.
Definition: vtkwriter.hh:672
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:758
virtual ~VTKWriter()
destructor
Definition: vtkwriter.hh:784
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:634
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:1041
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:833
base class for data array writers
Definition: dataarraywriter.hh:54
Descriptor struct for VTK fields.
Definition: common.hh:326
@ 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:350
Dump a .vtu/.vtp files contents to a stream.
Definition: pvtuwriter.hh:60
Dump a .vtu/.vtp files contents to a stream.
Definition: vtuwriter.hh:96
DataArrayWriter * makeArrayWriter(const std::string &name, unsigned ncomps, unsigned nitems, Precision prec)
acquire a DataArrayWriter
Definition: vtuwriter.hh:378
void endCellData()
finish CellData section
Definition: vtuwriter.hh:218
void beginCells()
start section for the grid cells/PolyData lines
Definition: vtuwriter.hh:272
void endPointData()
finish PointData section
Definition: vtuwriter.hh:180
void beginCellData(const std::string &scalars="", const std::string &vectors="")
start CellData section
Definition: vtuwriter.hh:203
void beginPointData(const std::string &scalars="", const std::string &vectors="")
start PointData section
Definition: vtuwriter.hh:165
void endPoints()
finish section for the point coordinates
Definition: vtuwriter.hh:247
void endCells()
start section for the grid cells/PolyData lines
Definition: vtuwriter.hh:283
void beginPoints()
start section for the point coordinates
Definition: vtuwriter.hh:236
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:269
OutputType
How the bulk data should be stored in the file.
Definition: common.hh:41
FileType
which type of VTK file to write
Definition: common.hh:250
DataMode
Whether to produce conforming or non-conforming output.
Definition: common.hh:65
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:38
#define DUNE_THROW(E, m)
Definition: exceptions.hh:216
void localFunction(DiscreteGlobalBasisFunction< TT... > &&t)=delete
Construction of local functions from a temporary DiscreteGlobalBasisFunction (forbidden)
const IndexSet & indexSet() const
obtain the index set
Definition: gridview.hh:175
Traits::Grid Grid
type of the grid
Definition: gridview.hh:80
const CollectiveCommunication & comm() const
obtain collective communication object
Definition: gridview.hh:249
Traits::IndexSet IndexSet
type of the index set
Definition: gridview.hh:83
Grid::ctype ctype
type used for coordinates in grid
Definition: gridview.hh:127
@ dimension
The dimension of the grid.
Definition: gridview.hh:130
@ dimensionworld
The dimension of the world the grid lives in.
Definition: gridview.hh:134
MCMGLayout mcmgVertexLayout()
layout for vertices (dim-0 entities)
Definition: mcmgmapper.hh:105
std::string relativePath(const std::string &newbase, const std::string &p)
compute a relative path between two paths
Definition: path.cc:151
std::string concatPaths(const std::string &base, const std::string &p)
concatenate two paths
Definition: path.cc:30
auto min(ADLTag< 0 >, const V &v1, const V &v2)
implements binary Simd::min()
Definition: defaults.hh:87
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:11
Utilities for handling filesystem paths.
Static tag representing a codimension.
Definition: dimension.hh:22
static const ReferenceElement & general(const GeometryType &type)
get general reference elements
Definition: referenceelements.hh:196
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:162
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:186
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:204
virtual void unbind()
Unbind data set from current grid entity - mostly here for performance and symmetry reasons.
Definition: vtkwriter.hh:199
virtual void bind(const Entity &e)
Bind data set to grid entity - must be called before evaluating (i.e. calling write())
Definition: vtkwriter.hh:194
Type erasure implementation for C++ functions, i.e., functions that can be evaluated in global coordi...
Definition: vtkwriter.hh:234
virtual void unbind()
Unbind data set from current grid entity - mostly here for performance and symmetry reasons.
Definition: vtkwriter.hh:248
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:253
virtual void bind(const Entity &e)
Bind data set to grid entity - must be called before evaluating (i.e. calling write())
Definition: vtkwriter.hh:243
Type erasure implementation for legacy VTKFunctions.
Definition: vtkwriter.hh:274
virtual void unbind()
Unbind data set from current grid entity - mostly here for performance and symmetry reasons.
Definition: vtkwriter.hh:285
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:290
virtual void bind(const Entity &e)
Bind data set to grid entity - must be called before evaluating (i.e. calling write())
Definition: vtkwriter.hh:280
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:26