Microflow 3D  v1.0
MF::RW::VTIWriter< T_VDBGridType > Class Template Reference

The VTIWriter class provides a write interface for data storage in .vti format. More...

#include <VTIWriter.h>

Public Member Functions

 VTIWriter (std::shared_ptr< T_VDBGridType > VDBGrid_Ptr, const std::string &FileName, float voxelSize)
 
 ~VTIWriter ()=default
 
template<typename T3 >
void addDataFromGrid (const std::string &DataName, unsigned int numComponents, std::shared_ptr< T_VDBGridType > VDBGrid_Ptr)
 
template<typename T3 >
void addDataFromGrid (const std::string &DataName, unsigned int numComponents)
 
void addNodeIDFromGrid ()
 
void addNodeCoordinateFromGrid ()
 
template<typename T3_VTKArrayType , typename T4 >
void addDataFromTable (const std::string &DataName, unsigned int numComponents, std::vector< T4 > &table_Ptr)
 
void addNodeCoordinateFromThread (const std::shared_ptr< MF::Database::ThreadArray > &ThreadArray_Ptr)
 
void addNodeTypeFromThread (const std::shared_ptr< MF::Database::ThreadArray > &ThreadArray_Ptr)
 
void addUidThreadFromThread (const std::shared_ptr< MF::Database::ThreadArray > &ThreadArray_Ptr)
 
void addNodeThreadsByNameThread (const std::shared_ptr< MF::Database::ThreadArray > &ThreadArray_Ptr)
 
void addVelocityLBFromThread (const std::shared_ptr< MF::Database::ThreadArray > &ThreadArray_Ptr)
 
void addInitialVelocityLBFromThread (const std::shared_ptr< MF::Database::ThreadArray > &ThreadArray_Ptr)
 
void addVelocityPhysFromThread (const std::shared_ptr< MF::Database::ThreadArray > &ThreadArray_Ptr)
 
void addRhoLBFromThread (const std::shared_ptr< MF::Database::ThreadArray > &ThreadArray_Ptr)
 
void addInitialRhoLBFromThread (const std::shared_ptr< MF::Database::ThreadArray > &ThreadArray_Ptr)
 
void addRhoPhysFromThread (const std::shared_ptr< MF::Database::ThreadArray > &ThreadArray_Ptr)
 
void addPressPhysFromThread (const std::shared_ptr< MF::Database::ThreadArray > &ThreadArray_Ptr)
 
void addFQ19FromThread (const std::shared_ptr< MF::Database::ThreadArray > &ThreadArray_Ptr)
 
void addPropagationDirFromThread (const std::shared_ptr< MF::Database::ThreadArray > &ThreadArray_Ptr)
 
void writeFile ()
 

Static Public Member Functions

static std::shared_ptr< MF::RW::VTIWriter< T_VDBGridType > > New (std::shared_ptr< T_VDBGridType > VDBGrid_Ptr, const std::string &FileName, float voxelSize)
 

Detailed Description

template<class T_VDBGridType>
class MF::RW::VTIWriter< T_VDBGridType >

The VTIWriter class provides a write interface for data storage in .vti format.

Parameters
template<class t_vdbgridtype>=""> = <openvdb::FloatGrid::Ptr>, <openvdb::Int32Grid::Ptr>, <openvdb::DoubleGrid::Ptr> etc.

Definition at line 44 of file VTIWriter.h.

Constructor & Destructor Documentation

◆ VTIWriter()

template<class T_VDBGridType >
MF::RW::VTIWriter< T_VDBGridType >::VTIWriter ( std::shared_ptr< T_VDBGridType >  VDBGrid_Ptr,
const std::string &  FileName,
float  voxelSize 
)
inline

Definition at line 46 of file VTIWriter.h.

References MF::RW::VTIWriter< T_VDBGridType >::~VTIWriter().

46  {
47  m_FileNameExt = FileName + ".vti";
48  // The bounding box for active voxels
49  VDBGrid_Ptr->treePtr()->evalActiveVoxelBoundingBox(m_bbox);
50  // Smart pointer to vtkImageData object
51  m_VTKimageData_Ptr = vtkSmartPointer<vtkImageData>::New();
52  m_VTKimageData_Ptr->SetDimensions(m_bbox.dim().x(), m_bbox.dim().y(), m_bbox.dim().z());
53  m_VTKimageData_Ptr->SetExtent(m_bbox.min().x(), m_bbox.max().x(), m_bbox.min().y(), m_bbox.max().y(), m_bbox.min().z(), m_bbox.max().z());
54  m_VTKimageData_Ptr->SetSpacing(voxelSize, voxelSize, voxelSize);
55  m_pointNumber = m_VTKimageData_Ptr->GetNumberOfPoints();
56  m_VDBGrid_Ptr = VDBGrid_Ptr;
57  };

◆ ~VTIWriter()

template<class T_VDBGridType >
MF::RW::VTIWriter< T_VDBGridType >::~VTIWriter ( )
default

Member Function Documentation

◆ addDataFromGrid() [1/2]

template<class T_VDBGridType >
template<typename T3 >
void MF::RW::VTIWriter< T_VDBGridType >::addDataFromGrid ( const std::string &  DataName,
unsigned int  numComponents,
std::shared_ptr< T_VDBGridType >  VDBGrid_Ptr 
)
inline

Definition at line 78 of file VTIWriter.h.

78  {
79  auto newData = vtkSmartPointer<T3>::New();
80  newData->SetNumberOfComponents(numComponents);
81  newData->SetNumberOfTuples(m_VTKimageData_Ptr->GetNumberOfPoints());
82  newData->SetName(DataName.c_str());
83 
84  // Get an accessor for coordinate-based access to voxels.
85  auto accessor = VDBGrid_Ptr->getAccessor();
86  openvdb::Coord ijk;
87  uint64_t nr = 0;
88  double tuple[numComponents];
89  for (unsigned int n=0; n<numComponents; ++n) tuple[n] = NAN;
90  int32_t &i = ijk[0], &j = ijk[1], &k = ijk[2];
91 
92  for (k = m_bbox.min().z(); k <= m_bbox.max().z(); ++k) {
93  for (j = m_bbox.min().y(); j <= m_bbox.max().y(); ++j) {
94  for (i = m_bbox.min().x(); i <= m_bbox.max().x(); ++i) {
95  for(unsigned int n=0; n<numComponents; ++n)
96  tuple[n] = accessor.getValue(ijk);
97  newData->SetTuple(nr, tuple);
98  nr ++;
99  }
100  }
101  }
102  m_VTKimageData_Ptr->GetPointData()->AddArray(newData);
103  };

◆ addDataFromGrid() [2/2]

template<class T_VDBGridType >
template<typename T3 >
void MF::RW::VTIWriter< T_VDBGridType >::addDataFromGrid ( const std::string &  DataName,
unsigned int  numComponents 
)
inline

Definition at line 106 of file VTIWriter.h.

106  {
107  auto newData = vtkSmartPointer<T3>::New();
108  newData->SetNumberOfComponents(numComponents);
109  newData->SetNumberOfTuples(m_pointNumber);
110  newData->SetName(DataName.c_str());
111 
112  // Get an accessor for coordinate-based access to voxels.
113  auto accessor = m_VDBGrid_Ptr->getAccessor();
114  openvdb::Coord ijk;
115  uint64_t nr = 0;
116  double tuple[numComponents]; //?! Why only double or float allowed?
117  for (unsigned int n=0; n<numComponents; ++n)
118  tuple[n] = NAN;
119 
120  int32_t &i = ijk[0], &j = ijk[1], &k = ijk[2];
121 
122  for (k = m_bbox.min().z(); k <= m_bbox.max().z(); ++k) {
123  for (j = m_bbox.min().y(); j <= m_bbox.max().y(); ++j) {
124  for (i = m_bbox.min().x(); i <= m_bbox.max().x(); ++i) {
125  for(unsigned int n=0; n<numComponents; ++n)
126  tuple[n] = accessor.getValue(ijk);
127  newData->SetTuple(nr, tuple);
128  nr ++;
129  }
130  }
131  }
132  m_VTKimageData_Ptr->GetPointData()->AddArray(newData);
133  };

◆ addDataFromTable()

template<class T_VDBGridType >
template<typename T3_VTKArrayType , typename T4 >
void MF::RW::VTIWriter< T_VDBGridType >::addDataFromTable ( const std::string &  DataName,
unsigned int  numComponents,
std::vector< T4 > &  table_Ptr 
)
inline

Definition at line 193 of file VTIWriter.h.

193  {
194  vtkSmartPointer<T3_VTKArrayType> newData = vtkSmartPointer<T3_VTKArrayType>::New();
195  newData->SetNumberOfComponents(numComponents);
196  newData->SetNumberOfTuples(m_pointNumber);
197  newData->SetName(DataName.c_str());
198 
199  double tuple[numComponents];
200  for (unsigned int n=0; n<numComponents; ++n)
201  tuple[n] = NAN;
202 
203  for (auto i = 0; i < m_pointNumber; ++i) {
204  for(unsigned int n=0; n<numComponents; ++n)
205  tuple[n] = table_Ptr[i*numComponents+n];
206  newData->SetTuple(i, tuple);
207  }
208  m_VTKimageData_Ptr->GetPointData()->AddArray(newData);
209  };

◆ addFQ19FromThread()

template<class T_VDBGridType >
void MF::RW::VTIWriter< T_VDBGridType >::addFQ19FromThread ( const std::shared_ptr< MF::Database::ThreadArray > &  ThreadArray_Ptr)
inline

Definition at line 517 of file VTIWriter.h.

References MFQ19.

517  {
518  vtkSmartPointer<vtkDoubleArray> newData = vtkSmartPointer<vtkDoubleArray>::New();
519  newData->SetNumberOfComponents(MFQ19);
520  newData->SetNumberOfTuples(m_pointNumber);
521  newData->SetName("FQ19");
522  const uint32_t Nx = m_bbox.dim().x();
523  const uint32_t Ny = m_bbox.dim().y();
524 
525  double tuple[MFQ19] = {NAN, NAN, NAN, NAN, NAN, NAN, NAN, NAN, NAN, NAN, NAN, NAN, NAN, NAN, NAN, NAN, NAN, NAN, NAN};
526  for (auto i = 0; i < m_pointNumber; ++i) {
527  newData->SetTuple(i, tuple);
528  }
529  for (auto & thread : *ThreadArray_Ptr->m_ThreadsTable_Ptr) {
530  for (size_t node = 0; node < thread->m_NodeArray_Ptr->size(); node++) {
531  for (unsigned int i = 0; i < MFQ19; i++)
532  tuple[i] = thread->m_NodeArray_Ptr->operator[](node).FQ19[i];
533  const uint32_t X = thread->m_Coordinate_Ptr->operator[](node).x - m_bbox.min().x();
534  const uint32_t Y = thread->m_Coordinate_Ptr->operator[](node).y - m_bbox.min().y();
535  const uint32_t Z = thread->m_Coordinate_Ptr->operator[](node).z - m_bbox.min().z();
536  const uint64_t node_nr = (Z * Ny + Y) * Nx + X;
537  newData->SetTuple(node_nr, tuple);
538  }
539  }
540  m_VTKimageData_Ptr->GetPointData()->AddArray(newData);
541  };
#define MFQ19
Number of lattice directions D3Q19.

◆ addInitialRhoLBFromThread()

template<class T_VDBGridType >
void MF::RW::VTIWriter< T_VDBGridType >::addInitialRhoLBFromThread ( const std::shared_ptr< MF::Database::ThreadArray > &  ThreadArray_Ptr)
inline

Definition at line 440 of file VTIWriter.h.

440  {
441  vtkSmartPointer<vtkDoubleArray> newData = vtkSmartPointer<vtkDoubleArray>::New();
442  newData->SetNumberOfComponents(1);
443  newData->SetNumberOfTuples(m_pointNumber);
444  newData->SetName("Initial rho LB [-]");
445  const uint32_t Nx = m_bbox.dim().x();
446  const uint32_t Ny = m_bbox.dim().y();
447 
448  double tuple[1] = {NAN};
449  for (size_t i = 0; i < m_pointNumber; ++i) {
450  newData->SetTuple(i, tuple);
451  }
452  for (auto & thread : *ThreadArray_Ptr->m_ThreadsTable_Ptr) {
453  tuple[0] = thread.get()->m_BoundaryValue_RhoLB;
454  for (size_t node = 0; node < thread->m_NodeArray_Ptr->size(); node++) {
455  const uint32_t X = thread->m_Coordinate_Ptr->operator[](node).x - m_bbox.min().x();
456  const uint32_t Y = thread->m_Coordinate_Ptr->operator[](node).y - m_bbox.min().y();
457  const uint32_t Z = thread->m_Coordinate_Ptr->operator[](node).z - m_bbox.min().z();
458  const uint64_t node_nr = (Z * Ny + Y) * Nx + X;
459  newData->SetTuple(node_nr, tuple);
460  }
461  }
462  m_VTKimageData_Ptr->GetPointData()->AddArray(newData);
463  };

◆ addInitialVelocityLBFromThread()

template<class T_VDBGridType >
void MF::RW::VTIWriter< T_VDBGridType >::addInitialVelocityLBFromThread ( const std::shared_ptr< MF::Database::ThreadArray > &  ThreadArray_Ptr)
inline

Definition at line 359 of file VTIWriter.h.

359  {
360  vtkSmartPointer<vtkDoubleArray> newData = vtkSmartPointer<vtkDoubleArray>::New();
361  newData->SetNumberOfComponents(3);
362  newData->SetNumberOfTuples(m_pointNumber);
363  newData->SetName("Initial velocity LB [-]");
364  const uint32_t Nx = m_bbox.dim().x();
365  const uint32_t Ny = m_bbox.dim().y();
366 
367  double tuple[3] = {NAN,NAN,NAN};
368  for (auto i = 0; i < m_pointNumber; ++i) {
369  newData->SetTuple(i, tuple);
370  }
371  for (auto & thread : *ThreadArray_Ptr->m_ThreadsTable_Ptr) {
372  tuple[0] = thread.get()->m_BoundaryValue_ULB.x;
373  tuple[1] = thread.get()->m_BoundaryValue_ULB.y;
374  tuple[2] = thread.get()->m_BoundaryValue_ULB.z;
375  for (size_t node = 0; node < thread->m_NodeArray_Ptr->size(); node++) {
376  const uint32_t X = thread->m_Coordinate_Ptr->operator[](node).x - m_bbox.min().x();
377  const uint32_t Y = thread->m_Coordinate_Ptr->operator[](node).y - m_bbox.min().y();
378  const uint32_t Z = thread->m_Coordinate_Ptr->operator[](node).z - m_bbox.min().z();
379  const uint64_t node_nr = (Z * Ny + Y) * Nx + X;
380  newData->SetTuple(node_nr, tuple);
381  }
382  }
383  m_VTKimageData_Ptr->GetPointData()->AddArray(newData);
384  };

◆ addNodeCoordinateFromGrid()

template<class T_VDBGridType >
void MF::RW::VTIWriter< T_VDBGridType >::addNodeCoordinateFromGrid ( )
inline

Definition at line 160 of file VTIWriter.h.

160  {
161  auto newData = vtkSmartPointer<vtkIntArray>::New();
162  newData->SetNumberOfComponents(3);
163  newData->SetNumberOfTuples(m_pointNumber);
164  newData->SetName("MF NodeCoord");
165 
166  // Get an accessor for coordinate-based access to voxels.
167  auto accessor = m_VDBGrid_Ptr->getAccessor();
168  openvdb::Coord ijk;
169  uint64_t nr = 0;
170  double tuple[3] ={NAN,NAN,NAN};
171  for (auto i = 0; i < m_pointNumber; ++i) {
172  newData->SetTuple(i, tuple);
173  }
174  int32_t &i = ijk[0], &j = ijk[1], &k = ijk[2];
175 
176  for (k = m_bbox.min().z(); k <= m_bbox.max().z(); ++k) {
177  for (j = m_bbox.min().y(); j <= m_bbox.max().y(); ++j) {
178  for (i = m_bbox.min().x(); i <= m_bbox.max().x(); ++i) {
179  if (accessor.getValue(ijk) != 0) {
180  tuple[0] = i;
181  tuple[1] = j;
182  tuple[2] = k;
183  newData->SetTuple(nr, tuple);
184  }
185  nr++;
186  }
187  }
188  }
189  m_VTKimageData_Ptr->GetPointData()->AddArray(newData);
190  };

◆ addNodeCoordinateFromThread()

template<class T_VDBGridType >
void MF::RW::VTIWriter< T_VDBGridType >::addNodeCoordinateFromThread ( const std::shared_ptr< MF::Database::ThreadArray > &  ThreadArray_Ptr)
inline

Definition at line 211 of file VTIWriter.h.

211  {
212  vtkSmartPointer<vtkIntArray> newData = vtkSmartPointer<vtkIntArray>::New();
213  newData->SetNumberOfComponents(3);
214  newData->SetNumberOfTuples(m_pointNumber);
215  newData->SetName("MF NodeCoord");
216  const uint32_t Nx = m_bbox.dim().x();
217  const uint32_t Ny = m_bbox.dim().y();
218 
219  double tuple[3] = {NAN,NAN,NAN};
220  for (auto i = 0; i < m_pointNumber; ++i) {
221  newData->SetTuple(i, tuple);
222  }
223  for (auto & thread : *ThreadArray_Ptr->m_ThreadsTable_Ptr) {
224  for (size_t node = 0; node < thread->m_NodeArray_Ptr->size(); node++) {
225  tuple[0] = thread->m_Coordinate_Ptr->operator[](node).x;
226  tuple[1] = thread->m_Coordinate_Ptr->operator[](node).y;
227  tuple[2] = thread->m_Coordinate_Ptr->operator[](node).z;
228  const uint32_t X = thread->m_Coordinate_Ptr->operator[](node).x - m_bbox.min().x();
229  const uint32_t Y = thread->m_Coordinate_Ptr->operator[](node).y - m_bbox.min().y();
230  const uint32_t Z = thread->m_Coordinate_Ptr->operator[](node).z - m_bbox.min().z();
231  const uint64_t node_nr = (Z * Ny + Y) * Nx + X;
232  newData->SetTuple(node_nr, tuple);
233  }
234  }
235  m_VTKimageData_Ptr->GetPointData()->AddArray(newData);
236  };

◆ addNodeIDFromGrid()

template<class T_VDBGridType >
void MF::RW::VTIWriter< T_VDBGridType >::addNodeIDFromGrid ( )
inline

Definition at line 135 of file VTIWriter.h.

135  {
136  auto newData = vtkSmartPointer<vtkTypeUInt64Array>::New();
137  newData->SetNumberOfComponents(1);
138  newData->SetNumberOfTuples(m_pointNumber);
139  newData->SetName("Node ID");
140 
141  // Get an accessor for coordinate-based access to voxels.
142  auto accessor = m_VDBGrid_Ptr->getAccessor();
143  openvdb::Coord ijk;
144  uint64_t nr = 0;
145  double tuple[1] = {0};
146  int32_t &i = ijk[0], &j = ijk[1], &k = ijk[2];
147 
148  for (k = m_bbox.min().z(); k <= m_bbox.max().z(); ++k) {
149  for (j = m_bbox.min().y(); j <= m_bbox.max().y(); ++j) {
150  for (i = m_bbox.min().x(); i <= m_bbox.max().x(); ++i) {
151  tuple[0] = accessor.getValue(ijk);
152  newData->SetTuple(nr, tuple);
153  nr ++;
154  }
155  }
156  }
157  m_VTKimageData_Ptr->GetPointData()->AddArray(newData);
158  };

◆ addNodeThreadsByNameThread()

template<class T_VDBGridType >
void MF::RW::VTIWriter< T_VDBGridType >::addNodeThreadsByNameThread ( const std::shared_ptr< MF::Database::ThreadArray > &  ThreadArray_Ptr)
inline

Definition at line 288 of file VTIWriter.h.

288  {
289 
290  auto threadVectorName = std::vector<std::string>();
291  for (auto & thread : *ThreadArray_Ptr->m_ThreadsTable_Ptr) {
292  std::string threadName = thread.get()->m_ThreadName;
293  if (!threadName.empty()) {
294  unsigned int number = 0;
295  for (auto & vec_it : threadVectorName) {
296  if (vec_it == threadName)
297  number++;
298  }
299  if (number == 0)
300  threadVectorName.push_back(threadName);
301  }
302  }
303 
304  for (auto & vec_it : threadVectorName) {
305  vtkSmartPointer<vtkIntArray> newData = vtkSmartPointer<vtkIntArray>::New();
306  newData->SetNumberOfComponents(1);
307  newData->SetNumberOfTuples(m_pointNumber);
308  newData->SetName(vec_it.c_str());
309  const uint32_t Nx = m_bbox.dim().x();
310  const uint32_t Ny = m_bbox.dim().y();
311 
312  double tuple[1] = {0};
313  for (auto i = 0; i < m_pointNumber; ++i)
314  newData->SetTuple(i, tuple);
315  tuple[0] = 1;
316  for (auto & thread : *ThreadArray_Ptr->m_ThreadsTable_Ptr) {
317  if (vec_it == thread.get()->m_ThreadName) {
318  for (size_t node = 0; node < thread->m_NodeArray_Ptr->size(); node++) {
319  const uint32_t X = thread->m_Coordinate_Ptr->operator[](node).x - m_bbox.min().x();
320  const uint32_t Y = thread->m_Coordinate_Ptr->operator[](node).y - m_bbox.min().y();
321  const uint32_t Z = thread->m_Coordinate_Ptr->operator[](node).z - m_bbox.min().z();
322  const uint64_t node_nr = (Z * Ny + Y) * Nx + X;
323  newData->SetTuple(node_nr, tuple);
324  }
325  }
326  }
327  m_VTKimageData_Ptr->GetPointData()->AddArray(newData);
328  }
329  };

◆ addNodeTypeFromThread()

template<class T_VDBGridType >
void MF::RW::VTIWriter< T_VDBGridType >::addNodeTypeFromThread ( const std::shared_ptr< MF::Database::ThreadArray > &  ThreadArray_Ptr)
inline

Definition at line 238 of file VTIWriter.h.

238  {
239  vtkSmartPointer<vtkIntArray> newData = vtkSmartPointer<vtkIntArray>::New();
240  newData->SetNumberOfComponents(1);
241  newData->SetNumberOfTuples(m_pointNumber);
242  newData->SetName("Node type");
243  const uint32_t Nx = m_bbox.dim().x();
244  const uint32_t Ny = m_bbox.dim().y();
245 
246  double tuple[1] = {0};
247  for (auto i = 0; i < m_pointNumber; ++i) {
248  newData->SetTuple(i, tuple);
249  }
250  for (auto & thread : *ThreadArray_Ptr->m_ThreadsTable_Ptr) {
251  tuple[0] = thread.get()->m_NodeType;
252  for (size_t node = 0; node < thread->m_NodeArray_Ptr->size(); node++) {
253  const uint32_t X = thread->m_Coordinate_Ptr->operator[](node).x - m_bbox.min().x();
254  const uint32_t Y = thread->m_Coordinate_Ptr->operator[](node).y - m_bbox.min().y();
255  const uint32_t Z = thread->m_Coordinate_Ptr->operator[](node).z - m_bbox.min().z();
256  const uint64_t node_nr = (Z * Ny + Y) * Nx + X;
257  newData->SetTuple(node_nr, tuple);
258  }
259  }
260  m_VTKimageData_Ptr->GetPointData()->AddArray(newData);
261  };

◆ addPressPhysFromThread()

template<class T_VDBGridType >
void MF::RW::VTIWriter< T_VDBGridType >::addPressPhysFromThread ( const std::shared_ptr< MF::Database::ThreadArray > &  ThreadArray_Ptr)
inline

Definition at line 491 of file VTIWriter.h.

References MF::Database::Vec4< T >::rho, and MF::LBPC::ParametersConversion::RhoLB_toPressurePhys().

491  {
492  vtkSmartPointer<vtkDoubleArray> newData = vtkSmartPointer<vtkDoubleArray>::New();
493  newData->SetNumberOfComponents(1);
494  newData->SetNumberOfTuples(m_pointNumber);
495  newData->SetName("Gauge pressure [Pa]");
496  const uint32_t Nx = m_bbox.dim().x();
497  const uint32_t Ny = m_bbox.dim().y();
498 
499  double tuple[1] = {NAN};
500  for (auto i = 0; i < m_pointNumber; ++i) {
501  newData->SetTuple(i, tuple);
502  }
503  for (auto & thread : *ThreadArray_Ptr->m_ThreadsTable_Ptr) {
504  for (size_t node = 0; node < thread->m_NodeArray_Ptr->size(); node++) {
505  const MF::Database::Vec4<double> Vector4 = (*thread.get()->m_pVRLBFunction)(&thread->m_NodeArray_Ptr->operator[](node));
507  const uint32_t X = thread->m_Coordinate_Ptr->operator[](node).x - m_bbox.min().x();
508  const uint32_t Y = thread->m_Coordinate_Ptr->operator[](node).y - m_bbox.min().y();
509  const uint32_t Z = thread->m_Coordinate_Ptr->operator[](node).z - m_bbox.min().z();
510  const uint64_t node_nr = (Z * Ny + Y) * Nx + X;
511  newData->SetTuple(node_nr, tuple);
512  }
513  }
514  m_VTKimageData_Ptr->GetPointData()->AddArray(newData);
515  };
T rho
Rho value.
Definition: Vec4.h:22
static double RhoLB_toPressurePhys(double rho_LB)
Converts rho LB to physical pressure.

◆ addPropagationDirFromThread()

template<class T_VDBGridType >
void MF::RW::VTIWriter< T_VDBGridType >::addPropagationDirFromThread ( const std::shared_ptr< MF::Database::ThreadArray > &  ThreadArray_Ptr)
inline

Definition at line 543 of file VTIWriter.h.

References MFQ19.

543  {
544  for (uint8_t k = 0; k < MFQ19; k++) {
545  vtkSmartPointer<vtkIntArray> newData = vtkSmartPointer<vtkIntArray>::New();
546  newData->SetNumberOfComponents(3);
547  newData->SetNumberOfTuples(m_pointNumber);
548  const std::string dataName = "PD" + std::to_string(k) + " [x,y,z]";
549  newData->SetName(dataName.c_str());
550  const uint32_t Nx = m_bbox.dim().x();
551  const uint32_t Ny = m_bbox.dim().y();
552  double tuple[3] = {NAN, NAN, NAN};
553  for (auto i = 0; i < m_pointNumber; ++i)
554  newData->SetTuple(i, tuple);
555 
556  for (auto & thread : *ThreadArray_Ptr->m_ThreadsTable_Ptr) {
557  for (size_t node = 0; node < thread->m_NodeArray_Ptr->size(); node++) {
558  if (thread->getLinkedNodePointer(node, k) != nullptr) {
559  auto neighbNodeNr = thread->getLinkedNodePointer(node, k)-> pMyThread_Ptr->operator*().getNodeNr(thread->getLinkedNodePointer(node, k));
560  tuple[0] = thread->getLinkedNodePointer(node, k)-> pMyThread_Ptr->operator*().getNodeCoordX(neighbNodeNr);
561  tuple[1] = thread->getLinkedNodePointer(node, k)-> pMyThread_Ptr->operator*().getNodeCoordY(neighbNodeNr);
562  tuple[2] = thread->getLinkedNodePointer(node, k)-> pMyThread_Ptr->operator*().getNodeCoordZ(neighbNodeNr);
563  const uint32_t X = thread->getNodeCoordX(node) - m_bbox.min().x();
564  const uint32_t Y = thread->getNodeCoordY(node) - m_bbox.min().y();
565  const uint32_t Z = thread->getNodeCoordZ(node) - m_bbox.min().z();
566  const uint64_t node_nr = (Z * Ny + Y) * Nx + X;
567  newData->SetTuple(node_nr, tuple);
568  }
569  }
570  }
571  m_VTKimageData_Ptr->GetPointData()->AddArray(newData);
572  }
573  };
#define MFQ19
Number of lattice directions D3Q19.

◆ addRhoLBFromThread()

template<class T_VDBGridType >
void MF::RW::VTIWriter< T_VDBGridType >::addRhoLBFromThread ( const std::shared_ptr< MF::Database::ThreadArray > &  ThreadArray_Ptr)
inline

Definition at line 414 of file VTIWriter.h.

References MF::Database::Vec4< T >::rho.

414  {
415  vtkSmartPointer<vtkDoubleArray> newData = vtkSmartPointer<vtkDoubleArray>::New();
416  newData->SetNumberOfComponents(1);
417  newData->SetNumberOfTuples(m_pointNumber);
418  newData->SetName("Rho LB [-]");
419  const uint32_t Nx = m_bbox.dim().x();
420  const uint32_t Ny = m_bbox.dim().y();
421 
422  double tuple[1] = {NAN};
423  for (size_t i = 0; i < m_pointNumber; ++i) {
424  newData->SetTuple(i, tuple);
425  }
426  for (auto & thread : *ThreadArray_Ptr->m_ThreadsTable_Ptr) {
427  for (size_t node = 0; node < thread->m_NodeArray_Ptr->size(); node++) {
428  const MF::Database::Vec4<double> Vector4 = (*thread.get()->m_pVRLBFunction)(&thread->m_NodeArray_Ptr->operator[](node));
429  tuple[0] = Vector4.rho;
430  const uint32_t X = thread->m_Coordinate_Ptr->operator[](node).x - m_bbox.min().x();
431  const uint32_t Y = thread->m_Coordinate_Ptr->operator[](node).y - m_bbox.min().y();
432  const uint32_t Z = thread->m_Coordinate_Ptr->operator[](node).z - m_bbox.min().z();
433  const uint64_t node_nr = (Z * Ny + Y) * Nx + X;
434  newData->SetTuple(node_nr, tuple);
435  }
436  }
437  m_VTKimageData_Ptr->GetPointData()->AddArray(newData);
438  };
T rho
Rho value.
Definition: Vec4.h:22

◆ addRhoPhysFromThread()

template<class T_VDBGridType >
void MF::RW::VTIWriter< T_VDBGridType >::addRhoPhysFromThread ( const std::shared_ptr< MF::Database::ThreadArray > &  ThreadArray_Ptr)
inline

Definition at line 465 of file VTIWriter.h.

References MF::Database::Vec4< T >::rho, and MF::LBPC::ParametersConversion::RhoLB_toRhoPhys().

465  {
466  vtkSmartPointer<vtkDoubleArray> newData = vtkSmartPointer<vtkDoubleArray>::New();
467  newData->SetNumberOfComponents(1);
468  newData->SetNumberOfTuples(m_pointNumber);
469  newData->SetName("Rho [kg/m3]");
470  const uint32_t Nx = m_bbox.dim().x();
471  const uint32_t Ny = m_bbox.dim().y();
472 
473  double tuple[1] = {NAN};
474  for (auto i = 0; i < m_pointNumber; ++i) {
475  newData->SetTuple(i, tuple);
476  }
477  for (auto & thread : *ThreadArray_Ptr->m_ThreadsTable_Ptr) {
478  for (size_t node = 0; node < thread->m_NodeArray_Ptr->size(); node++) {
479  const MF::Database::Vec4<double> Vector4 = (*thread.get()->m_pVRLBFunction)(&thread->m_NodeArray_Ptr->operator[](node));
481  const uint32_t X = thread->m_Coordinate_Ptr->operator[](node).x - m_bbox.min().x();
482  const uint32_t Y = thread->m_Coordinate_Ptr->operator[](node).y - m_bbox.min().y();
483  const uint32_t Z = thread->m_Coordinate_Ptr->operator[](node).z - m_bbox.min().z();
484  const uint64_t node_nr = (Z * Ny + Y) * Nx + X;
485  newData->SetTuple(node_nr, tuple);
486  }
487  }
488  m_VTKimageData_Ptr->GetPointData()->AddArray(newData);
489  };
T rho
Rho value.
Definition: Vec4.h:22
static double RhoLB_toRhoPhys(double rho_LB)
Converts rho LB to physical density.

◆ addUidThreadFromThread()

template<class T_VDBGridType >
void MF::RW::VTIWriter< T_VDBGridType >::addUidThreadFromThread ( const std::shared_ptr< MF::Database::ThreadArray > &  ThreadArray_Ptr)
inline

Definition at line 263 of file VTIWriter.h.

263  {
264  vtkSmartPointer<vtkIntArray> newData = vtkSmartPointer<vtkIntArray>::New();
265  newData->SetNumberOfComponents(1);
266  newData->SetNumberOfTuples(m_pointNumber);
267  newData->SetName("Node uidThreadNr");
268  const uint32_t Nx = m_bbox.dim().x();
269  const uint32_t Ny = m_bbox.dim().y();
270 
271  double tuple[1] = {NAN};
272  for (auto i = 0; i < m_pointNumber; ++i) {
273  newData->SetTuple(i, tuple);
274  }
275  for (auto & thread : *ThreadArray_Ptr->m_ThreadsTable_Ptr) {
276  tuple[0] = thread.get()->m_uidThreadNr;
277  for (size_t node = 0; node < thread->m_NodeArray_Ptr->size(); node++) {
278  const uint32_t X = thread->m_Coordinate_Ptr->operator[](node).x - m_bbox.min().x();
279  const uint32_t Y = thread->m_Coordinate_Ptr->operator[](node).y - m_bbox.min().y();
280  const uint32_t Z = thread->m_Coordinate_Ptr->operator[](node).z - m_bbox.min().z();
281  const uint64_t node_nr = (Z * Ny + Y) * Nx + X;
282  newData->SetTuple(node_nr, tuple);
283  }
284  }
285  m_VTKimageData_Ptr->GetPointData()->AddArray(newData);
286  };

◆ addVelocityLBFromThread()

template<class T_VDBGridType >
void MF::RW::VTIWriter< T_VDBGridType >::addVelocityLBFromThread ( const std::shared_ptr< MF::Database::ThreadArray > &  ThreadArray_Ptr)
inline

Definition at line 331 of file VTIWriter.h.

References MF::Database::Vec4< T >::x, MF::Database::Vec4< T >::y, and MF::Database::Vec4< T >::z.

331  {
332  vtkSmartPointer<vtkDoubleArray> newData = vtkSmartPointer<vtkDoubleArray>::New();
333  newData->SetNumberOfComponents(3);
334  newData->SetNumberOfTuples(m_pointNumber);
335  newData->SetName("Velocity LB [-]");
336  const uint32_t Nx = m_bbox.dim().x();
337  const uint32_t Ny = m_bbox.dim().y();
338 
339  double tuple[3] = {NAN,NAN,NAN};
340  for (auto i = 0; i < m_pointNumber; ++i) {
341  newData->SetTuple(i, tuple);
342  }
343  for (auto & thread : *ThreadArray_Ptr->m_ThreadsTable_Ptr) {
344  for (size_t node = 0; node < thread->m_NodeArray_Ptr->size(); node++) {
345  const MF::Database::Vec4<double> Vector4 = (*thread->m_pVRLBFunction)(&thread->m_NodeArray_Ptr->operator[](node));
346  tuple[0] = Vector4.x;
347  tuple[1] = Vector4.y;
348  tuple[2] = Vector4.z;
349  const uint32_t X = thread->m_Coordinate_Ptr->operator[](node).x - m_bbox.min().x();
350  const uint32_t Y = thread->m_Coordinate_Ptr->operator[](node).y - m_bbox.min().y();
351  const uint32_t Z = thread->m_Coordinate_Ptr->operator[](node).z - m_bbox.min().z();
352  const uint64_t node_nr = (Z * Ny + Y) * Nx + X;
353  newData->SetTuple(node_nr, tuple);
354  }
355  }
356  m_VTKimageData_Ptr->GetPointData()->AddArray(newData);
357  };
T z
Z direction value.
Definition: Vec4.h:21
T x
X direction value.
Definition: Vec4.h:19
T y
Y direction value.
Definition: Vec4.h:20

◆ addVelocityPhysFromThread()

template<class T_VDBGridType >
void MF::RW::VTIWriter< T_VDBGridType >::addVelocityPhysFromThread ( const std::shared_ptr< MF::Database::ThreadArray > &  ThreadArray_Ptr)
inline

Definition at line 386 of file VTIWriter.h.

References MF::LBPC::ParametersConversion::ULB_toPhys(), MF::Database::Vec4< T >::x, MF::Database::Vec4< T >::y, and MF::Database::Vec4< T >::z.

386  {
387  vtkSmartPointer<vtkDoubleArray> newData = vtkSmartPointer<vtkDoubleArray>::New();
388  newData->SetNumberOfComponents(3);
389  newData->SetNumberOfTuples(m_pointNumber);
390  newData->SetName("Velocity [m/s]");
391  const uint32_t Nx = m_bbox.dim().x();
392  const uint32_t Ny = m_bbox.dim().y();
393 
394  double tuple[3] = {NAN,NAN,NAN};
395  for (auto i = 0; i < m_pointNumber; ++i) {
396  newData->SetTuple(i, tuple);
397  }
398  for (auto & thread : *ThreadArray_Ptr->m_ThreadsTable_Ptr) {
399  for (size_t node = 0; node < thread->m_NodeArray_Ptr->size(); node++) {
400  const MF::Database::Vec4<double> Vector4 = (*thread.get()->m_pVRLBFunction)(&thread->m_NodeArray_Ptr->operator[](node));
401  tuple[0] = MF::LBPC::ParametersConversion::ULB_toPhys(Vector4.x);
402  tuple[1] = MF::LBPC::ParametersConversion::ULB_toPhys(Vector4.y);
403  tuple[2] = MF::LBPC::ParametersConversion::ULB_toPhys(Vector4.z);
404  const uint32_t X = thread->m_Coordinate_Ptr->operator[](node).x - m_bbox.min().x();
405  const uint32_t Y = thread->m_Coordinate_Ptr->operator[](node).y - m_bbox.min().y();
406  const uint32_t Z = thread->m_Coordinate_Ptr->operator[](node).z - m_bbox.min().z();
407  const uint64_t node_nr = (Z * Ny + Y) * Nx + X;
408  newData->SetTuple(node_nr, tuple);
409  }
410  }
411  m_VTKimageData_Ptr->GetPointData()->AddArray(newData);
412  };
T z
Z direction value.
Definition: Vec4.h:21
static double ULB_toPhys(double u_LB)
Converts LB velocity to physical velocity.
T x
X direction value.
Definition: Vec4.h:19
T y
Y direction value.
Definition: Vec4.h:20

◆ New()

template<class T_VDBGridType >
static std::shared_ptr<MF::RW::VTIWriter<T_VDBGridType> > MF::RW::VTIWriter< T_VDBGridType >::New ( std::shared_ptr< T_VDBGridType >  VDBGrid_Ptr,
const std::string &  FileName,
float  voxelSize 
)
inlinestatic

Definition at line 62 of file VTIWriter.h.

Referenced by MF::Database::GridVDB< T_VDBGridType >::saveAllDataToVTIFile(), and MF::Database::GridVDB< T_VDBGridType >::saveToVTIFile().

62  {
63  auto VTIWriter_Ptr = std::make_shared<MF::RW::VTIWriter<T_VDBGridType>>(VDBGrid_Ptr, FileName, voxelSize);
64  return VTIWriter_Ptr;
65  };

◆ writeFile()

template<class T_VDBGridType >
void MF::RW::VTIWriter< T_VDBGridType >::writeFile ( )
inline

Definition at line 575 of file VTIWriter.h.

575  {
576  vtkSmartPointer<vtkXMLImageDataWriter> writer = vtkSmartPointer<vtkXMLImageDataWriter>::New();
577  writer->SetFileName(m_FileNameExt.c_str());
578  writer->SetInputData(m_VTKimageData_Ptr);
579  writer->Write();
580  std::cout << "VTI file was saved successfully: ----> " << m_FileNameExt << std::endl;
581  }

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