OpenLB 1.7
Loading...
Searching...
No Matches
Public Member Functions | List of all members
olb::SuperVTMwriter3D< T, OUT_T, W > Class Template Reference

SuperVTMwriter3D writes any SuperF3D to vtk-based output files. More...

#include <superVtmWriter3D.h>

+ Collaboration diagram for olb::SuperVTMwriter3D< T, OUT_T, W >:

Public Member Functions

 SuperVTMwriter3D (const std::string &name, int overlap=1, bool binary=true, bool compress=true)
 Construct writer for functor output.
 
 SuperVTMwriter3D (CuboidGeometry3D< T > &cGeometry, const std::string &name, int overlap=1, bool binary=true, bool compress=true)
 Construct writer for CuboidGeometry3D debugging.
 
void write (int iT=0)
 writes functors stored in pointerVec every process writes a vti file with data for each of its cuboids the vti files are linked in a pvd file
 
void writePVD (int iT)
 writes only the linking pvd file for timestep iT, blocks must be written separately (e.g. asynchronously)
 
void writeGlobalVTI (int iT, int iC)
 writes the vti file for cuboid iC at timestep iT
 
void writeVTI (int iT, int iCloc)
 writes the vti file for cuboid iCloc at timestep iT
 
void write (SuperF3D< T, W > &f, int iT=0)
 writes functor instantaneously, same vti-pvd file structure as above
 
void write (std::shared_ptr< SuperF3D< T, W > > ptr_f, int iT=0)
 
void createMasterFile ()
 have to be called before calling write(int iT=0), since it creates
 
void addFunctor (SuperF3D< T, W > &f)
 put functor to _pointerVec to simplify writing process of several functors
 
void addFunctor (SuperF3D< T, W > &f, const std::string &functorName)
 put functor with specific name to _pointerVec to simplify writing process of several functors
 
void clearAddedFunctors ()
 to clear stored functors, not yet used due to lack of necessity
 
std::string getName () const
 getter for _name
 

Detailed Description

template<typename T, typename OUT_T = float, typename W = T>
class olb::SuperVTMwriter3D< T, OUT_T, W >

SuperVTMwriter3D writes any SuperF3D to vtk-based output files.

In .pvd files, there are only links/references to a VTKmultiblock file 'vtm'

.pvd file structure the time series is represented by different 'vtm' files.

.vtm file This file links cuboids ('vti') and represents the entire data of a single timestep.

Definition at line 51 of file superVtmWriter3D.h.

Constructor & Destructor Documentation

◆ SuperVTMwriter3D() [1/2]

template<typename T , typename OUT_T , typename W >
olb::SuperVTMwriter3D< T, OUT_T, W >::SuperVTMwriter3D ( const std::string & name,
int overlap = 1,
bool binary = true,
bool compress = true )

Construct writer for functor output.

Definition at line 54 of file superVtmWriter3D.hh.

55 : clout( std::cout,"SuperVTMwriter3D" ), _createFile(false), _name(name), _overlap(overlap), _binary(binary), _compress(compress)
56{
57 static_assert(std::is_same_v<OUT_T, float> || std::is_same_v<OUT_T, double>,
58 "OUT_T must be either float or double");
59}

◆ SuperVTMwriter3D() [2/2]

template<typename T , typename OUT_T , typename W >
olb::SuperVTMwriter3D< T, OUT_T, W >::SuperVTMwriter3D ( CuboidGeometry3D< T > & cGeometry,
const std::string & name,
int overlap = 1,
bool binary = true,
bool compress = true )

Construct writer for CuboidGeometry3D debugging.

Definition at line 62 of file superVtmWriter3D.hh.

64 : SuperVTMwriter3D(name, overlap, binary, compress)
65{
66 _cGeometry = &cGeometry;
67}
SuperVTMwriter3D(const std::string &name, int overlap=1, bool binary=true, bool compress=true)
Construct writer for functor output.

Member Function Documentation

◆ addFunctor() [1/2]

template<typename T , typename OUT_T , typename W >
void olb::SuperVTMwriter3D< T, OUT_T, W >::addFunctor ( SuperF3D< T, W > & f)

put functor to _pointerVec to simplify writing process of several functors

Definition at line 257 of file superVtmWriter3D.hh.

258{
259 _pointerVec.push_back(&f);
260}

◆ addFunctor() [2/2]

template<typename T , typename OUT_T , typename W >
void olb::SuperVTMwriter3D< T, OUT_T, W >::addFunctor ( SuperF3D< T, W > & f,
const std::string & functorName )

put functor with specific name to _pointerVec to simplify writing process of several functors

Definition at line 263 of file superVtmWriter3D.hh.

264{
265 f.getName() = functorName;
266 _pointerVec.push_back(&f);
267}
std::string & getName()
read and write access to name
Definition genericF.hh:51

References olb::GenericF< T, S >::getName().

+ Here is the call graph for this function:

◆ clearAddedFunctors()

template<typename T , typename OUT_T , typename W >
void olb::SuperVTMwriter3D< T, OUT_T, W >::clearAddedFunctors ( )

to clear stored functors, not yet used due to lack of necessity

Definition at line 270 of file superVtmWriter3D.hh.

271{
272 _pointerVec.clear();
273}

◆ createMasterFile()

template<typename T , typename OUT_T , typename W >
void olb::SuperVTMwriter3D< T, OUT_T, W >::createMasterFile ( )

have to be called before calling write(int iT=0), since it creates

Definition at line 241 of file superVtmWriter3D.hh.

242{
243 int rank = 0;
244#ifdef PARALLEL_MODE_MPI
245 rank = singleton::mpi().getRank();
246#endif
247 if ( rank == 0 ) {
248 const std::string fullNamePVDmaster = singleton::directories().getVtkOutDir()
249 + createFileName( _name ) + ".pvd";
250 preamblePVD(fullNamePVDmaster);
251 closePVD(fullNamePVDmaster);
252 _createFile = true;
253 }
254}
std::string getVtkOutDir() const
Definition singleton.h:97
int getRank() const
Returns the process ID.
MpiManager & mpi()
Directories & directories()
Definition singleton.h:150
std::string createFileName(std::string name)
for .pvd masterFile
Definition fileName.hh:34

References olb::createFileName(), olb::singleton::directories(), olb::singleton::MpiManager::getRank(), olb::singleton::Directories::getVtkOutDir(), and olb::singleton::mpi().

+ Here is the call graph for this function:

◆ getName()

template<typename T , typename OUT_T , typename W >
std::string olb::SuperVTMwriter3D< T, OUT_T, W >::getName ( ) const

getter for _name

Definition at line 276 of file superVtmWriter3D.hh.

277{
278 return _name;
279}

◆ write() [1/3]

template<typename T , typename OUT_T , typename W >
void olb::SuperVTMwriter3D< T, OUT_T, W >::write ( int iT = 0)

writes functors stored in pointerVec every process writes a vti file with data for each of its cuboids the vti files are linked in a pvd file

Definition at line 70 of file superVtmWriter3D.hh.

71{
72 // update to prevent gaps between vti files / cuboids
73 for (SuperF3D<T,W>* f : _pointerVec) {
74 f->getSuperStructure().communicate();
75 }
76
77 // to get first element _pointerVec
78 // problem if functors with different SuperStructure are stored
79 // since till now, there is only one origin
80 const auto it_begin = _pointerVec.cbegin();
81 if (!_cGeometry && it_begin == _pointerVec.end()) {
82 throw std::runtime_error("No functor to write");
83 }
84 CuboidGeometry3D<T> const& cGeometry = _cGeometry ? *_cGeometry : (**it_begin).getSuperStructure().getCuboidGeometry();
85
86 // PVD, owns all
87 writePVD(iT);
88 if (_cGeometry) {
89 // Write globally if cuboid geometry is provided
90 for (int iC = 0; iC < cGeometry.getNc(); ++iC) {
91 writeGlobalVTI(iT, iC);
92 }
93 } else {
94 LoadBalancer<T>& load = (**it_begin).getSuperStructure().getLoadBalancer();
95 // VTI, each process writes its cuboids
96 for (int iCloc = 0; iCloc < load.size(); iCloc++) {
97 writeVTI(iT, iCloc);
98 }
99 }
100}
void writeGlobalVTI(int iT, int iC)
writes the vti file for cuboid iC at timestep iT
void writeVTI(int iT, int iCloc)
writes the vti file for cuboid iCloc at timestep iT
void writePVD(int iT)
writes only the linking pvd file for timestep iT, blocks must be written separately (e....

References olb::CuboidGeometry3D< T >::getNc(), and olb::LoadBalancer< T >::size().

+ Here is the call graph for this function:

◆ write() [2/3]

template<typename T , typename OUT_T , typename W >
void olb::SuperVTMwriter3D< T, OUT_T, W >::write ( std::shared_ptr< SuperF3D< T, W > > ptr_f,
int iT = 0 )

Definition at line 235 of file superVtmWriter3D.hh.

236{
237 write(*ptr_f, iT);
238}
void write(int iT=0)
writes functors stored in pointerVec every process writes a vti file with data for each of its cuboid...

◆ write() [3/3]

template<typename T , typename OUT_T , typename W >
void olb::SuperVTMwriter3D< T, OUT_T, W >::write ( SuperF3D< T, W > & f,
int iT = 0 )

writes functor instantaneously, same vti-pvd file structure as above

Definition at line 183 of file superVtmWriter3D.hh.

184{
185 CuboidGeometry3D<T> const& cGeometry = f.getSuperStructure().getCuboidGeometry();
186 LoadBalancer<T>& load = f.getSuperStructure().getLoadBalancer();
187 // no gaps between vti files (cuboids)
189 const T delta = cGeometry.getMotherCuboid().getDeltaR();
190
191 int rank = 0;
192#ifdef PARALLEL_MODE_MPI
193 rank = singleton::mpi().getRank();
194#endif
195
196 // write a pvd file, which links all vti files
197 // each vti file is written by one thread, which may own severals cuboids
198 if ( rank == 0 ) {
199 // master only
200 const std::string pathVTM = singleton::directories().getVtkOutDir()
201 + createFileName( f.getName(), iT ) + ".vtm";
202
203 preambleVTM(pathVTM);
204 for (int iC = 0; iC < cGeometry.getNc(); iC++) {
205 const std::string nameVTI = "data/" + createFileName( f.getName(), iT, iC) + ".vti";
206 // puts name of .vti piece to a .pvd file [fullNamePVD]
207 dataVTM( iC, pathVTM, nameVTI );
208 }
209 closeVTM(pathVTM);
210 } // master only
211
212 for (int iCloc = 0; iCloc < load.size(); iCloc++) {
213 // get piece/whole extent
214 const Vector<int,3> extent0(-_overlap,-_overlap,-_overlap);
215 const Vector<int,3> extent1( cGeometry.get(load.glob(iCloc)).getExtent());
216
217 const std::string fullNameVTI = singleton::directories().getVtkOutDir() + "data/"
218 + createFileName( f.getName(), iT, load.glob(iCloc) ) + ".vti";
219
220 // get dimension/extent for each cuboid
221 const int originLatticeR[4] = {load.glob(iCloc),0,0,0};
222 T originPhysR[3] = {T()};
223 cGeometry.getPhysR(originPhysR,originLatticeR);
224
225 preambleVTI(fullNameVTI, extent0, (extent1+_overlap-1.), originPhysR, delta);
226
227 dataArray(fullNameVTI, f, load.glob(iCloc), extent1);
228 closePiece(fullNameVTI);
229 closeVTI(fullNameVTI);
230 } // cuboid
231}
SuperStructure< T, 3 > & getSuperStructure()
CuboidGeometry< T, D > & getCuboidGeometry()
Read and write access to cuboid geometry.
virtual void communicate()
LoadBalancer< T > & getLoadBalancer()
Read and write access to the load balancer.

References olb::SuperStructure< T, D >::communicate(), olb::createFileName(), olb::singleton::directories(), olb::CuboidGeometry3D< T >::get(), olb::SuperStructure< T, D >::getCuboidGeometry(), olb::SuperStructure< T, D >::getLoadBalancer(), olb::CuboidGeometry3D< T >::getMotherCuboid(), olb::GenericF< T, S >::getName(), olb::CuboidGeometry3D< T >::getNc(), olb::CuboidGeometry3D< T >::getPhysR(), olb::singleton::MpiManager::getRank(), olb::SuperF3D< T, W >::getSuperStructure(), olb::singleton::Directories::getVtkOutDir(), olb::LoadBalancer< T >::glob(), olb::singleton::mpi(), and olb::LoadBalancer< T >::size().

+ Here is the call graph for this function:

◆ writeGlobalVTI()

template<typename T , typename OUT_T , typename W >
void olb::SuperVTMwriter3D< T, OUT_T, W >::writeGlobalVTI ( int iT,
int iC )

writes the vti file for cuboid iC at timestep iT

Definition at line 128 of file superVtmWriter3D.hh.

129{
130 const auto it_begin = _pointerVec.cbegin();
131 CuboidGeometry3D<T> const& cGeometry = _cGeometry ? *_cGeometry : (**it_begin).getSuperStructure().getCuboidGeometry();
132 const T delta = cGeometry.getMotherCuboid().getDeltaR();
133
134 // get piece/whole extent
135 const Vector<int,3> extent0(-_overlap,-_overlap,-_overlap);
136 const Vector<int,3> extent1(cGeometry.get(iC).getExtent());
137
138 const std::string fullNameVTI = singleton::directories().getVtkOutDir() + "data/"
139 + createFileName(_name, iT, iC) + ".vti";
140
141 // get dimension/extent for each cuboid
142 const int originLatticeR[4] = {iC,0,0,0};
143 T originPhysR[3] = {T()};
144 cGeometry.getPhysR(originPhysR,originLatticeR);
145
146 preambleVTI(fullNameVTI, extent0, (extent1+_overlap-1), originPhysR, delta);
147 for (auto it : _pointerVec) {
148 dataArray(fullNameVTI, *it, iC, extent1);
149 }
150 closePiece(fullNameVTI);
151 closeVTI(fullNameVTI);
152}

References olb::createFileName(), olb::singleton::directories(), olb::CuboidGeometry3D< T >::get(), olb::CuboidGeometry3D< T >::getMotherCuboid(), olb::CuboidGeometry3D< T >::getPhysR(), and olb::singleton::Directories::getVtkOutDir().

+ Here is the call graph for this function:

◆ writePVD()

template<typename T , typename OUT_T , typename W >
void olb::SuperVTMwriter3D< T, OUT_T, W >::writePVD ( int iT)

writes only the linking pvd file for timestep iT, blocks must be written separately (e.g. asynchronously)

Definition at line 103 of file superVtmWriter3D.hh.

104{
105 // to get first element _pointerVec
106 // problem if functors with different SuperStructure are stored
107 // since till now, there is only one origin
108 const auto it_begin = _pointerVec.cbegin();
109 CuboidGeometry3D<T> const& cGeometry = _cGeometry ? *_cGeometry : (**it_begin).getSuperStructure().getCuboidGeometry();
110
111 // PVD, owns all
112 if (singleton::mpi().isMainProcessor()) {
113 const std::string pathPVD = singleton::directories().getVtkOutDir()
114 + createFileName(_name) + ".pvd";
115 dataPVDmaster(iT, pathPVD, "data/" + createFileName(_name, iT) + ".vtm");
116
117 const std::string pathVTM = singleton::directories().getVtkOutDir()
118 + "data/" + createFileName(_name, iT) + ".vtm";
119 preambleVTM(pathVTM);
120 for (int iC = 0; iC < cGeometry.getNc(); iC++) {
121 dataVTM(iC, pathVTM, createFileName(_name, iT, iC) + ".vti" );
122 }
123 closeVTM(pathVTM);
124 }
125}

References olb::createFileName(), olb::singleton::directories(), olb::CuboidGeometry3D< T >::getNc(), olb::singleton::Directories::getVtkOutDir(), and olb::singleton::mpi().

+ Here is the call graph for this function:

◆ writeVTI()

template<typename T , typename OUT_T , typename W >
void olb::SuperVTMwriter3D< T, OUT_T, W >::writeVTI ( int iT,
int iCloc )

writes the vti file for cuboid iCloc at timestep iT

Definition at line 155 of file superVtmWriter3D.hh.

156{
157 const auto it_begin = _pointerVec.cbegin();
158 CuboidGeometry3D<T> const& cGeometry = (**it_begin).getSuperStructure().getCuboidGeometry();
159 LoadBalancer<T>& load = (**it_begin).getSuperStructure().getLoadBalancer();
160 const T delta = cGeometry.getMotherCuboid().getDeltaR();
161
162 // get piece/whole extent
163 const Vector<int,3> extent0(-_overlap,-_overlap,-_overlap);
164 const Vector<int,3> extent1(cGeometry.get(load.glob(iCloc)).getExtent());
165
166 const std::string fullNameVTI = singleton::directories().getVtkOutDir() + "data/"
167 + createFileName(_name, iT, load.glob(iCloc)) + ".vti";
168
169 // get dimension/extent for each cuboid
170 const int originLatticeR[4] = {load.glob(iCloc),0,0,0};
171 T originPhysR[3] = {T()};
172 cGeometry.getPhysR(originPhysR,originLatticeR);
173
174 preambleVTI(fullNameVTI, extent0, (extent1+_overlap-1), originPhysR, delta);
175 for (auto it : _pointerVec) {
176 dataArray(fullNameVTI, *it, load.glob(iCloc), extent1);
177 }
178 closePiece(fullNameVTI);
179 closeVTI(fullNameVTI);
180}

References olb::createFileName(), olb::singleton::directories(), olb::CuboidGeometry3D< T >::get(), olb::CuboidGeometry3D< T >::getMotherCuboid(), olb::CuboidGeometry3D< T >::getPhysR(), olb::singleton::Directories::getVtkOutDir(), and olb::LoadBalancer< T >::glob().

+ Here is the call graph for this function:

The documentation for this class was generated from the following files: