Dune Core Modules (2.4.2)

vtksequencewriterbase.hh
1 // -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
2 // vi: set et ts=4 sw=2 sts=2:
3 
4 #ifndef DUNE_VTKSEQUENCEBASE_HH
5 #define DUNE_VTKSEQUENCEBASE_HH
6 
7 #include <vector>
8 #include <iostream>
9 #include <sstream>
10 #include <fstream>
11 #include <iomanip>
12 
14 #include <dune/common/path.hh>
16 
18 
19 namespace Dune {
20 
31  template<class GridView>
33  {
34  shared_ptr<VTKWriter<GridView> > vtkWriter_;
35  std::vector<double> timesteps_;
36  std::string name_,path_,extendpath_;
37  int rank_;
38  int size_;
39  public:
46  explicit VTKSequenceWriterBase( shared_ptr<VTKWriter<GridView> > vtkWriter,
47  const std::string& name,
48  const std::string& path,
49  const std::string& extendpath,
50  int rank,
51  int size)
52  : vtkWriter_(vtkWriter),
53  name_(name), path_(path),
54  extendpath_(extendpath),
55  rank_(rank),
56  size_(size)
57  {}
58 
60 
62  void addCellData (const shared_ptr<const typename VTKWriter<GridView>::VTKFunction> &p)
63  {
64  vtkWriter_->addCellData(p);
65  }
66 
69  {
70  vtkWriter_->addCellData(p);
71  }
72 
78  template<class V >
79  void addCellData (const V &v, const std::string &name, int ncomps=1)
80  {
81  vtkWriter_->addCellData(v, name, ncomps);
82  }
83 
86  {
87  vtkWriter_->addVertexData(p);
88  }
89 
91  void addVertexData (const typename VTKWriter<GridView>::VTKFunctionPtr &p)
92  {
93  vtkWriter_->addVertexData(p);
94  }
95 
101  template<class V >
102  void addVertexData (const V &v, const std::string &name, int ncomps=1)
103  {
104  vtkWriter_->addVertexData(v, name, ncomps);
105  }
106 
107 
113  void write (double time, VTK::OutputType ot = VTK::ascii)
114  {
115  /* remember current time step */
116  unsigned int count = timesteps_.size();
117  timesteps_.push_back(time);
118 
119  /* write VTK file */
120  if(size_==1)
121  vtkWriter_->write(concatPaths(path_,seqName(count)),ot);
122  else
123  vtkWriter_->pwrite(seqName(count), path_,extendpath_,ot);
124 
125  /* write pvd file ... only on rank 0 */
126  if (rank_==0) {
127  std::ofstream pvdFile;
128  pvdFile.exceptions(std::ios_base::badbit | std::ios_base::failbit |
129  std::ios_base::eofbit);
130  std::string pvdname = name_ + ".pvd";
131  pvdFile.open(pvdname.c_str());
132  pvdFile << "<?xml version=\"1.0\"?> \n"
133  << "<VTKFile type=\"Collection\" version=\"0.1\" byte_order=\"" << VTK::getEndiannessString() << "\"> \n"
134  << "<Collection> \n";
135  for (unsigned int i=0; i<=count; i++)
136  {
137  // filename
138  std::string piecepath;
139  std::string fullname;
140  if(size_==1) {
141  piecepath = path_;
142  fullname = vtkWriter_->getSerialPieceName(seqName(i), piecepath);
143  }
144  else {
145  piecepath = concatPaths(path_, extendpath_);
146  fullname = vtkWriter_->getParallelHeaderName(seqName(i), piecepath, size_);
147  }
148  pvdFile << "<DataSet timestep=\"" << timesteps_[i]
149  << "\" group=\"\" part=\"0\" name=\"\" file=\""
150  << fullname << "\"/> \n";
151  }
152  pvdFile << "</Collection> \n"
153  << "</VTKFile> \n" << std::flush;
154  pvdFile.close();
155  }
156  }
157  private:
158 
159  // create sequence name
160  std::string seqName(unsigned int count) const
161  {
162  std::stringstream n;
163  n.fill('0');
164  n << name_ << "-" << std::setw(5) << count;
165  return n.str();
166  }
167  };
168 
169 } // end namespace Dune
170 
171 #endif
A base class for grid functions with any return type and dimension.
Definition: function.hh:39
Base class to write pvd-files which contains a list of all collected vtk-files.
Definition: vtksequencewriterbase.hh:33
void addVertexData(const V &v, const std::string &name, int ncomps=1)
Adds a field of vertex data to the VTK file.
Definition: vtksequencewriterbase.hh:102
VTKSequenceWriterBase(shared_ptr< VTKWriter< GridView > > vtkWriter, const std::string &name, const std::string &path, const std::string &extendpath, int rank, int size)
Set up the VTKSequenceWriterBase class.
Definition: vtksequencewriterbase.hh:46
void addCellData(const shared_ptr< const typename VTKWriter< GridView >::VTKFunction > &p)
Adds a field of cell data to the VTK file.
Definition: vtksequencewriterbase.hh:62
void addVertexData(typename VTKWriter< GridView >::VTKFunction *p)
Adds a field of vertex data to the VTK file.
Definition: vtksequencewriterbase.hh:85
void addCellData(const V &v, const std::string &name, int ncomps=1)
Adds a field of cell data to the VTK file.
Definition: vtksequencewriterbase.hh:79
void addVertexData(const typename VTKWriter< GridView >::VTKFunctionPtr &p)
Adds a field of vertex data to the VTK file.
Definition: vtksequencewriterbase.hh:91
void addCellData(typename VTKWriter< GridView >::VTKFunction *p)
Adds a field of cell data to the VTK file.
Definition: vtksequencewriterbase.hh:68
void write(double time, VTK::OutputType ot=VTK::ascii)
Writes VTK data for the given time,.
Definition: vtksequencewriterbase.hh:113
Writer for the ouput of grid functions in the vtk format.
Definition: vtkwriter.hh:87
Common stuff for the VTKWriter.
OutputType
How the bulk data should be stored in the file.
Definition: common.hh:40
std::string concatPaths(const std::string &base, const std::string &p)
concatenate two paths
Definition: path.cc:30
Dune namespace.
Definition: alignment.hh:10
Utilities for handling filesystem paths.
This file implements the class shared_ptr (a reference counting pointer), for those systems that don'...
Provides file i/o for the visualization toolkit.
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.80.0 (May 16, 22:29, 2024)