34template<
typename T,
typename FUNCTOR, vtkType VTKTYPE>
36 : clout( std::cout,
"VTKwriter" ), _createFile(false), _name(name), _binary(binary), _compress(compress)
40template<
typename T,
typename FUNCTOR, vtkType VTKTYPE>
43 _pointerVec.push_back(&f);
46template<
typename T,
typename FUNCTOR, vtkType VTKTYPE>
49 f.getName() = functorName;
50 _pointerVec.push_back(&f);
54template<
typename T,
typename FUNCTOR, vtkType VTKTYPE>
55template<
bool vtmAsMaster>
60 std::string pathPrefixFileData;
61 if constexpr(!vtmAsMaster){
65 pathPrefixFileData =
"";
69 dataPVD(iT, fullNamePVD, nameVTM);
72 pathPrefixFileData =
"data/";
77 preambleVTM(fullNameVTM);
79 for (
int iC = 0; iC < rankSize; iC++) {
80 std::string nameFileData = pathPrefixFileData +
createFileName(nameCollection, iT, iC) + fileExtension;
82 dataVTM(iC, fullNameVTM, nameFileData);
86 closeVTM(fullNameVTM);
90template<
typename T,
typename FUNCTOR, vtkType VTKTYPE>
91template<
bool vtmAsMaster,
typename F>
94 std::string fileExtension;
98#ifdef PARALLEL_MODE_MPI
102 if ( _pointerVec.empty() ) {
103 clout <<
"Error: Did you add a Functor ?" << std::endl;
107 if constexpr (parallel) {
110 if constexpr(VTKTYPE==
VTI){
111 fileExtension =
".vti";
114 auto& superStructure =
dynamic_cast<SuperF<D,T>*
>(_pointerVec[0])->getSuperStructure();
115 auto const& cGeometry = superStructure.getCuboidDecomposition();
116 LoadBalancer<T> const& loadBalancer = superStructure.getLoadBalancer();
119 superStructure.communicate();
121 const T delta = cGeometry.getMotherCuboid().getDeltaR();
124 for (
int iCloc = 0; iCloc < loadBalancer.
size(); iCloc++) {
127 int globIC = loadBalancer.
glob(iCloc);
129 auto blockGeometry = cGeometry.get(globIC);
136 +
"data/" +
createFileName(nameCollection, iT, globIC ) + fileExtension;
140 auto originPhysR = cGeometry.getPhysR(originCuboid.
withPrefix(globIC));
142 preambleVTI(fullNameFileData, extent0, extent1, originPhysR.
data(), delta);
144 applyFunctors(fullNameFileData, extent1, globIC);
146 closePiece(fullNameFileData);
147 closeVTI(fullNameFileData);
151 }
else if constexpr(VTKTYPE==
VTU){
152 fileExtension =
".vtu";
155 LoadBalancer<T> const& loadBalancer = _pointerVec[0]->getLoadBalancer();
159 for (
int iCloc = 0; iCloc < loadBalancer.
size(); iCloc++) {
161 int globIC = loadBalancer.
glob(iCloc);
163 Vector<int,1> extent1 = _pointerVec[0]->getContainerF(iCloc).getContainerSize();
166 +
"data/" +
createFileName(nameCollection, iT, globIC ) + fileExtension;
168 preambleVTU(fullNameFileData,extent1);
170 applyFunctors(fullNameFileData, extent1, globIC);
172 this->cellDataVTU(fullNameFileData, extent1);
174 this->dataArrayPoints(fullNameFileData, extent1, globIC);
176 closeVTU(fullNameFileData);
181 fileExtension =
".vtp";
182 clout <<
"Error: VTP type not implemented yet" << std::endl;
186 if (rank == 0) { writeVTM<vtmAsMaster>(iT,rankSize,fileExtension,nameCollection); }
190 fileExtension =
".vtu";
200 std::string nameFileData =
"data/" +
createFileName(nameCollection, iT ) + fileExtension;
202 dataPVD(iT, fullNamePVD, nameFileData);
206 preambleVTU(fullNameFileData, extent1);
208 applyFunctors(fullNameFileData, extent1);
210 this->cellDataVTU(fullNameFileData, extent1);
212 this->dataArrayPoints(fullNameFileData, extent1);
214 closeVTU(fullNameFileData);
221template<
typename T,
typename FUNCTOR, vtkType VTKTYPE>
225 if constexpr(VTKTYPE==
VTI){
227 auto applyFunctors = [&](std::string fullNameFileData,
Vector<int,D> extent1,
int globIC=0 )
229 for (
int iF=0; iF<_pointerVec.size(); ++iF) {
230 this->dataArraySingleFunctor(fullNameFileData, *_pointerVec[iF], extent1, globIC);
234 write<false>(iT,_name,applyFunctors);
237 }
else if constexpr(VTKTYPE==
VTU){
239 auto applyFunctors = [&](std::string fullNameFileData,
Vector<int,1> extent1,
int globIC=0 )
241 for (std::size_t iF=1; iF<_pointerVec.size(); ++iF) {
242 this->dataArraySingleFunctor(fullNameFileData, *_pointerVec[iF], extent1, globIC);
246 write<false>(iT,_name, applyFunctors);
251template<
typename T,
typename FUNCTOR, vtkType VTKTYPE>
255 if constexpr(VTKTYPE==
VTI){
257 auto applyFunctors = [&](std::string fullNameFileData,
Vector<int,D> extent1,
int globIC=0 ){
258 this->dataArraySingleFunctor(fullNameFileData, f, extent1, globIC); };
260 write<true>(iT,f.getName(),applyFunctors);
263 }
else if constexpr(VTKTYPE==
VTU){
265 auto applyFunctors = [&](std::string fullNameFileData,
Vector<int,1> extent1,
int globIC=0 ){
266 this->dataArraySingleFunctor(fullNameFileData, f, extent1, globIC); };
268 write<true>(iT,f.getName(),applyFunctors);
273template<
typename T,
typename FUNCTOR, vtkType VTKTYPE>
280template<
typename T,
typename FUNCTOR, vtkType VTKTYPE>
284#ifdef PARALLEL_MODE_MPI
290 preamblePVD(fullNamePVD);
291 closePVD(fullNamePVD);
296template<
typename T,
typename FUNCTOR, vtkType VTKTYPE>
300 std::ofstream fout(fullName, std::ios::trunc);
302 clout <<
"Error: could not open " << fullName << std::endl;
304 fout <<
"<?xml version=\"1.0\"?>\n";
305 fout <<
"<VTKFile type=\"UnstructuredGrid\" version=\"0.1\" ";
307 fout <<
"byte_order=\"LittleEndian\" compressor=\"vtkZLibDataCompressor\">\n";
310 fout <<
"byte_order=\"LittleEndian\">\n";
312 fout <<
"<UnstructuredGrid>" << std::endl;
313 fout <<
"<Piece NumberOfPoints=\"" << extent1[0]
314 <<
"\" NumberOfCells=\"" << extent1[0] <<
"\">"
316 fout <<
"<PointData Vectors=\"Particles\">" << std::endl;
320template<
typename T,
typename FUNCTOR, vtkType VTKTYPE>
322 const std::string& fullNamePiece)
324 std::ofstream fout(fullNamePiece.c_str(), std::ios::app);
326 clout <<
"Error: could not open " << fullNamePiece << std::endl;
328 fout <<
"</Piece>" << std::endl;
329 fout <<
"</UnstructuredGrid>\n";
330 fout <<
"</VTKFile>\n";
338template<
typename T,
typename FUNCTOR, vtkType VTKTYPE>
340 const std::string& fullName,
344 std::ofstream fout(fullName, std::ios::trunc);
346 clout <<
"Error: could not open " << fullName << std::endl;
348 fout <<
"<?xml version=\"1.0\"?>\n";
349 fout <<
"<VTKFile type=\"ImageData\" version=\"0.1\" ";
351 fout <<
"byte_order=\"LittleEndian\" compressor=\"vtkZLibDataCompressor\">\n";
354 fout <<
"byte_order=\"LittleEndian\">\n";
357 fout <<
"<ImageData WholeExtent=\""
358 << extent0[0] <<
" "<< extent1[0];
359 for (
unsigned iDim=1; iDim<D; ++iDim){
360 fout <<
" " << extent0[iDim] <<
" " << extent1[iDim];
362 fout <<
"\" Origin=\"" << origin[0];
363 for (
unsigned iDim=1; iDim<D; ++iDim){
364 fout <<
" " << origin[iDim];
366 fout <<
"\" Spacing=\"" << delta <<
" " << delta <<
" " << delta <<
"\">\n";
367 fout <<
"<Piece Extent=\""
368 << extent0[0] <<
" "<< extent1[0];
369 for (
unsigned iDim=1; iDim<D; ++iDim){
370 fout <<
" " << extent0[iDim] <<
" "<< extent1[iDim];
373 fout <<
"<PointData>\n";
378template<
typename T,
typename FUNCTOR, vtkType VTKTYPE>
380 const std::string& fullNamePiece)
382 std::ofstream fout(fullNamePiece, std::ios::app );
384 clout <<
"Error: could not open " << fullNamePiece << std::endl;
386 fout <<
"</ImageData>\n";
387 fout <<
"</VTKFile>\n";
391template<
typename T,
typename FUNCTOR, vtkType VTKTYPE>
394 std::ofstream fout(fullNamePiece, std::ios::app );
396 clout <<
"Error: could not open " << fullNamePiece << std::endl;
398 fout <<
"</PointData>\n";
399 fout <<
"</Piece>\n";
404template<
typename T,
typename FUNCTOR, vtkType VTKTYPE>
406 const std::string& fullNamePVD)
408 std::ofstream fout(fullNamePVD.c_str(), std::ios::trunc);
410 clout <<
"Error: could not open " << fullNamePVD << std::endl;
413 fout <<
"<?xml version=\"1.0\"?>\n";
414 fout <<
"<VTKFile type=\"Collection\" version=\"0.1\" "
415 <<
"byte_order=\"LittleEndian\">\n" <<
"<Collection>\n";
419template<
typename T,
typename FUNCTOR, vtkType VTKTYPE>
421 const std::string& fullNamePVD)
423 std::ofstream fout(fullNamePVD.c_str(), std::ios::app);
425 clout <<
"Error: could not open " << fullNamePVD << std::endl;
427 fout <<
"</Collection>\n";
428 fout <<
"</VTKFile>\n";
435template<
typename T,
typename FUNCTOR, vtkType VTKTYPE>
438 std::ofstream fout(fullNameVTM, std::ios::trunc);
440 clout <<
"Error: could not open " << fullNameVTM << std::endl;
442 fout <<
"<?xml version=\"1.0\"?>\n";
443 fout <<
"<VTKFile type=\"vtkMultiBlockDataSet\" version=\"1.0\" "
444 <<
"byte_order=\"LittleEndian\">\n"
445 <<
"<vtkMultiBlockDataSet>\n" ;
449template<
typename T,
typename FUNCTOR, vtkType VTKTYPE>
452 std::ofstream fout(fullNameVTM, std::ios::app );
454 clout <<
"Error: could not open " << fullNameVTM << std::endl;
456 fout <<
"</vtkMultiBlockDataSet>\n";
457 fout <<
"</VTKFile>\n";
462template<
typename T,
typename FUNCTOR, vtkType VTKTYPE>
464 const std::string& namePiece)
466 std::ofstream fout(fullNameVTM, std::ios::app);
468 clout <<
"Error: could not open " << fullNameVTM << std::endl;
470 fout <<
"<Block index=\"" << iC <<
"\" >\n";
471 fout <<
"<DataSet index= \"0\" " <<
"file=\"" << namePiece <<
"\">\n"
473 fout <<
"</Block>\n";
477template<
typename T,
typename FUNCTOR, vtkType VTKTYPE>
479 const std::string& fullNamePVD,
const std::string& namePiece)
481 std::ofstream fout(fullNamePVD.c_str(),
482 std::ios::in | std::ios::out | std::ios::ate);
484 fout.seekp(-25, std::ios::end);
486 fout <<
"<DataSet timestep=\"" << iT <<
"\" "
487 <<
"group=\"\" part=\"\" "
488 <<
"file=\"" << namePiece <<
"\"/>\n";
490 closePVD(fullNamePVD);
493 clout <<
"Error: could not open " << fullNamePVD << std::endl;
497template<
typename T,
typename FUNCTOR, vtkType VTKTYPE>
498template<
unsigned sourceDim>
500 const std::string& fullName,
504 std::ofstream fout(fullName.c_str(), std::ios::app);
506 clout <<
"Error: could not open " << fullName << std::endl;
508 fout <<
"<Points>" << std::endl;
510 this->dataArraySingleFunctor(fullName, *_pointerVec[0], extent1, iC);
511 fout <<
"</Points>" << std::endl;
516template<
typename T,
typename FUNCTOR, vtkType VTKTYPE>
518 std::ofstream fout(fullName.c_str(), std::ios::app);
520 clout <<
"Error: could not open " << fullName << std::endl;
522 fout <<
"</PointData>" << std::endl;
523 fout <<
"<CellData /> " << std::endl;
524 fout <<
"<Cells>" << std::endl;
525 fout <<
"<DataArray type=\"Int32\" Name=\"connectivity\" format=\"ascii\">"
528 for (
int iTmp=0; iTmp<extent1[0]; ++iTmp){ fout << i32++ <<
" "; }
529 fout <<
"</DataArray>" << std::endl;
530 fout <<
"<DataArray type=\"Int32\" Name=\"offsets\" format=\"ascii\">"
533 for (
int iTmp=0; iTmp<extent1[0]; ++iTmp){ fout << i32++ <<
" "; }
534 fout <<
"</DataArray>" << std::endl;
535 fout <<
"<DataArray type=\"UInt8\" Name=\"types\" format=\"ascii\">"
537 for (
int iTmp=0; iTmp<extent1[0]; ++iTmp){ fout << 1 <<
" "; }
538 fout <<
"</DataArray>" << std::endl;
539 fout <<
"</Cells>" << std::endl;
544template<
typename T,
typename FUNCTOR, vtkType VTKTYPE>
545template<
unsigned sourceDim>
547 const std::string& fullName,
551 std::ofstream fout(fullName.c_str(), std::ios::app);
553 clout <<
"Error: could not open " << fullName << std::endl;
555 fout <<
"<DataArray type=\"Float32\" Name=\"" << f.getName() <<
"\" NumberOfComponents=\"" << f.getTargetDim() <<
"\" ";
556 if (_compress || _binary) {
557 fout <<
"format=\"binary\" encoding=\"base64\">\n";
565 std::vector<T> evaluated(f.getTargetDim());
568 if constexpr(parallel){ i[0] = iC;}
571 size_t numberOfFloats;
572 std::unique_ptr<float[]> streamFloat;
576 if constexpr(VTKTYPE==
VTI){
579 numberOfFloats = f.getTargetDim() * (extent1[0]+2) * (extent1[1]+2);
580 streamFloat = std::make_unique<float[]>(numberOfFloats);
581 for (i[2]=-1; i[2]<extent1[1]+1; ++i[2]) {
582 for (i[1]=-1; i[1]<extent1[0]+1; ++i[1]) {
583 f(evaluated.data(),i);
584 for (
int iDim=0; iDim < f.getTargetDim(); ++iDim){
585 streamFloat[itter] = float( evaluated[iDim] );
591 }
else if constexpr(D==3){
592 numberOfFloats = f.getTargetDim() * (extent1[0]+2) * (extent1[1]+2) * (extent1[2]+2);
593 streamFloat = std::make_unique<float[]>(numberOfFloats);
594 for (i[3]=-1; i[3]<extent1[2]+1; ++i[3]) {
595 for (i[2]=-1; i[2]<extent1[1]+1; ++i[2]) {
596 for (i[1]=-1; i[1]<extent1[0]+1; ++i[1]) {
597 f(evaluated.data(),i);
598 for (
int iDim=0; iDim < f.getTargetDim(); ++iDim){
599 streamFloat[itter] = float( evaluated[iDim] );
606 clout <<
"Error: only 2D and 3D supportet" << std::endl;
610 }
else if constexpr(VTKTYPE==
VTU){
611 numberOfFloats = f.getTargetDim() * extent1[0];
612 streamFloat = std::make_unique<float[]>(numberOfFloats);
614 if constexpr(parallel){
615 for (i[1]=0; i[1]<extent1[0]; ++i[1]) {
616 f(evaluated.data(),i);
617 for (
int iDim=0; iDim < f.getTargetDim(); ++iDim){
618 streamFloat[itter] = float( evaluated[iDim] );
624 for (i[0]=0; i[0]<extent1[0]; ++i[0]) {
625 f(evaluated.data(),i);
626 for (
int iDim=0; iDim < f.getTargetDim(); ++iDim){
627 streamFloat[itter] = float( evaluated[iDim] );
633 clout <<
"Error: VTP format not implemented yet" << std::endl;
637 uint32_t binarySize =
static_cast<uint32_t
>( numberOfFloats*
sizeof(float) );
642 const unsigned char* charData =
reinterpret_cast<unsigned char*
>(streamFloat.get());
644 std::unique_ptr<unsigned char[]> comprData(
new unsigned char[ binarySize ]);
647 uLongf sizeCompr = compressBound(binarySize);
648 compress2( comprData.get(), &sizeCompr, charData, binarySize, -1);
652 uint32_t prefix[4] = {1,binarySize,binarySize,
static_cast<uint32_t
>(sizeCompr)};
653 prefixEncoder.
encode(prefix, 4);
657 dataEncoder.
encode(comprData.get(), sizeCompr);
658 }
else if (_binary) {
661 prefixEncoder.
encode(&binarySize, 1);
664 dataEncoder.
encode(streamFloat.get(),numberOfFloats);
666 for (
size_t iOut = 0; iOut < numberOfFloats; ++iOut ) {
667 fout << streamFloat[iOut] <<
" ";
672 std::ofstream ffout( fullName, std::ios::out | std::ios::app );
674 clout <<
"Error: could not open " << fullName << std::endl;
676 ffout <<
"\n</DataArray>\n";
void encode(const T *data, size_t length)
Base class for all LoadBalancer.
void closeVTI(const std::string &fullNamePiece)
performes </ImageData> and </VTKFile>
void dataArraySingleFunctor(const std::string &fullName, FUNCTOR &f, Vector< int, sourceDim > extent1, int iC=0)
writes functors stored at pointerVec
void closeVTU(const std::string &fullNamePiece)
performes </ImageData> and </VTKFile>
void preambleVTM(const std::string &fullNameVTM)
performes <VTKFile ...> and <Collection>
VTKwriter(const std::string &name, bool binary=true, bool compress=true)
constructor
void preamblePVD(const std::string &fullNamePVD)
performes <VTKFile ...> and <Collection>
void createMasterFile()
have to be called before calling write(int iT=0), since it creates
void closeVTM(const std::string &fullNameVTM)
performes </Collection> and </VTKFile>
void dataVTM(int iC, const std::string &fullNameVTM, const std::string &namePiece)
performes <DataSet timestep= ... file=namePiece > used for linking vti into pvd files
void dataArrayPoints(const std::string &fullName, Vector< int, sourceDim > extent1, int iC=0)
writes points necessary for VTU
void write(int iT, std::string nameCollection, F applyFunctors)
to clear stored functors, not yet used due to lack of necessity
void closePiece(const std::string &fullNamePiece)
performes </PointData> and </Piece>
void preambleVTU(const std::string &fullName, Vector< int, 1 > extent1)
performes <VTKFile ...>, <ImageData ...> and <PieceExtent ...>
void preambleVTI(const std::string &fullName, const LatticeR< D > extent0, const LatticeR< D > extent1, PhysR< T, D > origin, T delta)
performes <VTKFile ...>, <ImageData ...>, <PieceExtent ...> and <PointData ...>
void closePVD(const std::string &fullNamePVD)
performes </Collection> and </VTKFile>
void cellDataVTU(const std::string &fullName, Vector< int, 1 > extent1)
TODO: add description: connectivity, offsete, type of unscructured nodes.
void addFunctor(FUNCTOR &f)
put functor to _pointerVec to simplify writing process of several functors
void writeVTM(int iT, int rankSize, std::string fileExtension, std::string nameCollection)
wrapper for VTM file creation
void dataPVD(int iT, const std::string &fullNamePVD, const std::string &namePiece)
performes <DataSet timestep= ... file=namePiece >
constexpr const T * data() const any_platform
constexpr Vector< T, D+1 > withPrefix(T prefix) const any_platform
Return new (prefix,this...) vector.
std::string getVtkOutDir() const
int getRank() const
Returns the process ID.
Directories & directories()
Top level namespace for all of OpenLB.
std::string createFileName(std::string name)
for .pvd masterFile
std::conditional_t< D==2, SuperF2D< T, U >, SuperF3D< T, U > > SuperF