Dune Core Modules (2.6.0)

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_GRID_IO_FILE_VTK_VTKSEQUENCEWRITERBASE_HH
5#define DUNE_GRID_IO_FILE_VTK_VTKSEQUENCEWRITERBASE_HH
6
7#include <vector>
8#include <iostream>
9#include <sstream>
10#include <fstream>
11#include <iomanip>
12#include <memory>
13
15#include <dune/common/path.hh>
16
18
19namespace Dune {
20
30 template<class GridView>
32 {
33 std::shared_ptr<VTKWriter<GridView> > vtkWriter_;
34 std::vector<double> timesteps_;
35 std::string name_,path_,extendpath_;
36 int rank_;
37 int size_;
38 public:
46 const std::string& name,
47 const std::string& path,
48 const std::string& extendpath,
49 int rank,
50 int size)
51 : vtkWriter_(vtkWriter),
52 name_(name), path_(path),
53 extendpath_(extendpath),
54 rank_(rank),
55 size_(size)
56 {}
57
61 const std::shared_ptr< VTKWriter<GridView> >& vtkWriter() const
62 {
63 return vtkWriter_;
64 }
65
67 void addCellData (const std::shared_ptr<const typename VTKWriter<GridView>::VTKFunction> &p)
68 {
69 vtkWriter_->addCellData(p);
70 }
71
77 template<class V >
78 void addCellData (const V &v, const std::string &name, int ncomps=1)
79 {
80 vtkWriter_->addCellData(v, name, ncomps);
81 }
82
84 void addVertexData (const std::shared_ptr<const typename VTKWriter<GridView>::VTKFunction> &p)
85 {
86 vtkWriter_->addVertexData(p);
87 }
88
94 template<class V >
95 void addVertexData (const V &v, const std::string &name, int ncomps=1)
96 {
97 vtkWriter_->addVertexData(v, name, ncomps);
98 }
99
100
106 void write (double time, VTK::OutputType type = VTK::ascii)
107 {
108 /* remember current time step */
109 unsigned int count = timesteps_.size();
110 timesteps_.push_back(time);
111
112 /* write VTK file */
113 if(size_==1)
114 vtkWriter_->write(concatPaths(path_,seqName(count)),type);
115 else
116 vtkWriter_->pwrite(seqName(count), path_,extendpath_,type);
117
118 /* write pvd file ... only on rank 0 */
119 if (rank_==0) {
120 std::ofstream pvdFile;
121 pvdFile.exceptions(std::ios_base::badbit | std::ios_base::failbit |
122 std::ios_base::eofbit);
123 std::string pvdname = name_ + ".pvd";
124 pvdFile.open(pvdname.c_str());
125 pvdFile << "<?xml version=\"1.0\"?> \n"
126 << "<VTKFile type=\"Collection\" version=\"0.1\" byte_order=\"" << VTK::getEndiannessString() << "\"> \n"
127 << "<Collection> \n";
128 for (unsigned int i=0; i<=count; i++)
129 {
130 // filename
131 std::string piecepath;
132 std::string fullname;
133 if(size_==1) {
134 piecepath = path_;
135 fullname = vtkWriter_->getSerialPieceName(seqName(i), piecepath);
136 }
137 else {
138 piecepath = concatPaths(path_, extendpath_);
139 fullname = vtkWriter_->getParallelHeaderName(seqName(i), piecepath, size_);
140 }
141 pvdFile << "<DataSet timestep=\"" << timesteps_[i]
142 << "\" group=\"\" part=\"0\" name=\"\" file=\""
143 << fullname << "\"/> \n";
144 }
145 pvdFile << "</Collection> \n"
146 << "</VTKFile> \n" << std::flush;
147 pvdFile.close();
148 }
149 }
150 private:
151
152 // create sequence name
153 std::string seqName(unsigned int count) const
154 {
155 std::stringstream n;
156 n.fill('0');
157 n << name_ << "-" << std::setw(5) << count;
158 return n.str();
159 }
160 };
161
162} // end namespace Dune
163
164#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:32
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:95
const std::shared_ptr< VTKWriter< GridView > > & vtkWriter() const
Definition: vtksequencewriterbase.hh:61
void write(double time, VTK::OutputType type=VTK::ascii)
Writes VTK data for the given time,.
Definition: vtksequencewriterbase.hh:106
void addVertexData(const std::shared_ptr< const typename VTKWriter< GridView >::VTKFunction > &p)
Adds a field of vertex data to the VTK file.
Definition: vtksequencewriterbase.hh:84
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:78
VTKSequenceWriterBase(std::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:45
void addCellData(const std::shared_ptr< const typename VTKWriter< GridView >::VTKFunction > &p)
Adds a field of cell data to the VTK file.
Definition: vtksequencewriterbase.hh:67
Writer for the ouput of grid functions in the vtk format.
Definition: vtkwriter.hh:88
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: alignedallocator.hh:10
Utilities for handling filesystem paths.
Provides file i/o for the visualization toolkit.
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Dec 26, 23:30, 2024)