OpenLB 1.7
Loading...
Searching...
No Matches
Public Member Functions | Static Public Attributes | Protected Member Functions | List of all members
olb::VTKwriter< T, FUNCTOR, VTKTYPE > Class Template Reference

#include <vtkWriter.h>

+ Collaboration diagram for olb::VTKwriter< T, FUNCTOR, VTKTYPE >:

Public Member Functions

 VTKwriter (const std::string &name, bool binary=true, bool compress=true)
 constructor
 
void addFunctor (FUNCTOR &f)
 put functor to _pointerVec to simplify writing process of several functors
 
void addFunctor (FUNCTOR &f, const std::string &functorName)
 put functor with specific name to _pointerVec to simplify writing process of several functors
 
template<bool vtmAsMaster, typename F >
void write (int iT, std::string nameCollection, F applyFunctors)
 to clear stored functors, not yet used due to lack of necessity
 
void write (int iT=0)
 
void write (FUNCTOR &f, int iT=0)
 
void write (std::shared_ptr< FUNCTOR > ptr_f, int iT=0)
 
void createMasterFile ()
 have to be called before calling write(int iT=0), since it creates
 

Static Public Attributes

static constexpr bool parallel = FUNCTOR::isSuper
 
static constexpr unsigned D = FUNCTOR::d
 

Protected Member Functions

void preamblePVD (const std::string &fullNamePVD)
 performes <VTKFile ...> and <Collection>
 
void dataPVD (int iT, const std::string &fullNamePVD, const std::string &namePiece)
 performes <DataSet timestep= ... file=namePiece >
 
void closePVD (const std::string &fullNamePVD)
 performes </Collection> and </VTKFile>
 
void preambleVTM (const std::string &fullNameVTM)
 performes <VTKFile ...> and <Collection>
 
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
 
template<bool vtmAsMaster = false>
void writeVTM (int iT, int rankSize, std::string fileExtension, std::string nameCollection)
 wrapper for VTM file creation
 
void preambleVTU (const std::string &fullName, Vector< int, 1 > extent1)
 performes <VTKFile ...>, <ImageData ...> and <PieceExtent ...>
 
void closeVTU (const std::string &fullNamePiece)
 performes </ImageData> and </VTKFile>
 
void cellDataVTU (const std::string &fullName, Vector< int, 1 > extent1)
 TODO: add description: connectivity, offsete, type of unscructured nodes.
 
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 closeVTI (const std::string &fullNamePiece)
 performes </ImageData> and </VTKFile>
 
void closePiece (const std::string &fullNamePiece)
 performes </PointData> and </Piece>
 
template<unsigned sourceDim>
void dataArrayPoints (const std::string &fullName, Vector< int, sourceDim > extent1, int iC=0)
 writes points necessary for VTU
 
template<unsigned sourceDim>
void dataArraySingleFunctor (const std::string &fullName, FUNCTOR &f, Vector< int, sourceDim > extent1, int iC=0)
 writes functors stored at pointerVec
 

Detailed Description

template<typename T, typename FUNCTOR, vtkType VTKTYPE>
class olb::VTKwriter< T, FUNCTOR, VTKTYPE >

Definition at line 56 of file vtkWriter.h.

Constructor & Destructor Documentation

◆ VTKwriter()

template<typename T , typename FUNCTOR , vtkType VTKTYPE>
olb::VTKwriter< T, FUNCTOR, VTKTYPE >::VTKwriter ( const std::string & name,
bool binary = true,
bool compress = true )

constructor

Definition at line 35 of file vtkWriter.hh.

36 : clout( std::cout, "VTKwriter" ), _createFile(false), _name(name), _binary(binary), _compress(compress)
37{}

Member Function Documentation

◆ addFunctor() [1/2]

template<typename T , typename FUNCTOR , vtkType VTKTYPE>
void olb::VTKwriter< T, FUNCTOR, VTKTYPE >::addFunctor ( FUNCTOR & f)

put functor to _pointerVec to simplify writing process of several functors

Definition at line 41 of file vtkWriter.hh.

42{
43 _pointerVec.push_back(&f);
44}

◆ addFunctor() [2/2]

template<typename T , typename FUNCTOR , vtkType VTKTYPE>
void olb::VTKwriter< T, FUNCTOR, VTKTYPE >::addFunctor ( FUNCTOR & f,
const std::string & functorName )

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

Definition at line 47 of file vtkWriter.hh.

48{
49 f.getName() = functorName;
50 _pointerVec.push_back(&f);
51}

◆ cellDataVTU()

template<typename T , typename FUNCTOR , vtkType VTKTYPE>
void olb::VTKwriter< T, FUNCTOR, VTKTYPE >::cellDataVTU ( const std::string & fullName,
Vector< int, 1 > extent1 )
protected

TODO: add description: connectivity, offsete, type of unscructured nodes.

Definition at line 519 of file vtkWriter.hh.

519 {
520 std::ofstream fout(fullName.c_str(), std::ios::app);
521 if (!fout) {
522 clout << "Error: could not open " << fullName << std::endl;
523 }
524 fout << "</PointData>" << std::endl;
525 fout << "<CellData /> " << std::endl; //TODO: open tag missing?
526 fout << "<Cells>" << std::endl;
527 fout << "<DataArray type=\"Int32\" Name=\"connectivity\" format=\"ascii\">"
528 << std::endl;
529 int32_t i32 = 0;
530 for (int iTmp=0; iTmp<extent1[0]; ++iTmp){ fout << i32++ << " "; }
531 fout << "</DataArray>" << std::endl;
532 fout << "<DataArray type=\"Int32\" Name=\"offsets\" format=\"ascii\">"
533 << std::endl;
534 i32 = 1;
535 for (int iTmp=0; iTmp<extent1[0]; ++iTmp){ fout << i32++ << " "; }
536 fout << "</DataArray>" << std::endl;
537 fout << "<DataArray type=\"UInt8\" Name=\"types\" format=\"ascii\">"
538 << std::endl;
539 for (int iTmp=0; iTmp<extent1[0]; ++iTmp){ fout << 1 << " "; }
540 fout << "</DataArray>" << std::endl;
541 fout << "</Cells>" << std::endl;
542}

◆ closePiece()

template<typename T , typename FUNCTOR , vtkType VTKTYPE>
void olb::VTKwriter< T, FUNCTOR, VTKTYPE >::closePiece ( const std::string & fullNamePiece)
protected

performes </PointData> and </Piece>

Definition at line 394 of file vtkWriter.hh.

395{
396 std::ofstream fout(fullNamePiece, std::ios::app );
397 if (!fout) {
398 clout << "Error: could not open " << fullNamePiece << std::endl;
399 }
400 fout << "</PointData>\n";
401 fout << "</Piece>\n";
402 fout.close();
403}

◆ closePVD()

template<typename T , typename FUNCTOR , vtkType VTKTYPE>
void olb::VTKwriter< T, FUNCTOR, VTKTYPE >::closePVD ( const std::string & fullNamePVD)
protected

performes </Collection> and </VTKFile>

Definition at line 422 of file vtkWriter.hh.

424{
425 std::ofstream fout(fullNamePVD.c_str(), std::ios::app);
426 if (!fout) {
427 clout << "Error: could not open " << fullNamePVD << std::endl;
428 }
429 fout << "</Collection>\n";
430 fout << "</VTKFile>\n";
431 fout.close();
432}

◆ closeVTI()

template<typename T , typename FUNCTOR , vtkType VTKTYPE>
void olb::VTKwriter< T, FUNCTOR, VTKTYPE >::closeVTI ( const std::string & fullNamePiece)
protected

performes </ImageData> and </VTKFile>

Definition at line 381 of file vtkWriter.hh.

383{
384 std::ofstream fout(fullNamePiece, std::ios::app );
385 if (!fout) {
386 clout << "Error: could not open " << fullNamePiece << std::endl;
387 }
388 fout << "</ImageData>\n";
389 fout << "</VTKFile>\n";
390 fout.close();
391}

◆ closeVTM()

template<typename T , typename FUNCTOR , vtkType VTKTYPE>
void olb::VTKwriter< T, FUNCTOR, VTKTYPE >::closeVTM ( const std::string & fullNameVTM)
protected

performes </Collection> and </VTKFile>

Definition at line 452 of file vtkWriter.hh.

453{
454 std::ofstream fout(fullNameVTM, std::ios::app );
455 if (!fout) {
456 clout << "Error: could not open " << fullNameVTM << std::endl;
457 }
458 fout << "</vtkMultiBlockDataSet>\n";
459 fout << "</VTKFile>\n";
460 fout.close();
461}

◆ closeVTU()

template<typename T , typename FUNCTOR , vtkType VTKTYPE>
void olb::VTKwriter< T, FUNCTOR, VTKTYPE >::closeVTU ( const std::string & fullNamePiece)
protected

performes </ImageData> and </VTKFile>

Definition at line 323 of file vtkWriter.hh.

325{
326 std::ofstream fout(fullNamePiece.c_str(), std::ios::app);
327 if (!fout) {
328 clout << "Error: could not open " << fullNamePiece << std::endl;
329 }
330 fout << "</Piece>" << std::endl;
331 fout << "</UnstructuredGrid>\n";
332 fout << "</VTKFile>\n";
333 fout.close();
334}

◆ createMasterFile()

template<typename T , typename FUNCTOR , vtkType VTKTYPE>
void olb::VTKwriter< T, FUNCTOR, VTKTYPE >::createMasterFile ( )

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

Definition at line 283 of file vtkWriter.hh.

284{
285 int rank = 0;
286#ifdef PARALLEL_MODE_MPI
287 rank = singleton::mpi().getRank();
288#endif
289 if (rank == 0) {
290 std::string fullNamePVD = singleton::directories().getVtkOutDir()
291 + createFileName(_name) + ".pvd";
292 preamblePVD(fullNamePVD);
293 closePVD(fullNamePVD);
294 _createFile = true;
295 }
296}
void preamblePVD(const std::string &fullNamePVD)
performes <VTKFile ...> and <Collection>
Definition vtkWriter.hh:407
void closePVD(const std::string &fullNamePVD)
performes </Collection> and </VTKFile>
Definition vtkWriter.hh:422
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:

◆ dataArrayPoints()

template<typename T , typename FUNCTOR , vtkType VTKTYPE>
template<unsigned sourceDim>
void olb::VTKwriter< T, FUNCTOR, VTKTYPE >::dataArrayPoints ( const std::string & fullName,
Vector< int, sourceDim > extent1,
int iC = 0 )
protected

writes points necessary for VTU

Definition at line 501 of file vtkWriter.hh.

505{
506 std::ofstream fout(fullName.c_str(), std::ios::app);
507 if (!fout) {
508 clout << "Error: could not open " << fullName << std::endl;
509 }
510 fout << "<Points>" << std::endl;
511 // Call dataArray with first functor (which should contain points)
512 this->dataArraySingleFunctor(fullName, *_pointerVec[0], extent1, iC);
513 fout << "</Points>" << std::endl;
514}
void dataArraySingleFunctor(const std::string &fullName, FUNCTOR &f, Vector< int, sourceDim > extent1, int iC=0)
writes functors stored at pointerVec
Definition vtkWriter.hh:548

◆ dataArraySingleFunctor()

template<typename T , typename FUNCTOR , vtkType VTKTYPE>
template<unsigned sourceDim>
void olb::VTKwriter< T, FUNCTOR, VTKTYPE >::dataArraySingleFunctor ( const std::string & fullName,
FUNCTOR & f,
Vector< int, sourceDim > extent1,
int iC = 0 )
protected

writes functors stored at pointerVec

Definition at line 548 of file vtkWriter.hh.

552{
553 std::ofstream fout(fullName.c_str(), std::ios::app);
554 if (!fout) {
555 clout << "Error: could not open " << fullName << std::endl;
556 }
557 fout << "<DataArray type=\"Float32\" Name=\"" << f.getName() << "\" NumberOfComponents=\"" << f.getTargetDim() << "\" ";
558 if (_compress || _binary) {
559 fout << "format=\"binary\" encoding=\"base64\">\n";
560 }
561 else {
562 fout << ">\n";
563 }
564
565 // Create functor input and output
566 int i[4] = {0}; //4 as maximum dimension
567 std::vector<T> evaluated(f.getTargetDim());
568
569 //Set ic if parallel
570 if constexpr(parallel){ i[0] = iC;}
571
572 //Create float array
573 size_t numberOfFloats;
574 std::unique_ptr<float[]> streamFloat;
575 int itter = 0;
576
577 //VTI
578 if constexpr(VTKTYPE==VTI){
579 //2D
580 if constexpr(D==2){
581 numberOfFloats = f.getTargetDim() * (extent1[0]+2) * (extent1[1]+2);
582 streamFloat = std::make_unique<float[]>(numberOfFloats);
583 for (i[2]=-1; i[2]<extent1[1]+1; ++i[2]) {
584 for (i[1]=-1; i[1]<extent1[0]+1; ++i[1]) {
585 f(evaluated.data(),i);
586 for (int iDim=0; iDim < f.getTargetDim(); ++iDim){
587 streamFloat[itter] = float( evaluated[iDim] );
588 ++itter;
589 }
590 }
591 }
592 //3D
593 } else if constexpr(D==3){
594 numberOfFloats = f.getTargetDim() * (extent1[0]+2) * (extent1[1]+2) * (extent1[2]+2);
595 streamFloat = std::make_unique<float[]>(numberOfFloats);
596 for (i[3]=-1; i[3]<extent1[2]+1; ++i[3]) {
597 for (i[2]=-1; i[2]<extent1[1]+1; ++i[2]) {
598 for (i[1]=-1; i[1]<extent1[0]+1; ++i[1]) {
599 f(evaluated.data(),i);
600 for (int iDim=0; iDim < f.getTargetDim(); ++iDim){
601 streamFloat[itter] = float( evaluated[iDim] );
602 ++itter;
603 }
604 }
605 }
606 }
607 } else {
608 clout << "Error: only 2D and 3D supportet" << std::endl;
609 }
610
611 //VTU
612 } else if constexpr(VTKTYPE==VTU){
613 numberOfFloats = f.getTargetDim() * extent1[0];
614 streamFloat = std::make_unique<float[]>(numberOfFloats);
615 //PARALLEL
616 if constexpr(parallel){
617 for (i[1]=0; i[1]<extent1[0]; ++i[1]) {
618 f(evaluated.data(),i);
619 for (int iDim=0; iDim < f.getTargetDim(); ++iDim){
620 streamFloat[itter] = float( evaluated[iDim] );
621 ++itter;
622 }
623 }
624 //NONPARALLEL
625 } else {
626 for (i[0]=0; i[0]<extent1[0]; ++i[0]) {
627 f(evaluated.data(),i);
628 for (int iDim=0; iDim < f.getTargetDim(); ++iDim){
629 streamFloat[itter] = float( evaluated[iDim] );
630 ++itter;
631 }
632 }
633 }
634 } else {
635 clout << "Error: VTP format not implemented yet" << std::endl;
636 }
637
638 //Define binarSize from number of floats
639 uint32_t binarySize = static_cast<uint32_t>( numberOfFloats*sizeof(float) );
640
641 //Write data (TODO: repair for non parallised vtu version)
642 if (_compress) {
643 // char buffer for functor data
644 const unsigned char* charData = reinterpret_cast<unsigned char*>(streamFloat.get());
645 // buffer for compression
646 std::unique_ptr<unsigned char[]> comprData(new unsigned char[ binarySize ]); // stack may be too small
647
648 // compress data (not yet decoded as base64) by zlib
649 uLongf sizeCompr = compressBound(binarySize);
650 compress2( comprData.get(), &sizeCompr, charData, binarySize, -1);
651
652 // encode prefix to base64 documented in http://www.earthmodels.org/software/vtk-and-paraview/vtk-file-formats
653 Base64Encoder<uint32_t> prefixEncoder(fout, 4);
654 uint32_t prefix[4] = {1,binarySize,binarySize,static_cast<uint32_t>(sizeCompr)};
655 prefixEncoder.encode(prefix, 4);
656
657 // encode compressed data to base64
658 Base64Encoder<unsigned char> dataEncoder( fout, sizeCompr );
659 dataEncoder.encode(comprData.get(), sizeCompr);
660 } else if (_binary) {
661 // encode prefix to base64 documented in http://www.earthmodels.org/software/vtk-and-paraview/vtk-file-formats
662 Base64Encoder<uint32_t> prefixEncoder(fout, 1);
663 prefixEncoder.encode(&binarySize, 1);
664 // write numbers from functor
665 Base64Encoder<float> dataEncoder(fout, numberOfFloats);
666 dataEncoder.encode(streamFloat.get(),numberOfFloats);
667 } else {
668 for ( size_t iOut = 0; iOut < numberOfFloats; ++iOut ) {
669 fout << streamFloat[iOut] << " ";
670 }
671 }
672 fout.close();
673
674 std::ofstream ffout( fullName, std::ios::out | std::ios::app );
675 if (!ffout) {
676 clout << "Error: could not open " << fullName << std::endl;
677 }
678 ffout << "\n</DataArray>\n";
679 ffout.close();
680}
static constexpr unsigned D
Definition vtkWriter.h:59
static constexpr bool parallel
Definition vtkWriter.h:58
@ VTU
Definition vtkWriter.h:52
@ VTI
Definition vtkWriter.h:52

References olb::Base64Encoder< T >::encode(), olb::VTI, and olb::VTU.

+ Here is the call graph for this function:

◆ dataPVD()

template<typename T , typename FUNCTOR , vtkType VTKTYPE>
void olb::VTKwriter< T, FUNCTOR, VTKTYPE >::dataPVD ( int iT,
const std::string & fullNamePVD,
const std::string & namePiece )
protected

performes <DataSet timestep= ... file=namePiece >

Definition at line 480 of file vtkWriter.hh.

482{
483 std::ofstream fout(fullNamePVD.c_str(),
484 std::ios::in | std::ios::out | std::ios::ate);
485 if (fout) {
486 fout.seekp(-25, std::ios::end); // jump -25 form the end of file to overwrite closePVD
487
488 fout << "<DataSet timestep=\"" << iT << "\" "
489 << "group=\"\" part=\"\" "
490 << "file=\"" << namePiece << "\"/>\n";
491 fout.close();
492 closePVD(fullNamePVD);
493 }
494 else {
495 clout << "Error: could not open " << fullNamePVD << std::endl;
496 }
497}

◆ dataVTM()

template<typename T , typename FUNCTOR , vtkType VTKTYPE>
void olb::VTKwriter< T, FUNCTOR, VTKTYPE >::dataVTM ( int iC,
const std::string & fullNameVTM,
const std::string & namePiece )
protected

performes <DataSet timestep= ... file=namePiece > used for linking vti into pvd files

Definition at line 465 of file vtkWriter.hh.

467{
468 std::ofstream fout(fullNameVTM, std::ios::app);
469 if (!fout) {
470 clout << "Error: could not open " << fullNameVTM << std::endl;
471 }
472 fout << "<Block index=\"" << iC << "\" >\n";
473 fout << "<DataSet index= \"0\" " << "file=\"" << namePiece << "\">\n"
474 << "</DataSet>\n";
475 fout << "</Block>\n";
476 fout.close();
477}

◆ preamblePVD()

template<typename T , typename FUNCTOR , vtkType VTKTYPE>
void olb::VTKwriter< T, FUNCTOR, VTKTYPE >::preamblePVD ( const std::string & fullNamePVD)
protected

performes <VTKFile ...> and <Collection>

Definition at line 407 of file vtkWriter.hh.

409{
410 std::ofstream fout(fullNamePVD.c_str(), std::ios::trunc);
411 if (!fout) {
412 clout << "Error: could not open " << fullNamePVD << std::endl;
413 }
414
415 fout << "<?xml version=\"1.0\"?>\n";
416 fout << "<VTKFile type=\"Collection\" version=\"0.1\" "
417 << "byte_order=\"LittleEndian\">\n" << "<Collection>\n";
418 fout.close();
419}

◆ preambleVTI()

template<typename T , typename FUNCTOR , vtkType VTKTYPE>
void olb::VTKwriter< T, FUNCTOR, VTKTYPE >::preambleVTI ( const std::string & fullName,
const LatticeR< D > extent0,
const LatticeR< D > extent1,
PhysR< T, D > origin,
T delta )
protected

performes <VTKFile ...>, <ImageData ...>, <PieceExtent ...> and <PointData ...>

Definition at line 341 of file vtkWriter.hh.

345{
346 std::ofstream fout(fullName, std::ios::trunc);
347 if (!fout) {
348 clout << "Error: could not open " << fullName << std::endl;
349 }
350 fout << "<?xml version=\"1.0\"?>\n";
351 fout << "<VTKFile type=\"ImageData\" version=\"0.1\" ";
352 if (_compress) {
353 fout << "byte_order=\"LittleEndian\" compressor=\"vtkZLibDataCompressor\">\n";
354 }
355 else {
356 fout << "byte_order=\"LittleEndian\">\n";
357 }
358
359 fout << "<ImageData WholeExtent=\""
360 << extent0[0] <<" "<< extent1[0];
361 for (unsigned iDim=1; iDim<D; ++iDim){
362 fout << " " << extent0[iDim] << " " << extent1[iDim];
363 }
364 fout << "\" Origin=\"" << origin[0];
365 for (unsigned iDim=1; iDim<D; ++iDim){
366 fout << " " << origin[iDim];
367 }
368 fout << "\" Spacing=\"" << delta << " " << delta << " " << delta << "\">\n";
369 fout << "<Piece Extent=\""
370 << extent0[0] <<" "<< extent1[0];
371 for (unsigned iDim=1; iDim<D; ++iDim){
372 fout << " " << extent0[iDim] <<" "<< extent1[iDim];
373 }
374 fout <<"\">\n";
375 fout << "<PointData>\n";
376 fout.close();
377}

◆ preambleVTM()

template<typename T , typename FUNCTOR , vtkType VTKTYPE>
void olb::VTKwriter< T, FUNCTOR, VTKTYPE >::preambleVTM ( const std::string & fullNameVTM)
protected

performes <VTKFile ...> and <Collection>

Definition at line 438 of file vtkWriter.hh.

439{
440 std::ofstream fout(fullNameVTM, std::ios::trunc);
441 if (!fout) {
442 clout << "Error: could not open " << fullNameVTM << std::endl;
443 }
444 fout << "<?xml version=\"1.0\"?>\n";
445 fout << "<VTKFile type=\"vtkMultiBlockDataSet\" version=\"1.0\" "
446 << "byte_order=\"LittleEndian\">\n"
447 << "<vtkMultiBlockDataSet>\n" ;
448 fout.close();
449}

◆ preambleVTU()

template<typename T , typename FUNCTOR , vtkType VTKTYPE>
void olb::VTKwriter< T, FUNCTOR, VTKTYPE >::preambleVTU ( const std::string & fullName,
Vector< int, 1 > extent1 )
protected

performes <VTKFile ...>, <ImageData ...> and <PieceExtent ...>

Definition at line 299 of file vtkWriter.hh.

301{
302 std::ofstream fout(fullName, std::ios::trunc);
303 if (!fout) {
304 clout << "Error: could not open " << fullName << std::endl;
305 }
306 fout << "<?xml version=\"1.0\"?>\n";
307 fout << "<VTKFile type=\"UnstructuredGrid\" version=\"0.1\" ";
308 if (_compress) {
309 fout << "byte_order=\"LittleEndian\" compressor=\"vtkZLibDataCompressor\">\n";
310 }
311 else {
312 fout << "byte_order=\"LittleEndian\">\n";
313 }
314 fout << "<UnstructuredGrid>" << std::endl;
315 fout << "<Piece NumberOfPoints=\"" << extent1[0]
316 << "\" NumberOfCells=\"" << extent1[0] << "\">"
317 << std::endl;
318 fout << "<PointData Vectors=\"Particles\">" << std::endl;
319 fout.close();
320}

◆ write() [1/4]

template<typename T , typename FUNCTOR , vtkType VTKTYPE>
void olb::VTKwriter< T, FUNCTOR, VTKTYPE >::write ( FUNCTOR & f,
int iT = 0 )

Definition at line 254 of file vtkWriter.hh.

255{
256 //VTI
257 if constexpr(VTKTYPE==VTI){
258 //Define functor
259 auto applyFunctors = [&](std::string fullNameFileData, Vector<int,D> extent1, int globIC=0 ){
260 this->dataArraySingleFunctor(fullNameFileData, f, extent1, globIC); };
261 //Call write
262 write<true>(iT,f.getName(),applyFunctors);
263
264 //VTU
265 } else if constexpr(VTKTYPE==VTU){
266 //Define functor
267 auto applyFunctors = [&](std::string fullNameFileData, Vector<int,1> extent1, int globIC=0 ){
268 this->dataArraySingleFunctor(fullNameFileData, f, extent1, globIC); };
269 //Call write
270 write<true>(iT,f.getName(),applyFunctors);
271 }
272
273}

References olb::VTI, and olb::VTU.

◆ write() [2/4]

template<typename T , typename FUNCTOR , vtkType VTKTYPE>
template<bool vtmAsMaster, typename F >
void olb::VTKwriter< T, FUNCTOR, VTKTYPE >::write ( int iT,
std::string nameCollection,
F applyFunctors )

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

Definition at line 92 of file vtkWriter.hh.

93{
94
95 std::string fileExtension;
96
97 int rank = 0;
98 int rankSize;
99#ifdef PARALLEL_MODE_MPI
100 rank = singleton::mpi().getRank();
101#endif
102
103 if ( _pointerVec.empty() ) {
104 clout << "Error: Did you add a Functor ?" << std::endl;
105 } else {
106
107 //PARALLEL FUNCTORTYPE
108 if constexpr (parallel) {
109
110 //VTI
111 if constexpr(VTKTYPE==VTI){
112 fileExtension = ".vti";
113
114 //Retrieve cuboidGeometry from first functor (TODO: include into common functortype)
115 auto& superStructure = dynamic_cast<SuperF<D,T>*>(_pointerVec[0])->getSuperStructure();
116 auto const& cGeometry = superStructure.getCuboidGeometry();
117 LoadBalancer<T> const& loadBalancer = superStructure.getLoadBalancer();
118 rankSize = loadBalancer.getRankSize();
119 //Ensure no gaps between vti files (TODO: still necessary?)
120 superStructure.communicate();
121 //Retrieve delta
122 const T delta = cGeometry.getMotherCuboid().getDeltaR();
123
124 //Cor each cuboid
125 for (int iCloc = 0; iCloc < loadBalancer.size(); iCloc++) {
126
127 //Retrieve global IC
128 int globIC = loadBalancer.glob(iCloc);
129 //Retrieve BlockGeometry
130 auto blockGeometry = cGeometry.get(globIC);
131 //Retrieve extent
132 Vector<int,D> extent0(-1);
133 Vector<int,D> extent1( blockGeometry.getExtent() );
134
135 //Name of Data file
136 std::string fullNameFileData = singleton::directories().getVtkOutDir()
137 + "data/" + createFileName(nameCollection, iT, globIC ) + fileExtension;
138
139 // get dimension/extent for each cuboid
140 LatticeR<D> originCuboid(0.);
141 T originPhysR[D] = {T()};
142 cGeometry.getPhysR(originPhysR, globIC, originCuboid);
143 //Write preamble VTU
144 preambleVTI(fullNameFileData, extent0, extent1, originPhysR, delta);
145 //Loop over functors or call individual one
146 applyFunctors(fullNameFileData, extent1, globIC);
147 //Close piece and vti
148 closePiece(fullNameFileData);
149 closeVTI(fullNameFileData);
150 }
151
152 //VTU
153 } else if constexpr(VTKTYPE==VTU){
154 fileExtension = ".vtu";
155
156 //Retrive loadBalancer from first functor (TODO: include into common functortype)
157 LoadBalancer<T> const& loadBalancer = _pointerVec[0]->getLoadBalancer();
158 rankSize = loadBalancer.getRankSize();
159
160 //Cor each cuboid
161 for (int iCloc = 0; iCloc < loadBalancer.size(); iCloc++) {
162 //Retrieve global IC
163 int globIC = loadBalancer.glob(iCloc);
164 // Retrieve container size (number of particles) (TODO: include into common functortype)
165 Vector<int,1> extent1 = _pointerVec[0]->getContainerF(iCloc).getContainerSize();
166 //Name of VTU
167 std::string fullNameFileData = singleton::directories().getVtkOutDir()
168 + "data/" + createFileName(nameCollection, iT, globIC ) + fileExtension;
169 //Write preamble VTU
170 preambleVTU(fullNameFileData,extent1);
171 //Loop over functors or call individual one
172 applyFunctors(fullNameFileData, extent1, globIC);
173 // Write celldata (connectivity, offset, types)
174 this->cellDataVTU(fullNameFileData, extent1);
175 // Write points
176 this->dataArrayPoints(fullNameFileData, extent1, globIC);
177 // CloseVTU
178 closeVTU(fullNameFileData);
179 }
180
181 //VTP
182 } else {
183 fileExtension = ".vtp";
184 clout << "Error: VTP type not implemented yet" << std::endl;
185 }
186
187 //Write VTM (only master)
188 if (rank == 0) { writeVTM<vtmAsMaster>(iT,rankSize,fileExtension,nameCollection); }
189
190 //UNPARRALLIZED FUNCTORTYPE (for now only relevant for VTU)
191 } else {
192 fileExtension = ".vtu";
193
194 //if master
195 if (rank == 0) {
196 // Full name pvd
197 std::string fullNamePVD = singleton::directories().getVtkOutDir()
198 + createFileName(nameCollection) + ".pvd";
199 // Retrieve container size (number of particles)(TODO: include into common functortype)
200 Vector<int,1> extent1(_pointerVec[0]->getContainerSize());
201 // Name VTU
202 std::string nameFileData = "data/" + createFileName(nameCollection, iT ) + fileExtension;
203 // Create PVD
204 dataPVD(iT, fullNamePVD, nameFileData);
205 // Full Name VTU
206 std::string fullNameFileData = singleton::directories().getVtkOutDir() + nameFileData;
207 // Write preample VTU
208 preambleVTU(fullNameFileData, extent1);
209 //Loop over functors or call individual one
210 applyFunctors(fullNameFileData, extent1);
211 // Write celldata (connectivity, offset, types)
212 this->cellDataVTU(fullNameFileData, extent1);
213 // Write points
214 this->dataArrayPoints(fullNameFileData, extent1);
215 // CloseVTU
216 closeVTU(fullNameFileData);
217 } // master only
218 } // if parralized type
219 } //if ( _pointerVec.empty() )
220}
void closeVTI(const std::string &fullNamePiece)
performes </ImageData> and </VTKFile>
Definition vtkWriter.hh:381
void closeVTU(const std::string &fullNamePiece)
performes </ImageData> and </VTKFile>
Definition vtkWriter.hh:323
void dataArrayPoints(const std::string &fullName, Vector< int, sourceDim > extent1, int iC=0)
writes points necessary for VTU
Definition vtkWriter.hh:501
void closePiece(const std::string &fullNamePiece)
performes </PointData> and </Piece>
Definition vtkWriter.hh:394
void preambleVTU(const std::string &fullName, Vector< int, 1 > extent1)
performes <VTKFile ...>, <ImageData ...> and <PieceExtent ...>
Definition vtkWriter.hh:299
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 ...>
Definition vtkWriter.hh:341
void cellDataVTU(const std::string &fullName, Vector< int, 1 > extent1)
TODO: add description: connectivity, offsete, type of unscructured nodes.
Definition vtkWriter.hh:519
void dataPVD(int iT, const std::string &fullNamePVD, const std::string &namePiece)
performes <DataSet timestep= ... file=namePiece >
Definition vtkWriter.hh:480

References olb::createFileName(), olb::singleton::directories(), olb::singleton::MpiManager::getRank(), olb::LoadBalancer< T >::getRankSize(), olb::singleton::Directories::getVtkOutDir(), olb::LoadBalancer< T >::glob(), olb::singleton::mpi(), olb::LoadBalancer< T >::size(), olb::VTI, and olb::VTU.

+ Here is the call graph for this function:

◆ write() [3/4]

template<typename T , typename FUNCTOR , vtkType VTKTYPE>
void olb::VTKwriter< T, FUNCTOR, VTKTYPE >::write ( int iT = 0)

Definition at line 224 of file vtkWriter.hh.

225{
226 //VTI
227 if constexpr(VTKTYPE==VTI){
228 //Define functors to be applied
229 auto applyFunctors = [&](std::string fullNameFileData, Vector<int,D> extent1, int globIC=0 )
230 {
231 for (int iF=0; iF<_pointerVec.size(); ++iF) {
232 this->dataArraySingleFunctor(fullNameFileData, *_pointerVec[iF], extent1, globIC);
233 }
234 };
235 //Call write
236 write<false>(iT,_name,applyFunctors);
237
238 //VTU
239 } else if constexpr(VTKTYPE==VTU){
240 //Define functors to be applied
241 auto applyFunctors = [&](std::string fullNameFileData, Vector<int,1> extent1, int globIC=0 )
242 {
243 for (std::size_t iF=1; iF<_pointerVec.size(); ++iF) {
244 this->dataArraySingleFunctor(fullNameFileData, *_pointerVec[iF], extent1, globIC);
245 }
246 };
247 //Call write
248 write<false>(iT,_name, applyFunctors);
249 }
250}

References olb::VTI, and olb::VTU.

◆ write() [4/4]

template<typename T , typename FUNCTOR , vtkType VTKTYPE>
void olb::VTKwriter< T, FUNCTOR, VTKTYPE >::write ( std::shared_ptr< FUNCTOR > ptr_f,
int iT = 0 )

Definition at line 276 of file vtkWriter.hh.

277{
278 write(*ptr_f, iT);
279}
void write(int iT, std::string nameCollection, F applyFunctors)
to clear stored functors, not yet used due to lack of necessity
Definition vtkWriter.hh:92

◆ writeVTM()

template<typename T , typename FUNCTOR , vtkType VTKTYPE>
template<bool vtmAsMaster>
void olb::VTKwriter< T, FUNCTOR, VTKTYPE >::writeVTM ( int iT,
int rankSize,
std::string fileExtension,
std::string nameCollection )
protected

wrapper for VTM file creation

Definition at line 56 of file vtkWriter.hh.

57{
58 //Evaluate vtmAsMaster
59 std::string nameVTM;
60 std::string pathPrefixFileData;
61 if constexpr(!vtmAsMaster){
62 //Name of vtm
63 nameVTM = "data/" + createFileName(nameCollection, iT ) + ".vtm";
64 //path prefix
65 pathPrefixFileData = "";
66 //Create PVD
67 std::string fullNamePVD = singleton::directories().getVtkOutDir()
68 + createFileName(nameCollection) + ".pvd";
69 dataPVD(iT, fullNamePVD, nameVTM);
70 } else {
71 nameVTM = createFileName(nameCollection, iT) + ".vtm";
72 pathPrefixFileData = "data/";
73 }
74 //Create full name VTM
75 std::string fullNameVTM = singleton::directories().getVtkOutDir() + nameVTM;
76 //Write preamble VTM
77 preambleVTM(fullNameVTM);
78 //Loop over all cuboids (globIcs)
79 for (int iC = 0; iC < rankSize; iC++) {
80 std::string nameFileData = pathPrefixFileData + createFileName(nameCollection, iT, iC) + fileExtension;
81 // puts name of .vti piece to a .pvd file [fullNameVTM]
82 dataVTM(iC, fullNameVTM, nameFileData);
83 // adds a namePiece to master.pvd file.
84 }
85 // To do so we overwrite closePVD() and add new entry.
86 closeVTM(fullNameVTM); // timestep
87}
void preambleVTM(const std::string &fullNameVTM)
performes <VTKFile ...> and <Collection>
Definition vtkWriter.hh:438
void closeVTM(const std::string &fullNameVTM)
performes </Collection> and </VTKFile>
Definition vtkWriter.hh:452
void dataVTM(int iC, const std::string &fullNameVTM, const std::string &namePiece)
performes <DataSet timestep= ... file=namePiece > used for linking vti into pvd files
Definition vtkWriter.hh:465

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

+ Here is the call graph for this function:

Member Data Documentation

◆ D

template<typename T , typename FUNCTOR , vtkType VTKTYPE>
constexpr unsigned olb::VTKwriter< T, FUNCTOR, VTKTYPE >::D = FUNCTOR::d
staticconstexpr

Definition at line 59 of file vtkWriter.h.

◆ parallel

template<typename T , typename FUNCTOR , vtkType VTKTYPE>
constexpr bool olb::VTKwriter< T, FUNCTOR, VTKTYPE >::parallel = FUNCTOR::isSuper
staticconstexpr

Definition at line 58 of file vtkWriter.h.


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