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