00001
00002
00003
00004 #ifndef DUNE_VTKWRITER_HH
00005 #define DUNE_VTKWRITER_HH
00006
00007 #include <cstring>
00008 #include <iostream>
00009 #include <fstream>
00010 #include <sstream>
00011 #include <iomanip>
00012
00013 #include <vector>
00014 #include <list>
00015
00016 #include <dune/common/deprecated.hh>
00017 #include <dune/common/exceptions.hh>
00018 #include <dune/common/iteratorfacades.hh>
00019 #include <dune/grid/common/mcmgmapper.hh>
00020 #include <dune/grid/common/referenceelements.hh>
00021
00022
00023
00024
00025
00026
00027
00038 namespace Dune
00039 {
00042 struct VTKOptions
00043 {
00044 enum OutputType {
00046 ascii,
00048 binary,
00050 binaryappended
00051 };
00052 enum DataMode {
00054 conforming,
00056 nonconforming
00057 };
00058 };
00059
00060
00061
00062 template<class T>
00063 struct VTKTypeNameTraits {
00064 std::string operator () (){
00065 return "";
00066 }
00067 };
00068
00069 template<>
00070 struct VTKTypeNameTraits<char> {
00071 std::string operator () () {
00072 return "Int8";
00073 }
00074 typedef int PrintType;
00075 };
00076
00077 template<>
00078 struct VTKTypeNameTraits<unsigned char> {
00079 std::string operator () () {
00080 return "UInt8";
00081 }
00082 typedef int PrintType;
00083 };
00084
00085 template<>
00086 struct VTKTypeNameTraits<short> {
00087 std::string operator () () {
00088 return "Int16";
00089 }
00090 typedef short PrintType;
00091 };
00092
00093 template<>
00094 struct VTKTypeNameTraits<unsigned short> {
00095 std::string operator () () {
00096 return "UInt16";
00097 }
00098 typedef unsigned short PrintType;
00099 };
00100
00101 template<>
00102 struct VTKTypeNameTraits<int> {
00103 std::string operator () () {
00104 return "Int32";
00105 }
00106 typedef int PrintType;
00107 };
00108
00109 template<>
00110 struct VTKTypeNameTraits<unsigned int> {
00111 std::string operator () () {
00112 return "UInt32";
00113 }
00114 typedef unsigned int PrintType;
00115 };
00116
00117 template<>
00118 struct VTKTypeNameTraits<float> {
00119 std::string operator () () {
00120 return "Float32";
00121 }
00122 typedef float PrintType;
00123 };
00124
00125 template<>
00126 struct VTKTypeNameTraits<double> {
00127 std::string operator () () {
00128 return "Float64";
00129 }
00130 typedef double PrintType;
00131 };
00132
00133
00142 template< class GridView >
00143 class VTKWriter {
00144 template<int dim>
00145 struct P0Layout
00146 {
00147 bool contains (Dune::GeometryType gt)
00148 {
00149 return gt.dim()==dim;
00150 }
00151 };
00152
00153 template<int dim>
00154 struct P1Layout
00155 {
00156 bool contains (Dune::GeometryType gt)
00157 {
00158 return gt.dim()==0;
00159 }
00160 };
00161
00162
00163 typedef typename GridView::Grid Grid;
00164 typedef typename Grid::ctype DT;
00165 enum { n = GridView::dimension };
00166 enum { w = GridView::dimensionworld };
00167
00168 typedef typename GridView::template Codim< 0 >::Entity Cell;
00169 typedef typename GridView::template Codim< n >::Entity Vertex;
00170 typedef Cell Entity;
00171
00172 typedef typename GridView::IndexSet IndexSet;
00173
00174 static const PartitionIteratorType VTK_Partition = InteriorBorder_Partition;
00175
00176 typedef typename GridView::template Codim< 0 >
00177 ::template Partition< VTK_Partition >::Iterator
00178 GridCellIterator;
00179 typedef typename GridView::template Codim< n >
00180 ::template Partition< VTK_Partition >::Iterator
00181 GridVertexIterator;
00182
00183 typedef MultipleCodimMultipleGeomTypeMapper< Grid, IndexSet, P1Layout > VertexMapper;
00184
00185 public:
00191 class VTKFunction
00192 {
00193 public:
00195 virtual int ncomps () const = 0;
00196
00198
00204 virtual double evaluate (int comp, const Entity& e, const Dune::FieldVector<DT,n>& xi) const = 0;
00205
00207 virtual std::string name () const = 0;
00208
00210 virtual ~VTKFunction () {}
00211 };
00212
00213 protected:
00214 typedef typename std::list<VTKFunction*>::iterator FunctionIterator;
00215
00216 class CellIterator :
00217 public ForwardIteratorFacade<CellIterator, Entity, Entity&, int>
00218 {
00219 GridCellIterator git;
00220 GridCellIterator gend;
00221 public:
00222 CellIterator(const GridCellIterator & x, const GridCellIterator & end) : git(x), gend(end) {};
00223 void increment ()
00224 {
00225 ++git;
00226 while (git!=gend && git->partitionType()!=InteriorEntity) ++git;
00227 }
00228 bool equals (const CellIterator & cit) const
00229 {
00230 return git == cit.git;
00231 }
00232 Entity& dereference() const
00233 {
00234 return *git;
00235 }
00236 const FieldVector<DT,n> position() const
00237 {
00238 return ReferenceElements<DT,n>::general(git->type()).position(0,0);
00239 }
00240 };
00241
00242 CellIterator cellBegin() const
00243 {
00244 return CellIterator( gridView_.template begin< 0, VTK_Partition >(),
00245 gridView_.template end< 0, VTK_Partition >() );
00246 }
00247
00248 CellIterator cellEnd() const
00249 {
00250 return CellIterator( gridView_.template end< 0, VTK_Partition >(),
00251 gridView_.template end< 0, VTK_Partition >() );
00252 }
00253
00254 class VertexIterator :
00255 public ForwardIteratorFacade<VertexIterator, Entity, Entity&, int>
00256 {
00257 GridCellIterator git;
00258 GridCellIterator gend;
00259 VTKOptions::DataMode datamode;
00260 int index;
00261 const VertexMapper & vertexmapper;
00262 std::vector<bool> visited;
00263 const std::vector<int> & number;
00264 int offset;
00265 protected:
00266 void basicIncrement ()
00267 {
00268 if( git == gend )
00269 return;
00270 ++index;
00271 const int numCorners = git->template count< n >();
00272 if( index == numCorners )
00273 {
00274 offset += numCorners;
00275 index = 0;
00276
00277 ++git;
00278 while( (git != gend) && (git->partitionType() != InteriorEntity) )
00279 ++git;
00280 }
00281 }
00282 public:
00283 VertexIterator(const GridCellIterator & x,
00284 const GridCellIterator & end,
00285 const VTKOptions::DataMode & dm,
00286 const VertexMapper & vm,
00287 const std::vector<int> & num) :
00288 git(x), gend(end), datamode(dm), index(0),
00289 vertexmapper(vm), visited(vm.size(), false),
00290 number(num), offset(0)
00291 {
00292 if (datamode == VTKOptions::conforming && git != gend)
00293 visited[vertexmapper.template map<n>(*git,index)] = true;
00294 };
00295 void increment ()
00296 {
00297 switch (datamode)
00298 {
00299 case VTKOptions::conforming:
00300 while(visited[vertexmapper.template map<n>(*git,index)])
00301 {
00302 basicIncrement();
00303 if (git == gend) return;
00304 }
00305 visited[vertexmapper.template map<n>(*git,index)] = true;
00306 break;
00307 case VTKOptions::nonconforming:
00308 basicIncrement();
00309 break;
00310 }
00311 }
00312 bool equals (const VertexIterator & cit) const
00313 {
00314 return git == cit.git
00315 && index == cit.index && datamode == cit.datamode;
00316 }
00317 Entity& dereference() const
00318 {
00319 return *git;
00320 }
00321 int id () const
00322 {
00323 switch (datamode)
00324 {
00325 case VTKOptions::conforming:
00326 return
00327 number[vertexmapper.template map<n>(*git,renumber(*git,index))];
00328 case VTKOptions::nonconforming:
00329 return offset + renumber(*git,index);
00330 default:
00331 DUNE_THROW(IOError,"VTKWriter: unsupported DataMode" << datamode);
00332 }
00333 }
00334 int localindex () const
00335 {
00336 return index;
00337 }
00338 const FieldVector<DT,n> & position () const
00339 {
00340 return ReferenceElements<DT,n>::general(git->type()).position(index,n);
00341 }
00342 };
00343
00344 VertexIterator vertexBegin () const
00345 {
00346 return VertexIterator( gridView_.template begin< 0, VTK_Partition >(),
00347 gridView_.template end< 0, VTK_Partition >(),
00348 datamode, *vertexmapper, number );
00349 }
00350
00351 VertexIterator vertexEnd () const
00352 {
00353 return VertexIterator( gridView_.template end< 0, VTK_Partition >(),
00354 gridView_.template end< 0, VTK_Partition >(),
00355 datamode, *vertexmapper, number );
00356 }
00357
00358 class CornerIterator :
00359 public ForwardIteratorFacade<CornerIterator, Entity, Entity&, int>
00360 {
00361 GridCellIterator git;
00362 GridCellIterator gend;
00363 VTKOptions::DataMode datamode;
00364 int index;
00365 const VertexMapper & vertexmapper;
00366 std::vector<bool> visited;
00367 const std::vector<int> & number;
00368 int offset;
00369 protected:
00370 void basicIncrement ()
00371 {
00372 if( git == gend )
00373 return;
00374 ++index;
00375 const int numCorners = git->template count< n >();
00376 if( index == numCorners )
00377 {
00378 offset += numCorners;
00379 index = 0;
00380
00381 ++git;
00382 while( (git != gend) && (git->partitionType() != InteriorEntity) )
00383 ++git;
00384 }
00385 }
00386 public:
00387 CornerIterator(const GridCellIterator & x,
00388 const GridCellIterator & end,
00389 const VTKOptions::DataMode & dm,
00390 const VertexMapper & vm,
00391 const std::vector<int> & num) :
00392 git(x), gend(end), datamode(dm), index(0),
00393 vertexmapper(vm),
00394 number(num), offset(0) {};
00395 void increment ()
00396 {
00397 basicIncrement();
00398 }
00399 bool equals (const CornerIterator & cit) const
00400 {
00401 return git == cit.git
00402 && index == cit.index && datamode == cit.datamode;
00403 }
00404 Entity& dereference() const
00405 {
00406 return *git;
00407 }
00408 int id () const
00409 {
00410 switch (datamode)
00411 {
00412 case VTKOptions::conforming:
00413 return
00414 number[vertexmapper.template map<n>(*git,renumber(*git,index))];
00415 case VTKOptions::nonconforming:
00416 return offset + renumber(*git,index);
00417 default:
00418 DUNE_THROW(IOError,"VTKWriter: unsupported DataMode" << datamode);
00419 }
00420 }
00421 int localindex () const
00422 {
00423 return index;
00424 }
00425 };
00426
00427 CornerIterator cornerBegin () const
00428 {
00429 return CornerIterator( gridView_.template begin< 0, VTK_Partition >(),
00430 gridView_.template end< 0, VTK_Partition >(),
00431 datamode, *vertexmapper, number );
00432 }
00433
00434 CornerIterator cornerEnd () const
00435 {
00436 return CornerIterator( gridView_.template end< 0, VTK_Partition >(),
00437 gridView_.template end< 0, VTK_Partition >(),
00438 datamode, *vertexmapper, number );
00439 }
00440
00441 private:
00445 template<class V>
00446 class P0VectorWrapper : public VTKFunction
00447 {
00448 typedef MultipleCodimMultipleGeomTypeMapper< Grid, IndexSet, P0Layout > VM0;
00449 public:
00451 virtual int ncomps () const
00452 {
00453 return 1;
00454 }
00455
00457 virtual double evaluate (int comp, const Entity& e, const Dune::FieldVector<DT,n>& xi) const
00458 {
00459 return v[mapper.map(e)];
00460 }
00461
00463 virtual std::string name () const
00464 {
00465 return s;
00466 }
00467
00469 P0VectorWrapper ( const Grid &g_, const IndexSet &is_, const V &v_, std::string s_)
00470 : g(g_), is(is_), v(v_), s(s_), mapper(g_,is_)
00471 {
00472 if (v.size()!=(unsigned int)mapper.size())
00473 DUNE_THROW(IOError,"VTKWriter::P0VectorWrapper: size mismatch");
00474 }
00475
00476 virtual ~P0VectorWrapper() {}
00477
00478 private:
00479 const Grid& g;
00480 const IndexSet &is;
00481 const V& v;
00482 std::string s;
00483 VM0 mapper;
00484 };
00485
00489 template<class V>
00490 class P1VectorWrapper : public VTKFunction
00491 {
00492 typedef MultipleCodimMultipleGeomTypeMapper< Grid, IndexSet, P1Layout > VM1;
00493
00494 public:
00496 virtual int ncomps () const
00497 {
00498 return 1;
00499 }
00500
00502 virtual double evaluate (int comp, const Entity& e, const Dune::FieldVector<DT,n>& xi) const
00503 {
00504 double min=1E100;
00505 int imin=-1;
00506 Dune::GeometryType gt = e.type();
00507 for (int i=0; i<e.template count<n>(); ++i)
00508 {
00509 Dune::FieldVector<DT,n>
00510 local = Dune::ReferenceElements<DT,n>::general(gt).position(i,n);
00511 local -= xi;
00512 if (local.infinity_norm()<min)
00513 {
00514 min = local.infinity_norm();
00515 imin = i;
00516 }
00517 }
00518 return v[mapper.template map<n>(e,imin)];
00519 }
00520
00522 virtual std::string name () const
00523 {
00524 return s;
00525 }
00526
00528 P1VectorWrapper ( const Grid &g_, const IndexSet &is_, const V &v_, std::string s_ )
00529 : g(g_), is(is_), v(v_), s(s_), mapper(g_,is_)
00530 {
00531 if (v.size()!=(unsigned int)mapper.size())
00532 DUNE_THROW(IOError,"VTKWriter::P1VectorWrapper: size mismatch");
00533 }
00534
00535 virtual ~P1VectorWrapper() {}
00536
00537 private:
00538 const Grid& g;
00539 const IndexSet &is;
00540 const V& v;
00541 std::string s;
00542 VM1 mapper;
00543 };
00544
00545 public:
00553 explicit VTKWriter ( const GridView &gridView,
00554 VTKOptions::DataMode dm = VTKOptions::conforming )
00555 : gridView_( gridView ),
00556 grid( gridView.grid() ),
00557 is( gridView_.indexSet() ),
00558 datamode( dm )
00559 {
00560 indentCount = 0;
00561 numPerLine = 4*3;
00562 }
00563
00569 void addCellData (VTKFunction* p)
00570 {
00571 celldata.push_back(p);
00572 }
00573
00585 template<class V>
00586 void addCellData (const V& v, std::string name)
00587 {
00588 VTKFunction* p = new P0VectorWrapper<V>(grid,is,v,name);
00589 celldata.push_back(p);
00590 }
00591
00597 void addVertexData (VTKFunction* p)
00598 {
00599 vertexdata.push_back(p);
00600 }
00601
00613 template<class V>
00614 void addVertexData (const V& v, std::string name)
00615 {
00616 VTKFunction* p = new P1VectorWrapper<V>(grid,is,v,name);
00617 vertexdata.push_back(p);
00618 }
00619
00621 void clear ()
00622 {
00623 for (FunctionIterator it=celldata.begin();
00624 it!=celldata.end(); ++it)
00625 delete *it;
00626 celldata.clear();
00627 for (FunctionIterator it=vertexdata.begin();
00628 it!=vertexdata.end(); ++it)
00629 delete *it;
00630 vertexdata.clear();
00631 }
00632
00634 virtual ~VTKWriter ()
00635 {
00636 this->clear();
00637 }
00638
00644 std::string write ( const std::string &name,
00645 VTKOptions::OutputType type = VTKOptions::ascii )
00646 {
00647 return write( name, type, gridView_.comm().rank(), gridView_.comm().size() );
00648 }
00649
00657 std::string pwrite ( const char* name, const char* path, const char* extendpath,
00658 VTKOptions::OutputType type = VTKOptions::ascii )
00659 {
00660 return pwrite( name, path, extendpath, type, gridView_.comm().rank(), gridView_.comm().size() );
00661 }
00662
00663 protected:
00664 std::string write ( const std::string &name,
00665 VTKOptions::OutputType type,
00666 const int commRank,
00667 const int commSize )
00668 {
00669
00670 outputtype = type;
00671
00672
00673 bytecount = 0;
00674
00675
00676 std::ostringstream pieceName;
00677 if( commSize > 1 )
00678 {
00679 pieceName << "s" << std::setfill( '0' ) << std::setw( 4 ) << commSize << ":";
00680 pieceName << "p" << std::setfill( '0' ) << std::setw( 4 ) << commRank << ":";
00681 }
00682 pieceName << name << (GridView::dimension > 1 ? ".vtu" : ".vtp");
00683
00684
00685 std::ofstream file;
00686 if( outputtype == VTKOptions::binaryappended )
00687 file.open( pieceName.str().c_str(), std::ios::binary );
00688 else
00689 file.open( pieceName.str().c_str() );
00690 writeDataFile( file );
00691 file.close();
00692
00693
00694 if( commSize == 1 )
00695 return pieceName.str();
00696
00697
00698 gridView_.comm().barrier();
00699
00700
00701 std::ostringstream parallelName;
00702 parallelName << "s" << std::setfill( '0' ) << std::setw( 4 ) << commSize << ":";
00703 parallelName << name << (GridView::dimension > 1 ? ".pvtu" : ".pvtp");
00704
00705
00706 if( commRank == 0 )
00707 {
00708 file.open( parallelName.str().c_str() );
00709 writeParallelHeader( file, name.c_str(), "", commSize );
00710 file.close();
00711 }
00712
00713
00714 gridView_.comm().barrier();
00715 return parallelName.str();
00716 }
00717
00719 std::string pwrite ( const char* name, const char* path, const char* extendpath,
00720 VTKOptions::OutputType ot,
00721 const int commRank,
00722 const int commSize )
00723 {
00724
00725 outputtype=ot;
00726
00727
00728 bytecount = 0;
00729
00730
00731 std::ofstream file;
00732 char piecepath[ 4096 ];
00733 char relpiecepath[ 4096 ];
00734 int n=strlen(path);
00735 int m=strlen(extendpath);
00736 if (n>0 && path[0]=='/' && path[n-1]=='/')
00737 {
00738
00739
00740
00741 if (m==0)
00742 {
00743
00744 piecepath[0] = '/';
00745 piecepath[1] = '\0';
00746 }
00747 else
00748 {
00749
00750 char *p=piecepath;
00751 if (extendpath[0]!='/')
00752 {
00753 *p = '/';
00754 p++;
00755 }
00756 for (int i=0; i<m; i++)
00757 {
00758 *p = extendpath[i];
00759 p++;
00760 }
00761 if (*(p-1)!='/')
00762 {
00763 *p = '/';
00764 p++;
00765 }
00766 *p = '\0';
00767 }
00768
00769
00770 int k=0;
00771 const char *p=path;
00772 while (*p!='\0')
00773 {
00774 if (*p=='/') k++;
00775 p++;
00776 }
00777 char *pp = relpiecepath;
00778 if (k>1)
00779 {
00780 for (int i=0; i<k; i++)
00781 {
00782 *pp='.'; pp++; *pp='.'; pp++; *pp='/'; pp++;
00783 }
00784 }
00785
00786 for (int i=0; i<m; i++)
00787 {
00788 if (i==0 && extendpath[i]=='/') continue;
00789 *pp = extendpath[i];
00790 pp++;
00791 }
00792 if ( pp!=relpiecepath && (*(pp-1)!='/') )
00793 {
00794 *pp = '/';
00795 pp++;
00796 }
00797 *pp = '\0';
00798 }
00799 else
00800 {
00801
00802
00803 if (n==0 || m==0)
00804 sprintf(piecepath,"%s%s",path,extendpath);
00805 else
00806 {
00807
00808 if (path[n-1]!='/' && extendpath[0]!='/')
00809 sprintf(piecepath,"%s/%s",path,extendpath);
00810 else
00811 sprintf(piecepath,"%s%s",path,extendpath);
00812 }
00813
00814 sprintf(relpiecepath,"%s",extendpath);
00815 }
00816 char fullname[ 8192 ];
00817 if (GridView::dimension>1)
00818 sprintf(fullname,"%s/s%04d:p%04d:%s.vtu",piecepath, commSize, commRank, name);
00819 else
00820 sprintf(fullname,"%s/s%04d:p%04d:%s.vtp",piecepath, commSize, commRank, name);
00821 if (outputtype==VTKOptions::binaryappended)
00822 file.open(fullname,std::ios::binary);
00823 else
00824 file.open(fullname);
00825 writeDataFile(file);
00826 file.close();
00827 gridView_.comm().barrier();
00828 if( commRank ==0 )
00829 {
00830 if (GridView::dimension>1)
00831 sprintf(fullname,"%s/s%04d:%s.pvtu",path, commSize, name);
00832 else
00833 sprintf(fullname,"%s/s%04d:%s.pvtp",path, commSize, name);
00834 file.open(fullname);
00835 writeParallelHeader(file,name,relpiecepath, commSize );
00836 file.close();
00837 }
00838 gridView_.comm().barrier();
00839 return fullname;
00840 }
00841
00842 private:
00843 enum VTKGeometryType
00844 {
00845 vtkLine = 3,
00846 vtkTriangle = 5,
00847 vtkQuadrilateral = 9,
00848 vtkTetrahedron = 10,
00849 vtkHexahedron = 12,
00850 vtkPrism = 13,
00851 vtkPyramid = 14
00852 };
00853
00854 protected:
00856 static VTKGeometryType vtkType(const Dune::GeometryType & t)
00857 {
00858 if (t.isLine())
00859 return vtkLine;
00860 if (t.isTriangle())
00861 return vtkTriangle;
00862 if (t.isQuadrilateral())
00863 return vtkQuadrilateral;
00864 if (t.isTetrahedron())
00865 return vtkTetrahedron;
00866 if (t.isPyramid())
00867 return vtkPyramid;
00868 if (t.isPrism())
00869 return vtkPrism;
00870 if (t.isHexahedron())
00871 return vtkHexahedron;
00872 DUNE_THROW(IOError,"VTKWriter: unsupported GeometryType " << t);
00873 }
00874
00875 private:
00877 void writeParallelHeader ( std::ostream& s, const char* piecename, const char* piecepath,
00878 const int commSize )
00879 {
00880
00881 s << "<?xml version=\"1.0\"?>" << std::endl;
00882
00883
00884 if (n>1)
00885 s << "<VTKFile type=\"PUnstructuredGrid\" version=\"0.1\" byte_order=\"LittleEndian\">" << std::endl;
00886 else
00887 s << "<VTKFile type=\"PPolyData\" version=\"0.1\" byte_order=\"LittleEndian\">" << std::endl;
00888 indentUp();
00889
00890
00891 indent(s);
00892 if (n>1)
00893 s << "<PUnstructuredGrid GhostLevel=\"0\">" << std::endl;
00894 else
00895 s << "<PPolyData GhostLevel=\"0\">" << std::endl;
00896 indentUp();
00897
00898
00899 indent(s); s << "<PPointData";
00900 for (FunctionIterator it=vertexdata.begin(); it!=vertexdata.end(); ++it)
00901 if ((*it)->ncomps()==1)
00902 {
00903 s << " Scalars=\"" << (*it)->name() << "\"" ;
00904 break;
00905 }
00906 for (FunctionIterator it=vertexdata.begin(); it!=vertexdata.end(); ++it)
00907 if ((*it)->ncomps()>1)
00908 {
00909 s << " Vectors=\"" << (*it)->name() << "\"" ;
00910 break;
00911 }
00912 s << ">" << std::endl;
00913 indentUp();
00914 for (FunctionIterator it=vertexdata.begin(); it!=vertexdata.end(); ++it)
00915 {
00916 indent(s); s << "<PDataArray type=\"Float32\" Name=\"" << (*it)->name() << "\" ";
00917 s << "NumberOfComponents=\"" << ((*it)->ncomps()>1?3:1) << "\" ";
00918 if (outputtype==VTKOptions::ascii)
00919 s << "format=\"ascii\"/>" << std::endl;
00920 if (outputtype==VTKOptions::binary)
00921 s << "format=\"binary\"/>" << std::endl;
00922 if (outputtype==VTKOptions::binaryappended)
00923 s << "format=\"appended\"/>" << std::endl;
00924 }
00925 indentDown();
00926 indent(s); s << "</PPointData>" << std::endl;
00927
00928
00929 indent(s); s << "<PCellData";
00930 for (FunctionIterator it=celldata.begin(); it!=celldata.end(); ++it)
00931 if ((*it)->ncomps()==1)
00932 {
00933 s << " Scalars=\"" << (*it)->name() << "\"" ;
00934 break;
00935 }
00936 for (FunctionIterator it=celldata.begin(); it!=celldata.end(); ++it)
00937 if ((*it)->ncomps()>1)
00938 {
00939 s << " Vectors=\"" << (*it)->name() << "\"" ;
00940 break;
00941 }
00942 s << ">" << std::endl;
00943 indentUp();
00944 for (FunctionIterator it=celldata.begin(); it!=celldata.end(); ++it)
00945 {
00946 indent(s); s << "<PDataArray type=\"Float32\" Name=\"" << (*it)->name() << "\" ";
00947 s << "NumberOfComponents=\"" << ((*it)->ncomps()>1?3:1) << "\" ";
00948 if (outputtype==VTKOptions::ascii)
00949 s << "format=\"ascii\"/>" << std::endl;
00950 if (outputtype==VTKOptions::binary)
00951 s << "format=\"binary\"/>" << std::endl;
00952 if (outputtype==VTKOptions::binaryappended)
00953 s << "format=\"appended\"/>" << std::endl;
00954 }
00955 indentDown();
00956 indent(s); s << "</PCellData>" << std::endl;
00957
00958
00959 indent(s); s << "<PPoints>" << std::endl;
00960 indentUp();
00961 indent(s); s << "<PDataArray type=\"Float32\" Name=\"Coordinates\" NumberOfComponents=\"" << "3" << "\" ";
00962 if (outputtype==VTKOptions::ascii)
00963 s << "format=\"ascii\"/>" << std::endl;
00964 if (outputtype==VTKOptions::binary)
00965 s << "format=\"binary\"/>" << std::endl;
00966 if (outputtype==VTKOptions::binaryappended)
00967 s << "format=\"appended\"/>" << std::endl;
00968 indentDown();
00969 indent(s); s << "</PPoints>" << std::endl;
00970
00971
00972 for( int i = 0; i < commSize; ++i )
00973 {
00974 char fullname[ 4096 ];
00975 if (GridView::dimension>1)
00976 sprintf(fullname,"%s/s%04d:p%04d:%s.vtu",piecepath, commSize, i,piecename);
00977 else
00978 sprintf(fullname,"%s/s%04d:p%04d:%s.vtp",piecepath, commSize, i,piecename);
00979 indent(s); s << "<Piece Source=\"" << fullname << "\"/>" << std::endl;
00980 }
00981
00982
00983 indentDown();
00984 indent(s);
00985 if (n>1)
00986 s << "</PUnstructuredGrid>" << std::endl;
00987 else
00988 s << "</PPolyData>" << std::endl;
00989
00990
00991 indentDown();
00992 s << "</VTKFile>" << std::endl;
00993 }
00994
00996 void writeDataFile (std::ostream& s)
00997 {
00998
00999 s << "<?xml version=\"1.0\"?>" << std::endl;
01000
01001
01002 if (n>1)
01003 s << "<VTKFile type=\"UnstructuredGrid\" version=\"0.1\" byte_order=\"LittleEndian\">" << std::endl;
01004 else
01005 s << "<VTKFile type=\"PolyData\" version=\"0.1\" byte_order=\"LittleEndian\">" << std::endl;
01006 indentUp();
01007
01008
01009 vertexmapper = new VertexMapper(grid,is);
01010 if (datamode == VTKOptions::conforming)
01011 {
01012 number.resize(vertexmapper->size());
01013 for (std::vector<int>::size_type i=0; i<number.size(); i++) number[i] = -1;
01014 }
01015 countEntities(nvertices, ncells, ncorners);
01016
01017
01018 indent(s);
01019 if (n>1)
01020 s << "<UnstructuredGrid>" << std::endl;
01021 else
01022 s << "<PolyData>" << std::endl;
01023 indentUp();
01024
01025
01026 indent(s);
01027 if (n>1)
01028 s << "<Piece NumberOfPoints=\"" << nvertices << "\" NumberOfCells=\"" << ncells << "\">" << std::endl;
01029 else
01030 s << "<Piece NumberOfPoints=\"" << nvertices << "\""
01031 << " NumberOfVerts=\"0\""
01032 << " NumberOfLines=\"" << ncells << "\">"
01033 << " NumberOfPolys=\"0\"" << std::endl;
01034 indentUp();
01035
01036
01037 writeVertexData(s);
01038
01039
01040 writeCellData(s);
01041
01042
01043 writeGridPoints(s);
01044
01045
01046 writeGridCells(s);
01047
01048
01049 indentDown();
01050 indent(s); s << "</Piece>" << std::endl;
01051
01052
01053 indentDown();
01054 indent(s);
01055 if (n>1)
01056 s << "</UnstructuredGrid>" << std::endl;
01057 else
01058 s << "</PolyData>" << std::endl;
01059
01060
01061 if (outputtype==VTKOptions::binaryappended)
01062 writeAppendedData(s);
01063
01064
01065 indentDown();
01066 s << "</VTKFile>" << std::endl;
01067
01068 delete vertexmapper; number.clear();
01069 }
01070
01071 protected:
01073 virtual void countEntities(int &nvertices, int &ncells, int &ncorners)
01074 {
01075 nvertices = 0;
01076 ncells = 0;
01077 ncorners = 0;
01078 for (CellIterator it=cellBegin(); it!=cellEnd(); ++it)
01079 {
01080 ncells++;
01081 for (int i=0; i<it->template count<n>(); ++i)
01082 {
01083 ncorners++;
01084 if (datamode == VTKOptions::conforming)
01085 {
01086 int alpha = vertexmapper->template map<n>(*it,i);
01087 if (number[alpha]<0)
01088 number[alpha] = nvertices++;
01089 }
01090 else
01091 {
01092 nvertices++;
01093 }
01094 }
01095 }
01096 }
01097
01099 virtual void writeCellData (std::ostream& s)
01100 {
01101 indent(s); s << "<CellData";
01102 for (FunctionIterator it=celldata.begin(); it!=celldata.end(); ++it)
01103 if ((*it)->ncomps()==1)
01104 {
01105 s << " Scalars=\"" << (*it)->name() << "\"" ;
01106 break;
01107 }
01108 for (FunctionIterator it=celldata.begin(); it!=celldata.end(); ++it)
01109 if ((*it)->ncomps()>1)
01110 {
01111 s << " Vectors=\"" << (*it)->name() << "\"" ;
01112 break;
01113 }
01114 s << ">" << std::endl;
01115 indentUp();
01116 for (FunctionIterator it=celldata.begin(); it!=celldata.end(); ++it)
01117 {
01118 VTKDataArrayWriter<float> *p = makeVTKDataArrayWriter<float>(s, (*it)->name().c_str(), (*it)->ncomps(), (*it)->ncomps()*ncells);
01119 for (CellIterator i=cellBegin(); i!=cellEnd(); ++i)
01120 {
01121 for (int j=0; j<(*it)->ncomps(); j++)
01122 p->write((*it)->evaluate(j,*i,i.position()));
01123
01124 if((*it)->ncomps()==2)
01125 p->write(0.0);
01126 }
01127 delete p;
01128 }
01129 indentDown();
01130 indent(s); s << "</CellData>" << std::endl;
01131 }
01132
01134 virtual void writeVertexData (std::ostream& s)
01135 {
01136 indent(s); s << "<PointData";
01137 for (FunctionIterator it=vertexdata.begin(); it!=vertexdata.end(); ++it)
01138 if ((*it)->ncomps()==1)
01139 {
01140 s << " Scalars=\"" << (*it)->name() << "\"" ;
01141 break;
01142 }
01143 for (FunctionIterator it=vertexdata.begin(); it!=vertexdata.end(); ++it)
01144 if ((*it)->ncomps()>1)
01145 {
01146 s << " Vectors=\"" << (*it)->name() << "\"" ;
01147 break;
01148 }
01149 s << ">" << std::endl;
01150 indentUp();
01151 for (FunctionIterator it=vertexdata.begin(); it!=vertexdata.end(); ++it)
01152 {
01153 VTKDataArrayWriter<float> *p = makeVTKDataArrayWriter<float>(s, (*it)->name().c_str(), (*it)->ncomps(), (*it)->ncomps()*nvertices);
01154 for (VertexIterator vit=vertexBegin(); vit!=vertexEnd(); ++vit)
01155 {
01156 for (int j=0; j<(*it)->ncomps(); j++)
01157 p->write((*it)->evaluate(j,*vit,vit.position()));
01158
01159 if((*it)->ncomps()==2)
01160 p->write(0.0);
01161 }
01162 delete p;
01163 }
01164 indentDown();
01165 indent(s); s << "</PointData>" << std::endl;
01166 }
01167
01169 virtual void writeGridPoints (std::ostream& s)
01170 {
01171 indent(s); s << "<Points>" << std::endl;
01172 indentUp();
01173
01174 VTKDataArrayWriter<float> *p = makeVTKDataArrayWriter<float>(s, "Coordinates", 3, 3*nvertices);
01175 VertexIterator vEnd = vertexEnd();
01176 for (VertexIterator vit=vertexBegin(); vit!=vEnd; ++vit)
01177 {
01178 int dimw=w;
01179 for (int j=0; j<std::min(dimw,3); j++)
01180 p->write(vit->geometry().corner(vit.localindex())[j]);
01181 for (int j=std::min(dimw,3); j<3; j++)
01182 p->write(0.0);
01183 }
01184 delete p;
01185
01186 indentDown();
01187 indent(s); s << "</Points>" << std::endl;
01188 }
01189
01191 virtual void writeGridCells (std::ostream& s)
01192 {
01193 indent(s);
01194 if (n>1)
01195 s << "<Cells>" << std::endl;
01196 else
01197 s << "<Lines>" << std::endl;
01198 indentUp();
01199
01200
01201 VTKDataArrayWriter<int> *p1 = makeVTKDataArrayWriter<int>(s, "connectivity", 1, ncorners);
01202 for (CornerIterator it=cornerBegin(); it!=cornerEnd(); ++it)
01203 p1->write(it.id());
01204 delete p1;
01205
01206
01207 VTKDataArrayWriter<int> *p2 = makeVTKDataArrayWriter<int>(s, "offsets", 1, ncells);
01208 {
01209 int offset = 0;
01210 for (CellIterator it=cellBegin(); it!=cellEnd(); ++it)
01211 {
01212 offset += it->template count<n>();
01213 p2->write(offset);
01214 }
01215 }
01216 delete p2;
01217
01218
01219 if (n>1)
01220 {
01221 VTKDataArrayWriter<unsigned char> *p3 = makeVTKDataArrayWriter<unsigned char>(s, "types", 1, ncells);
01222 for (CellIterator it=cellBegin(); it!=cellEnd(); ++it)
01223 {
01224 int vtktype = vtkType(it->type());
01225 p3->write(vtktype);
01226 }
01227 delete p3;
01228 }
01229
01230 indentDown();
01231 indent(s);
01232 if (n>1)
01233 s << "</Cells>" << std::endl;
01234 else
01235 s << "</Lines>" << std::endl;
01236 }
01237
01239 virtual void writeAppendedData (std::ostream& s)
01240 {
01241 indent(s); s << "<AppendedData encoding=\"raw\">" << std::endl;
01242 indentUp();
01243 indent(s); s << "_";
01244
01245 SimpleStream stream(s);
01246
01247
01248 unsigned int blocklength;
01249
01250
01251 for (FunctionIterator it=vertexdata.begin(); it!=vertexdata.end(); ++it)
01252 {
01253
01254 blocklength = nvertices * (*it)->ncomps() * sizeof(float);
01255
01256 if((*it)->ncomps()==2)
01257 blocklength = nvertices * (3) * sizeof(float);
01258 stream.write(blocklength);
01259 std::vector<bool> visited(vertexmapper->size(), false);
01260 for (VertexIterator vit=vertexBegin(); vit!=vertexEnd(); ++vit)
01261 {
01262 for (int j=0; j<(*it)->ncomps(); j++)
01263 {
01264 float data = (*it)->evaluate(j,*vit,vit.position());
01265 stream.write(data);
01266 }
01267
01268 if((*it)->ncomps()==2)
01269 {
01270 float data=0.0;
01271 stream.write(data);
01272 }
01273 }
01274 }
01275
01276
01277 for (FunctionIterator it=celldata.begin(); it!=celldata.end(); ++it)
01278 {
01279 blocklength = ncells * (*it)->ncomps() * sizeof(float);
01280 stream.write(blocklength);
01281 for (CellIterator i=cellBegin(); i!=cellEnd(); ++i)
01282 {
01283 for (int j=0; j<(*it)->ncomps(); j++)
01284 {
01285 float data = (*it)->evaluate(j,*i,i.position());
01286 stream.write(data);
01287 }
01288
01289 if((*it)->ncomps()==2)
01290 {
01291 float data=0.0;
01292 stream.write(data);
01293 }
01294 }
01295 }
01296
01297
01298 blocklength = nvertices * 3 * sizeof(float);
01299 stream.write(blocklength);
01300 std::vector<bool> visited(vertexmapper->size(), false);
01301 for (VertexIterator vit=vertexBegin(); vit!=vertexEnd(); ++vit)
01302 {
01303 int dimw=w;
01304 float data;
01305 for (int j=0; j<std::min(dimw,3); j++)
01306 {
01307 data = vit->geometry().corner(vit.localindex())[j];
01308 stream.write(data);
01309 }
01310 data = 0;
01311 for (int j=std::min(dimw,3); j<3; j++)
01312 stream.write(data);
01313 }
01314
01315
01316 blocklength = ncorners * sizeof(unsigned int);
01317 stream.write(blocklength);
01318 for (CornerIterator it=cornerBegin(); it!=cornerEnd(); ++it)
01319 {
01320 stream.write(it.id());
01321 }
01322
01323
01324 blocklength = ncells * sizeof(unsigned int);
01325 stream.write(blocklength);
01326 {
01327 int offset = 0;
01328 for (CellIterator it=cellBegin(); it!=cellEnd(); ++it)
01329 {
01330 offset += it->template count<n>();
01331 stream.write(offset);
01332 }
01333 }
01334
01335
01336 if (n>1)
01337 {
01338 blocklength = ncells * sizeof(unsigned char);
01339 stream.write(blocklength);
01340 for (CellIterator it=cellBegin(); it!=cellEnd(); ++it)
01341 {
01342 unsigned char vtktype = vtkType(it->type());
01343 stream.write(vtktype);
01344 }
01345 }
01346
01347 s << std::endl;
01348 indentDown();
01349 indent(s); s << "</AppendedData>" << std::endl;
01350 }
01351
01353 template<class T>
01354 class VTKDataArrayWriter
01355 {
01356 public:
01358 virtual void write (T data) = 0;
01360 virtual ~VTKDataArrayWriter () {}
01361 };
01362
01363 private:
01365 template<class T>
01366 class VTKAsciiDataArrayWriter : public VTKDataArrayWriter<T>
01367 {
01368 public:
01370 VTKAsciiDataArrayWriter (std::ostream& theStream, std::string name, int ncomps)
01371 : s(theStream), counter(0), numPerLine(12)
01372 {
01373 VTKTypeNameTraits<T> tn;
01374 s << "<DataArray type=\"" << tn() << "\" Name=\"" << name << "\" ";
01375
01376 if (ncomps>3)
01377 DUNE_THROW(IOError, "VTKWriter does not support more than 3 components");
01378 s << "NumberOfComponents=\"" << (ncomps>1?3:1) << "\" ";
01379 s << "format=\"ascii\">" << std::endl;
01380 }
01381
01383 void write (T data)
01384 {
01385 typedef typename VTKTypeNameTraits<T>::PrintType PT;
01386 s << (PT) data << " ";
01387 counter++;
01388 if (counter%numPerLine==0) s << std::endl;
01389 }
01390
01392 ~VTKAsciiDataArrayWriter ()
01393 {
01394 if (counter%numPerLine!=0) s << std::endl;
01395 s << "</DataArray>" << std::endl;
01396 }
01397
01398 private:
01399 std::ostream& s;
01400 int counter;
01401 int numPerLine;
01402 };
01403
01404
01405 template<class T>
01406 class VTKBinaryDataArrayWriter : public VTKDataArrayWriter<T>
01407 {
01408 public:
01410 VTKBinaryDataArrayWriter (std::ostream& theStream, std::string name, int ncomps, int nitems)
01411 : s(theStream),bufsize(4096),n(0)
01412 {
01413 DUNE_THROW(IOError, "binary does not work yet, use binaryappended!");
01414 VTKTypeNameTraits<T> tn;
01415 s << "<DataArray type=\"" << tn() << "\" Name=\"" << name << "\" ";
01416
01417 if (ncomps>3)
01418 DUNE_THROW(IOError, "VTKWriter does not support more than 3 components");
01419 s << "NumberOfComponents=\"" << (ncomps>1?3:1) << "\" ";
01420 s << "format=\"binary\">" << std::endl;
01421 buffer = new char[bufsize*sizeof(T)];
01422 code = new char[2*bufsize*sizeof(T)];
01423 unsigned int size = nitems*sizeof(T);
01424 char* p = reinterpret_cast<char*>(&size);
01425 memcpy(buffer+n,p,sizeof(int));
01426 n += sizeof(int);
01427
01428 }
01429
01431 void write (T data)
01432 {
01433 if (n+sizeof(T)>bufsize)
01434 {
01435
01436
01437
01438 n=0;
01439 }
01440 char* p = reinterpret_cast<char*>(&data);
01441 memcpy(buffer+n,p,sizeof(T));
01442 n += sizeof(T);
01443 }
01444
01446 ~VTKBinaryDataArrayWriter ()
01447 {
01448
01449 if (n>0)
01450 {
01451
01452
01453 }
01454
01455
01456
01457 s << std::endl;
01458 s << "</DataArray>" << std::endl;
01459 delete [] code;
01460 delete [] buffer;
01461 }
01462
01463 private:
01464 std::ostream& s;
01465
01466 size_t bufsize;
01467 char* buffer;
01468 char* code;
01469 int n;
01470 };
01471
01473 template<class T>
01474 class VTKBinaryAppendedDataArrayWriter : public VTKDataArrayWriter<T>
01475 {
01476 public:
01478 VTKBinaryAppendedDataArrayWriter (std::ostream& theStream, std::string name, int ncomps, unsigned int& bc)
01479 : s(theStream),bytecount(bc)
01480 {
01481 VTKTypeNameTraits<T> tn;
01482 s << "<DataArray type=\"" << tn() << "\" Name=\"" << name << "\" ";
01483
01484 if (ncomps>3)
01485 DUNE_THROW(IOError, "VTKWriter does not support more than 3 components");
01486 s << "NumberOfComponents=\"" << (ncomps>1?3:1) << "\" ";
01487 s << "format=\"appended\" offset=\""<< bytecount << "\" />" << std::endl;
01488 bytecount += 4;
01489 }
01490
01492 void write (T data)
01493 {
01494 bytecount += sizeof(T);
01495 }
01496
01497 private:
01498 std::ostream& s;
01499 unsigned int& bytecount;
01500 };
01501
01502 protected:
01510 template<class T>
01511 VTKDataArrayWriter<T> *makeVTKDataArrayWriter(std::ostream &s,
01512 const char *name,
01513 unsigned int components,
01514 unsigned int totallength)
01515 {
01516 switch(outputtype) {
01517 case VTKOptions::ascii:
01518 return new VTKAsciiDataArrayWriter<T>(s, name, components);
01519 case VTKOptions::binary:
01520 return new VTKBinaryDataArrayWriter<T>(s, name, components, totallength);
01521 case VTKOptions::binaryappended:
01522 return new VTKBinaryAppendedDataArrayWriter<T>(s, name, components, bytecount);
01523 }
01524 DUNE_THROW(IOError, "VTKWriter: unsupported OutputType" << outputtype);
01525 }
01526
01528 class SimpleStream
01529 {
01530 public:
01532 SimpleStream (std::ostream& theStream)
01533 : s(theStream)
01534 {}
01536 template<class T>
01537 void write (T data)
01538 {
01539 char* p = reinterpret_cast<char*>(&data);
01540 s.write(p,sizeof(T));
01541 }
01542 private:
01543 std::ostream& s;
01544 };
01545
01547 void indentUp ()
01548 {
01549 indentCount++;
01550 }
01551
01553 void indentDown ()
01554 {
01555 indentCount--;
01556 }
01557
01559 void indent (std::ostream& s)
01560 {
01561 for (int i=0; i<indentCount; i++)
01562 s << " ";
01563 }
01564
01565
01566 static int renumber (const GeometryType &t, int i)
01567 {
01568 static const int quadRenumbering[4] = {0,1,3,2};
01569 static const int cubeRenumbering[8] = {0,1,3,2,4,5,7,6};
01570 static const int prismRenumbering[6] = {0,2,1,3,5,4};
01571 switch (vtkType(t))
01572 {
01573 case vtkQuadrilateral:
01574 return quadRenumbering[i];
01575 case vtkHexahedron:
01576 return cubeRenumbering[i];
01577 case vtkPrism:
01578 return prismRenumbering[i];
01579 default:
01580 return i;
01581 }
01582 }
01583 static int renumber (const Entity& e, int i)
01584 {
01585 return renumber(e.type(), i);
01586 }
01587
01588
01589 std::list<VTKFunction*> celldata;
01590 std::list<VTKFunction*> vertexdata;
01591
01592 private:
01593
01594 GridView gridView_;
01595 const Grid& grid;
01596
01597
01598 const IndexSet& is;
01599
01600
01601 int indentCount;
01602 int numPerLine;
01603
01604
01605 protected:
01606 int ncells;
01607 int nvertices;
01608 int ncorners;
01609 private:
01610 VertexMapper* vertexmapper;
01611 std::vector<int> number;
01612 VTKOptions::DataMode datamode;
01613 protected:
01614 VTKOptions::OutputType outputtype;
01615 private:
01616 unsigned int bytecount;
01617 };
01618
01622 template< class Grid >
01623 class LeafVTKWriter
01624 : public VTKWriter< typename Grid::LeafGridView >
01625 {
01626 typedef VTKWriter< typename Grid::LeafGridView > Base;
01627
01628 public:
01630 explicit LeafVTKWriter ( const Grid &grid,
01631 VTKOptions::DataMode dm = VTKOptions::conforming ) DUNE_DEPRECATED
01632 : Base( grid.leafView(), dm )
01633 {}
01634 };
01635
01639 template< class Grid >
01640 class LevelVTKWriter
01641 : public VTKWriter< typename Grid::LevelGridView >
01642 {
01643 typedef VTKWriter< typename Grid::LevelGridView > Base;
01644
01645 public:
01647 LevelVTKWriter ( const Grid &grid, int level,
01648 VTKOptions::DataMode dm = VTKOptions::conforming ) DUNE_DEPRECATED
01649 : Base( grid.levelView( level ), dm )
01650 {}
01651 };
01652
01653 }
01654
01655 #endif