NEMoSys  0.63.0
A modular, extensible resource with robust automated mesh generation, mesh quality analysis, adaptive mesh refinement, and data transfer between arbitrary meshes.
blockMeshGen Class Reference

blockMeshGen <– meshGen <– meshBase This class incorporates mesh generating method of blockMesh utility. More...

Detailed Description

Currently, we are supporting automatic full-hex mesh generation Box, Sphere, Cylinder, and Tapered Cone geometries. Output mesh is written in OpenFOAM polyMesh format.

Definition at line 54 of file blockMeshGen.H.

Public Member Functions

 blockMeshGen ()
 blockMeshGen standard constructor More...
 
 blockMeshGen (blockMeshParams *params)
 blockMeshGen alternate constructor. More...
 
 ~blockMeshGen ()
 blockMeshGen standard desctructor More...
 
int createMeshFromSTL (const char *fname)
 Generates mesh and returns VTK database. More...
 
vtkSmartPointer< vtkDataSet > getDataSet () const
 

Public Attributes

std::unique_ptr< NEM::MSH::geoMeshBasegmData
 

Protected Attributes

vtkSmartPointer< vtkDataSet > dataSet
 

Private Member Functions

void initialize ()
 Initializes OpenFOAM args and runtime. More...
 
void createBlockMshDict (const bool &write)
 Creates blockMeshDict from user arguments. More...
 

Private Attributes

bool defaults
 
blockMeshParamsparams_
 Definition of parameters pointer. More...
 
std::unique_ptr< Foam::Time > runTime_
 
std::unique_ptr< Foam::fvMesh > fmesh_
 
std::unique_ptr< Foam::dictionary > blockMshDict_
 
std::unique_ptr< Foam::dictionary > controlDict_
 
std::unique_ptr< Foam::dictionary > fvSchemes_
 
std::unique_ptr< Foam::dictionary > fvSolution_
 

Inherits meshGen.

Constructor & Destructor Documentation

◆ blockMeshGen() [1/2]

blockMeshGen::blockMeshGen ( )

Definition at line 56 of file blockMeshGen.C.

References defaults, initialize(), and params_.

57 {
58  // booleans
59  params_ = new blockMeshParams();
60  defaults = true;
61 
62  // initializing the Foam environment
63  initialize();
64 }
blockMeshParams * params_
Definition of parameters pointer.
Definition: blockMeshGen.H:97
blockMeshParams contains the parameters important for automatic meshing using blockMeshGen class...
void initialize()
Initializes OpenFOAM args and runtime.
Definition: blockMeshGen.C:75

◆ blockMeshGen() [2/2]

blockMeshGen::blockMeshGen ( blockMeshParams params)

Uses user-defined parameters to perform requested meshing operation.

Parameters
paramsblockMeshParams object

Definition at line 67 of file blockMeshGen.C.

References initialize().

68  : defaults(false), params_(params) {
69  // Initialize foam environment
70  initialize();
71 }
blockMeshParams * params_
Definition of parameters pointer.
Definition: blockMeshGen.H:97
void initialize()
Initializes OpenFOAM args and runtime.
Definition: blockMeshGen.C:75

◆ ~blockMeshGen()

blockMeshGen::~blockMeshGen ( )

Definition at line 73 of file blockMeshGen.C.

73 {}

Member Function Documentation

◆ createBlockMshDict()

void blockMeshGen::createBlockMshDict ( const bool &  write)
private

Definition at line 124 of file blockMeshGen.C.

References blockMshDict_, blockMeshParams::cellSize, blockMeshParams::cnvrtToMeters, blockMeshParams::nCells, params_, points, blockMeshParams::shape, sin45, and tan30.

Referenced by initialize().

124  {
125  // header
126  std::string contText =
127  "\
128 /*--------------------------------*- C++ -*----------------------------------*\n\
129 | ========= | |\n\
130 | \\\\ / F ield | NEMoSys: blockMesh interface |\n\
131 | \\\\ / O peration | |\n\
132 | \\\\ / A nd | |\n\
133 | \\\\/ M anipulation | |\n\
134 \\*---------------------------------------------------------------------------*/\n\
135 \n\
136 FoamFile\n\
137 {\n\
138  version 2.0;\n\
139  format ascii;\n\
140  class dictionary;\n\
141  location \"system\";\n\
142  object blockMeshDict;\n\
143 }\n\n";
144 
145  contText =
146  contText +
147  "// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //\n\n";
148 
149  if (auto box = std::dynamic_pointer_cast<bmBox>(params_->shape)) {
150  double finalX;
151  double finalY;
152  double finalZ;
153 
154  // Automatically Generates Box Around Packs and Meshes that Hexahedrally.
155  if (box->autoGenerate.has_value()) {
156  const auto &autoGenerate = box->autoGenerate.value();
157  Foam::fileName inFileName((autoGenerate.packFileName));
158 
159  Foam::Module::triSurf origSurface(inFileName);
160  Foam::Module::triSurfModifier sMod(origSurface);
161  Foam::pointField &points = sMod.pointsAccess();
162 
163  const Foam::boundBox bb(points);
164 
165  Foam::vector negOffset, posOffset;
166 
167  negOffset[0] = (autoGenerate.offset[0]);
168  negOffset[1] = (autoGenerate.offset[1]);
169  negOffset[2] = (autoGenerate.offset[2]);
170  posOffset[0] = (autoGenerate.offset[0]);
171  posOffset[1] = (autoGenerate.offset[1]);
172  posOffset[2] = (autoGenerate.offset[2]);
173 
174  const Foam::boundBox newBB(bb.min() - negOffset, bb.max() + posOffset);
175 
176  (box->init[0]) = newBB.min()[0];
177  (box->init[1]) = newBB.min()[1];
178  (box->init[2]) = newBB.min()[2];
179 
180  finalX = newBB.max()[0];
181  finalY = newBB.max()[1];
182  finalZ = newBB.max()[2];
183 
184  std::array<double, 3> minPointsBox{newBB.min()[0], newBB.min()[1],
185  newBB.min()[2]};
186  std::array<double, 3> maxPointsBox{newBB.max()[0], newBB.max()[1],
187  newBB.max()[2]};
188 
189  box->coordsBox = std::make_pair(minPointsBox, maxPointsBox);
190 
191  if (params_->cellSize.has_value()) {
192  double xLength = std::sqrt(((box->init[0]) - (finalX)) *
193  ((box->init[0]) - (finalX)));
194  double yLength = std::sqrt(((box->init[1]) - (finalY)) *
195  ((box->init[1]) - (finalY)));
196  double zLength = std::sqrt(((box->init[2]) - (finalZ)) *
197  ((box->init[2]) - (finalZ)));
198 
199  (params_->nCells[0]) =
200  static_cast<int>(std::round(xLength / (params_->cellSize.value())));
201  (params_->nCells[1]) =
202  static_cast<int>(std::round(yLength / (params_->cellSize.value())));
203  (params_->nCells[2]) =
204  static_cast<int>(std::round(zLength / (params_->cellSize.value())));
205  }
206  } else {
207  finalX = (box->init[0]) + (box->len[0]);
208  finalY = (box->init[1]) + (box->len[1]);
209  finalZ = (box->init[2]) + (box->len[2]);
210 
211  std::array<double, 3> minPointsBox{box->init[0], box->init[1],
212  box->init[2]};
213  std::array<double, 3> maxPointsBox{finalX, finalY, finalZ};
214  box->coordsBox = std::make_pair(minPointsBox, maxPointsBox);
215  }
216 
217  // Box data
218  contText = contText + "convertToMeters " +
219  std::to_string(params_->cnvrtToMeters) + ";\n";
220  contText = contText + "\nvertices\n";
221  contText = contText + "(\n\n";
222 
223  contText = contText + "\t(" + std::to_string(box->init[0]) + " " +
224  std::to_string(box->init[1]) + " " +
225  std::to_string(box->init[2]) + ")\n";
226 
227  contText = contText + "\t(" + std::to_string(finalX) + " " +
228  std::to_string(box->init[1]) + " " +
229  std::to_string(box->init[2]) + ")\n";
230 
231  contText = contText + "\t(" + std::to_string(finalX) + " " +
232  std::to_string(finalY) + " " + std::to_string(box->init[2]) +
233  ")\n";
234 
235  contText = contText + "\t(" + std::to_string(box->init[0]) + " " +
236  std::to_string(finalY) + " " + std::to_string(box->init[2]) +
237  ")\n";
238 
239  contText = contText + "\t(" + std::to_string(box->init[0]) + " " +
240  std::to_string(box->init[1]) + " " + std::to_string(finalZ) +
241  ")\n";
242 
243  contText = contText + "\t(" + std::to_string(finalX) + " " +
244  std::to_string(box->init[1]) + " " + std::to_string(finalZ) +
245  ")\n";
246 
247  contText = contText + "\t(" + std::to_string(finalX) + " " +
248  std::to_string(finalY) + " " + std::to_string(finalZ) + ")\n";
249 
250  contText = contText + "\t(" + std::to_string(box->init[0]) + " " +
251  std::to_string(finalY) + " " + std::to_string(finalZ) + ")\n";
252  contText = contText + "\n);\n";
253 
254  contText = contText + "\nblocks\n(";
255  contText = contText + "\n\thex (0 1 2 3 4 5 6 7) (" +
256  std::to_string(params_->nCells[0]) + " " +
257  std::to_string(params_->nCells[1]) + " " +
258  std::to_string(params_->nCells[2]) + ") simpleGrading (" +
259  std::to_string(box->smplGrading[0]) + " " +
260  std::to_string(box->smplGrading[1]) + " " +
261  std::to_string(box->smplGrading[2]) + ")\n\n";
262  contText = contText + ");\n";
263 
264  contText = contText + "\nedges\n";
265  contText = contText + "(\n\n);\n";
266 
267  contText = contText + "\npatches";
268  contText = contText + "\n(\n";
269  contText = contText + "\n\tpatch up\n";
270  contText = contText + "\t(\n";
271  contText = contText + "\t\t(3 7 6 2)\n";
272  contText = contText + "\t)\n";
273 
274  contText = contText + "\tpatch down\n";
275  contText = contText + "\t(\n";
276  contText = contText + "\t\t(0 4 5 1)\n";
277  contText = contText + "\t)\n";
278 
279  contText = contText + "\tpatch in\n";
280  contText = contText + "\t(\n";
281  contText = contText + "\t\t(3 7 4 0)\n";
282  contText = contText + "\t)\n";
283 
284  contText = contText + "\tpatch out\n";
285  contText = contText + "\t(\n";
286  contText = contText + "\t\t(2 6 5 1)\n";
287  contText = contText + "\t)\n";
288 
289  contText = contText + "\tpatch front\n";
290  contText = contText + "\t(\n";
291  contText = contText + "\t\t(0 1 2 3)\n";
292  contText = contText + "\t)\n";
293 
294  contText = contText + "\tpatch back\n";
295  contText = contText + "\t(\n";
296  contText = contText + "\t\t(4 5 6 7)\n";
297  contText = contText + "\t)\n";
298 
299  contText = contText + "\n);\n";
300 
301  contText = contText + "\nmergePatchPairs\n(\n);\n";
302  contText = contText +
303  "\n //* * * * * * * * * * * * * * * * * * * * * * * * * * * * * "
304  "* * *//";
305 
306  }
307 
308  else if (auto sphere = std::dynamic_pointer_cast<bmSphere>(params_->shape)) {
309  // Put sphere implementation here
310  contText = contText + "\ncx " + std::to_string(sphere->center[0]) + ";";
311  contText = contText + "\ncy " + std::to_string(sphere->center[1]) + ";";
312  contText = contText + "\ncz " + std::to_string(sphere->center[2]) + ";";
313  contText = contText + "\nrad " + std::to_string(sphere->radius) + ";";
314 
315  contText = contText + "\ngeometry\n{\n";
316  contText = contText + "\tsphere\n\t{\n";
317  contText = contText + "\t\ttype searchableSphere;\n";
318  contText = contText + "\t\tcentre ($cx $cy $cz);\n";
319  contText = contText + "\t\tradius $rad;\n";
320  contText = contText + "\t}\n}\n";
321 
322  contText = contText + "\nconvertToMeters " +
323  std::to_string(params_->cnvrtToMeters) + ";\n";
324 
325  double vx = (sphere->center[0]) + tan30 * (sphere->radius);
326  double mvx = (sphere->center[0]) - tan30 * (sphere->radius);
327  contText = contText + "\nvx " + std::to_string(vx) + ";\n";
328  contText = contText + "mvx " + std::to_string(mvx) + ";\n";
329 
330  double vy = (sphere->center[1]) + tan30 * (sphere->radius);
331  double mvy = (sphere->center[1]) - tan30 * (sphere->radius);
332  contText = contText + "\nvy " + std::to_string(vy) + ";\n";
333  contText = contText + "mvy " + std::to_string(mvy) + ";\n";
334 
335  double vz = (sphere->center[2]) + tan30 * (sphere->radius);
336  double mvz = (sphere->center[2]) - tan30 * (sphere->radius);
337  contText = contText + "\nvz " + std::to_string(vz) + ";\n";
338  contText = contText + "mvz " + std::to_string(mvz) + ";\n";
339 
340  double ax = (sphere->center[0]) + sin45 * (sphere->radius);
341  double max = (sphere->center[0]) - sin45 * (sphere->radius);
342  contText = contText + "\nax " + std::to_string(ax) + ";\n";
343  contText = contText + "max " + std::to_string(max) + ";\n";
344 
345  double ay = (sphere->center[1]) + sin45 * (sphere->radius);
346  double may = (sphere->center[1]) - sin45 * (sphere->radius);
347  contText = contText + "\nay " + std::to_string(ay) + ";\n";
348  contText = contText + "may " + std::to_string(may) + ";\n";
349 
350  double az = (sphere->center[2]) + sin45 * (sphere->radius);
351  double naz = (sphere->center[2]) - sin45 * (sphere->radius);
352  contText = contText + "\naz " + std::to_string(az) + ";\n";
353  contText = contText + "maz " + std::to_string(naz) + ";\n";
354 
355  contText = contText + "\nvertices\n";
356  contText = contText + "(\n";
357  contText = contText + "\t($mvx $mvy $mvz)\n";
358  contText = contText + "\t( $vx $mvy $mvz)\n";
359  contText = contText + "\t( $vx $vy $mvz)\n";
360  contText = contText + "\t($mvx $vy $mvz)\n";
361  contText = contText + "\t($mvx $mvy $vz)\n";
362  contText = contText + "\t( $vx $mvy $vz)\n";
363  contText = contText + "\t( $vx $vy $vz)\n";
364  contText = contText + "\t($mvx $vy $vz)\n";
365  contText = contText + ");\n";
366 
367  contText = contText + "\nblocks\n";
368  contText = contText + "(\n";
369  contText = contText + "\thex (0 1 2 3 4 5 6 7) (" +
370  std::to_string(params_->nCells[0]) + " " +
371  std::to_string(params_->nCells[1]) + " " +
372  std::to_string(params_->nCells[2]) + ") simpleGrading (" +
373  std::to_string(sphere->sphrGrading[0]) + " " +
374  std::to_string(sphere->sphrGrading[1]) + " " +
375  std::to_string(sphere->sphrGrading[2]) + ")\n\n";
376 
377  contText = contText + ");\n";
378 
379  contText = contText + "\nedges\n";
380  contText = contText + "(\n";
381  contText = contText + "\tarc 0 1 ($cx $may $maz)\n";
382  contText = contText + "\tarc 2 3 ($cx $ay $maz)\n";
383  contText = contText + "\tarc 6 7 ($cx $ay $az)\n";
384  contText = contText + "\tarc 4 5 ($cx $may $az)\n";
385 
386  contText = contText + "\tarc 0 3 ($max $cy $maz)\n";
387  contText = contText + "\tarc 1 2 ($ax $cy $maz)\n";
388  contText = contText + "\tarc 5 6 ($ax $cy $az)\n";
389  contText = contText + "\tarc 4 7 ($max $cy $az)\n";
390 
391  contText = contText + "\tarc 0 4 ($max $may $cz)\n";
392  contText = contText + "\tarc 1 5 ($ax $may $cz)\n";
393  contText = contText + "\tarc 2 6 ($ax $ay $cz)\n";
394  contText = contText + "\tarc 3 7 ($max $ay $cz)\n";
395 
396  contText = contText + ");\n";
397 
398  contText = contText + "\nfaces\n(\n";
399  contText = contText + "\tproject (0 4 7 3) sphere\n";
400  contText = contText + "\tproject (2 6 5 1) sphere\n";
401  contText = contText + "\tproject (1 5 4 0) sphere\n";
402  contText = contText + "\tproject (3 7 6 2) sphere\n";
403  contText = contText + "\tproject (0 3 2 1) sphere\n";
404  contText = contText + "\tproject (4 5 6 7) sphere\n";
405  contText = contText + ");\n";
406 
407  contText = contText + "\nboundary\n(\n";
408  contText = contText + "\twalls\n\t{\n";
409  contText = contText + "\t\ttype wall;\n";
410  contText = contText + "\t\tfaces\n\t\t(\n";
411  contText = contText + "\t\t\t(0 4 7 3)\n";
412  contText = contText + "\t\t\t(2 6 5 1)\n";
413  contText = contText + "\t\t\t(1 5 4 0)\n";
414  contText = contText + "\t\t\t(3 7 6 2)\n";
415  contText = contText + "\t\t\t(0 3 2 1)\n";
416  contText = contText + "\t\t\t(4 5 6 7)\n";
417  contText = contText + "\t\t);\n";
418  contText = contText + "\t}\n";
419  contText = contText + ");\n";
420  contText = contText +
421  "\n //* * * * * * * * * * * * * * * * * * * * * * * * * * * * * "
422  "* * *//";
423  } else if (auto cylTaperedCone =
424  std::dynamic_pointer_cast<bmCylTaperedCone>(params_->shape)) {
425  contText = contText + "\n\ncx " +
426  std::to_string(cylTaperedCone->centerCyl[0]) + ";\n";
427  contText =
428  contText + "cy " + std::to_string(cylTaperedCone->centerCyl[1]) + ";\n";
429  contText =
430  contText + "cz " + std::to_string(cylTaperedCone->centerCyl[2]) + ";\n";
431  contText =
432  contText + "rad1 " + std::to_string(cylTaperedCone->radius1) + ";\n";
433  contText = contText + "rad2 " +
434  std::to_string(
435  cylTaperedCone->radius2.value_or(cylTaperedCone->radius1)) +
436  ";\n";
437  contText = contText + "h " + std::to_string(cylTaperedCone->height) + ";\n";
438  contText = contText + "cellX " + std::to_string(params_->nCells[0]) + ";\n";
439  contText = contText + "cellY " + std::to_string(params_->nCells[1]) + ";\n";
440  contText = contText + "cellZ " + std::to_string(params_->nCells[2]) + ";\n";
441  contText = contText + "grdX " +
442  std::to_string(cylTaperedCone->cylGrading[0]) + ";\n";
443  contText = contText + "grdY " +
444  std::to_string(cylTaperedCone->cylGrading[1]) + ";\n";
445  contText = contText + "grdZ " +
446  std::to_string(cylTaperedCone->cylGrading[2]) + ";\n\n";
447 
448  double rad1 = cylTaperedCone->radius1;
449  double rad2 = cylTaperedCone->radius2.value_or(cylTaperedCone->radius1);
450  double X0 = (cylTaperedCone->centerCyl[0]) - sin45 * rad1;
451  contText = contText + "\nX0 " + std::to_string(X0) + ";";
452 
453  double Y0 = (cylTaperedCone->centerCyl[1]) - sin45 * rad1;
454  contText = contText + "\nY0 " + std::to_string(Y0) + ";";
455 
456  double X1 = (cylTaperedCone->centerCyl[0]) + sin45 * rad1;
457  contText = contText + "\nX1 " + std::to_string(X1) + ";";
458 
459  double Y1 = (cylTaperedCone->centerCyl[1]) - sin45 * rad1;
460  contText = contText + "\nY1 " + std::to_string(Y1) + ";";
461 
462  double X2 = (cylTaperedCone->centerCyl[0]) - sin45 * 0.3 * rad1;
463  contText = contText + "\nX2 " + std::to_string(X2) + ";";
464 
465  double Y2 = (cylTaperedCone->centerCyl[1]) - sin45 * 0.3 * rad1;
466  contText = contText + "\nY2 " + std::to_string(Y2) + ";";
467 
468  double X3 = (cylTaperedCone->centerCyl[0]) + sin45 * 0.3 * rad1;
469  contText = contText + "\nX3 " + std::to_string(X3) + ";";
470 
471  double Y3 = (cylTaperedCone->centerCyl[1]) - sin45 * 0.3 * rad1;
472  contText = contText + "\nY3 " + std::to_string(Y3) + ";";
473 
474  double X4 = (cylTaperedCone->centerCyl[0]) + sin45 * 0.3 * rad1;
475  contText = contText + "\nX4 " + std::to_string(X4) + ";";
476 
477  double Y4 = (cylTaperedCone->centerCyl[1]) + sin45 * 0.3 * rad1;
478  contText = contText + "\nY4 " + std::to_string(Y4) + ";";
479 
480  double X5 = (cylTaperedCone->centerCyl[0]) + sin45 * rad1;
481  contText = contText + "\nX5 " + std::to_string(X5) + ";";
482 
483  double Y5 = (cylTaperedCone->centerCyl[1]) + sin45 * rad1;
484  contText = contText + "\nY5 " + std::to_string(Y5) + ";";
485 
486  double X6 = (cylTaperedCone->centerCyl[0]) - sin45 * 0.3 * rad1;
487  contText = contText + "\nX6 " + std::to_string(X6) + ";";
488 
489  double Y6 = (cylTaperedCone->centerCyl[1]) + sin45 * 0.3 * rad1;
490  contText = contText + "\nY6 " + std::to_string(Y6) + ";";
491 
492  double X7 = (cylTaperedCone->centerCyl[0]) - sin45 * rad1;
493  contText = contText + "\nX7 " + std::to_string(X7) + ";";
494 
495  double Y7 = (cylTaperedCone->centerCyl[1]) + sin45 * rad1;
496  contText = contText + "\nY7 " + std::to_string(Y7) + ";";
497 
498  double X8 = (cylTaperedCone->centerCyl[0]) - sin45 * rad2;
499  contText = contText + "\nX8 " + std::to_string(X8) + ";";
500 
501  double Y8 = (cylTaperedCone->centerCyl[1]) - sin45 * rad2;
502  contText = contText + "\nY8 " + std::to_string(Y8) + ";";
503 
504  double X9 = (cylTaperedCone->centerCyl[0]) + sin45 * rad2;
505  contText = contText + "\nX9 " + std::to_string(X9) + ";";
506 
507  double Y9 = (cylTaperedCone->centerCyl[1]) - sin45 * rad2;
508  contText = contText + "\nY9 " + std::to_string(Y9) + ";";
509 
510  double X10 = (cylTaperedCone->centerCyl[0]) - sin45 * 0.3 * rad2;
511  contText = contText + "\nX10 " + std::to_string(X10) + ";";
512 
513  double Y10 = (cylTaperedCone->centerCyl[1]) - sin45 * 0.3 * rad2;
514  contText = contText + "\nY10 " + std::to_string(Y10) + ";";
515 
516  double X11 = (cylTaperedCone->centerCyl[0]) + sin45 * 0.3 * rad2;
517  contText = contText + "\nX11 " + std::to_string(X11) + ";";
518 
519  double Y11 = (cylTaperedCone->centerCyl[1]) - sin45 * 0.3 * rad2;
520  contText = contText + "\nY11 " + std::to_string(Y11) + ";";
521 
522  double X12 = (cylTaperedCone->centerCyl[0]) + sin45 * 0.3 * rad2;
523  contText = contText + "\nX12 " + std::to_string(X12) + ";";
524 
525  double Y12 = (cylTaperedCone->centerCyl[1]) + sin45 * 0.3 * rad2;
526  contText = contText + "\nY12 " + std::to_string(Y12) + ";";
527 
528  double X13 = (cylTaperedCone->centerCyl[0]) + sin45 * rad2;
529  contText = contText + "\nX13 " + std::to_string(X13) + ";";
530 
531  double Y13 = (cylTaperedCone->centerCyl[1]) + sin45 * rad2;
532  contText = contText + "\nY13 " + std::to_string(Y13) + ";";
533 
534  double X14 = (cylTaperedCone->centerCyl[0]) - sin45 * 0.3 * rad2;
535  contText = contText + "\nX14 " + std::to_string(X14) + ";";
536 
537  double Y14 = (cylTaperedCone->centerCyl[1]) + sin45 * 0.3 * rad2;
538  contText = contText + "\nY14 " + std::to_string(Y14) + ";";
539 
540  double X15 = (cylTaperedCone->centerCyl[0]) - sin45 * rad2;
541  contText = contText + "\nX15 " + std::to_string(X15) + ";";
542 
543  double Y15 = (cylTaperedCone->centerCyl[1]) + sin45 * rad2;
544  contText = contText + "\nY15 " + std::to_string(Y15) + ";";
545 
546  double arcX1 = cylTaperedCone->centerCyl[0];
547  contText = contText + "\narcX1 " + std::to_string(arcX1) + ";";
548 
549  double arcY1 = (cylTaperedCone->centerCyl[1] - rad1);
550  contText = contText + "\narcY1 " + std::to_string(arcY1) + ";";
551 
552  double arcX2 = (cylTaperedCone->centerCyl[0] + rad1);
553  contText = contText + "\narcX2 " + std::to_string(arcX2) + ";";
554 
555  double arcY2 = cylTaperedCone->centerCyl[1];
556  contText = contText + "\narcY2 " + std::to_string(arcY2) + ";";
557 
558  double arcX3 = cylTaperedCone->centerCyl[0];
559  contText = contText + "\narcX3 " + std::to_string(arcX3) + ";";
560 
561  double arcY3 = (cylTaperedCone->centerCyl[1] + rad1);
562  contText = contText + "\narcY3 " + std::to_string(arcY3) + ";";
563 
564  double arcX4 = (cylTaperedCone->centerCyl[0] - rad1);
565  contText = contText + "\narcX4 " + std::to_string(arcX4) + ";";
566 
567  double arcY4 = cylTaperedCone->centerCyl[1];
568  contText = contText + "\narcY4 " + std::to_string(arcY4) + ";";
569 
570  double arcX5 = cylTaperedCone->centerCyl[0];
571  contText = contText + "\narcX5 " + std::to_string(arcX5) + ";";
572 
573  double arcY5 = (cylTaperedCone->centerCyl[1] - rad2);
574  contText = contText + "\narcY5 " + std::to_string(arcY5) + ";";
575 
576  double arcX6 = (cylTaperedCone->centerCyl[0] + rad2);
577  contText = contText + "\narcX6 " + std::to_string(arcX6) + ";";
578 
579  double arcY6 = cylTaperedCone->centerCyl[1];
580  contText = contText + "\narcY6 " + std::to_string(arcY6) + ";";
581 
582  double arcX7 = cylTaperedCone->centerCyl[0];
583  contText = contText + "\narcX7 " + std::to_string(arcX7) + ";";
584 
585  double arcY7 = (cylTaperedCone->centerCyl[1] + rad2);
586  contText = contText + "\narcY7 " + std::to_string(arcY7) + ";";
587 
588  double arcX8 = (cylTaperedCone->centerCyl[0] - rad2);
589  contText = contText + "\narcX8 " + std::to_string(arcX8) + ";";
590 
591  double arcY8 = cylTaperedCone->centerCyl[1];
592  contText = contText + "\narcY8 " + std::to_string(arcY8) + ";";
593 
594  double fz = (cylTaperedCone->centerCyl[1] + (cylTaperedCone->height));
595  contText = contText + "\nfz " + std::to_string(fz) + ";";
596 
597  contText = contText + "\nconvertToMeters " +
598  std::to_string(params_->cnvrtToMeters) + ";\n";
599 
600  contText = contText + "\nvertices\n(\n\n";
601  contText = contText + "\t($X0 $Y0 $cz)\n";
602  contText = contText + "\t($X1 $Y1 $cz)\n";
603  contText = contText + "\t($X2 $Y2 $cz)\n";
604  contText = contText + "\t($X3 $Y3 $cz)\n";
605  contText = contText + "\t($X4 $X4 $cz)\n";
606  contText = contText + "\t($X5 $Y5 $cz)\n";
607  contText = contText + "\t($X6 $Y6 $cz)\n";
608  contText = contText + "\t($X7 $Y7 $cz)\n";
609  contText = contText + "\t($X8 $Y8 $fz)\n";
610  contText = contText + "\t($X9 $Y9 $fz)\n";
611  contText = contText + "\t($X10 $Y10 $fz)\n";
612  contText = contText + "\t($X11 $Y11 $fz)\n";
613  contText = contText + "\t($X12 $X12 $fz)\n";
614  contText = contText + "\t($X13 $Y13 $fz)\n";
615  contText = contText + "\t($X14 $Y14 $fz)\n";
616  contText = contText + "\t($X15 $Y15 $fz)\n";
617  contText = contText + ");\n\n";
618 
619  contText = contText + "\nblocks\n(\n";
620  contText = contText +
621  "\thex (2 3 4 6 10 11 12 14) ($cellX $cellY $cellZ) " +
622  "simpleGrading ($grdX $grdY $grdZ)\n";
623  contText = contText + "\thex (0 1 3 2 8 9 11 10) ($cellX $cellY $cellZ) " +
624  "simpleGrading ($grdX $grdY $grdZ)\n";
625  contText = contText + "\thex (3 1 5 4 11 9 13 12) ($cellX $cellY $cellZ) " +
626  "simpleGrading ($grdX $grdY $grdZ)\n";
627  contText = contText +
628  "\thex (4 5 7 6 12 13 15 14) ($cellX $cellY $cellZ) " +
629  "simpleGrading ($grdX $grdY $grdZ)\n";
630  contText = contText + "\thex (6 7 0 2 14 15 8 10) ($cellX $cellY $cellZ) " +
631  "simpleGrading ($grdX $grdY $grdZ)\n";
632  contText = contText + ");\n";
633 
634  contText = contText + "\nedges\n(\n";
635  contText = contText + "\tarc 0 1 ($arcX1 $arcY1 $cz)\n";
636  contText = contText + "\tarc 1 5 ($arcX2 $arcY2 $cz)\n";
637  contText = contText + "\tarc 5 7 ($arcX3 $arcY3 $cz)\n";
638  contText = contText + "\tarc 7 0 ($arcX4 $arcY4 $cz)\n";
639  contText = contText + "\tarc 8 9 ($arcX5 $arcY5 $fz)\n";
640  contText = contText + "\tarc 9 13 ($arcX6 $arcY6 $fz)\n";
641  contText = contText + "\tarc 13 15 ($arcX7 $arcY7 $fz)\n";
642  contText = contText + "\tarc 15 8 ($arcX8 $arcY8 $fz)\n";
643  contText = contText + ");\n\n";
644 
645  contText = contText + "boundary\n(\n";
646  contText = contText + "\tinlet\n\t{\n";
647  contText = contText + "\t\ttype patch;\n";
648  contText = contText + "\t\tfaces\n\t\t(\n";
649  contText = contText + "\t\t\t(0 1 2 3)\n";
650  contText = contText + "\t\t\t(1 3 4 5)\n";
651  contText = contText + "\t\t\t(4 5 7 6)\n";
652  contText = contText + "\t\t\t(0 2 6 7)\n";
653  contText = contText + "\t\t\t(2 3 4 6)\n";
654  contText = contText + "\t\t);\n";
655  contText = contText + "\t}\n";
656 
657  contText = contText + "\toutlet\n\t{\n";
658  contText = contText + "\t\ttype patch;\n";
659  contText = contText + "\t\tfaces\n\t\t(\n";
660  contText = contText + "\t\t\t(8 9 10 11)\n";
661  contText = contText + "\t\t\t(11 9 12 13)\n";
662  contText = contText + "\t\t\t(12 13 14 15)\n";
663  contText = contText + "\t\t\t(14 15 10 8)\n";
664  contText = contText + "\t\t\t(10 11 14 12)\n";
665  contText = contText + "\t\t);\n";
666  contText = contText + "\t}\n";
667 
668  contText = contText + "\twalls\n\t{\n";
669  contText = contText + "\t\ttype patch;\n";
670  contText = contText + "\t\tfaces\n\t\t(\n";
671  contText = contText + "\t\t\t(1 5 13 9)\n";
672  contText = contText + "\t\t\t(5 7 15 13)\n";
673  contText = contText + "\t\t\t(7 0 8 15)\n";
674  contText = contText + "\t\t\t(0 1 9 8)\n";
675  contText = contText + "\t\t);\n";
676  contText = contText + "\t}\n";
677  contText = contText + ");\n";
678 
679  contText =
680  contText +
681  "\n //* * * * * * * * * * * * * * * * * * * * * * * * * * * * *//";
682  } else {
683  std::cerr << "User cannot select multiple geometries in single run!\n"
684  << std::endl;
685  exit(1);
686  }
687 
688  // Keep this dictionary in memory
689  Foam::dictionary tmptmpDuc =
690  new Foam::dictionary(Foam::IStringStream(contText)(), true);
691  blockMshDict_ =
692  std::unique_ptr<Foam::dictionary>(new Foam::dictionary("blockMeshDict"));
693  blockMshDict_->merge(tmptmpDuc);
694 
695  // Write in traditional way due to formatting issues
696  if (write) {
697  // creating a base system directory
698  const char dir_path[] = "./system";
699  boost::filesystem::path dir(dir_path);
700  try {
701  boost::filesystem::create_directory(dir);
702  } catch (boost::filesystem::filesystem_error &e) {
703  std::cerr << "Problem in creating system directory for the cfMesh"
704  << "\n";
705  std::cerr << e.what() << std::endl;
706  throw;
707  }
708  std::ofstream contDict;
709  contDict.open(std::string(dir_path) + "/blockMeshDict");
710  contDict << contText;
711  contDict.close();
712  }
713 }
std::array< int, 3 > nCells
Defines number of cells in X/Y/Z direction; unused if cellSize has a value.
double cnvrtToMeters
Defines mesh scale to meters.
std::unique_ptr< Foam::dictionary > blockMshDict_
Definition: blockMeshGen.H:103
blockMeshParams * params_
Definition of parameters pointer.
Definition: blockMeshGen.H:97
jsoncons::optional< double > cellSize
user can define mesh cell size here (if not set, nCells used instead).
constexpr double sin45
Definition: blockMeshGen.C:54
std::shared_ptr< bmShape > shape
3D block, sphere, cylinder, or tapered cone
std::vector< vtkIdType > points
points given by id in .inp file
Definition: inpGeoMesh.C:133
constexpr double tan30
Definition: blockMeshGen.C:53

◆ createMeshFromSTL()

int blockMeshGen::createMeshFromSTL ( const char *  fname)
virtual
Parameters
fnameInput surface file name
Returns
meshBase mesh dataset

Implements meshGen.

Definition at line 98 of file blockMeshGen.C.

References blockMshDict_, fmesh_, meshGen::gmData, blockMeshParams::isPackMesh, params_, NEM::MSH::Read(), and runTime_.

Referenced by NEM::DRV::BlockMeshMeshGenDriver::execute().

98  {
99  auto mshObj = blockMsh();
100 
101  bool writeMsh;
102  if (params_->isPackMesh)
103  writeMsh = true;
104  else
105  writeMsh = false;
106 
107  fmesh_ = std::unique_ptr<Foam::fvMesh>(
108  mshObj.generate(runTime_, blockMshDict_, writeMsh));
109 
110  // Create foamGeoMesh
111  if (writeMsh) {
112  gmData.reset();
113  gmData = std::unique_ptr<NEM::MSH::geoMeshBase>(NEM::MSH::Read(".foam"));
114  } else {
115  auto fgm_ = std::unique_ptr<NEM::MSH::foamGeoMesh>(
116  new NEM::MSH::foamGeoMesh(fmesh_.get(), ""));
117  gmData = std::unique_ptr<NEM::MSH::geoMeshBase>(
118  dynamic_cast<NEM::MSH::geoMeshBase *>(fgm_.get()));
119  }
120 
121  return 0;
122 }
std::unique_ptr< Foam::fvMesh > fmesh_
Definition: blockMeshGen.H:102
std::unique_ptr< NEM::MSH::geoMeshBase > gmData
Definition: meshGen.H:54
std::unique_ptr< Foam::Time > runTime_
Definition: blockMeshGen.H:101
geoMeshBase * Read(const std::string &fileName)
Read a mesh from file.
std::unique_ptr< Foam::dictionary > blockMshDict_
Definition: blockMeshGen.H:103
blockMeshParams * params_
Definition of parameters pointer.
Definition: blockMeshGen.H:97
bool isPackMesh
Boolean for if the operation is packmesh.
A concrete implementation of geoMeshBase representing a mesh in a fvMesh.
Definition: foamGeoMesh.H:52
abstract class to specify geometry and mesh data
Definition: geoMeshBase.H:102

◆ getDataSet()

vtkSmartPointer<vtkDataSet> meshGen::getDataSet ( ) const
inlineinherited

Definition at line 53 of file meshGen.H.

Referenced by NEM::DRV::MeshGenDriver::MeshGenDriver().

53 { return dataSet; }
vtkSmartPointer< vtkDataSet > dataSet
Definition: meshGen.H:57

◆ initialize()

void blockMeshGen::initialize ( )
private

Definition at line 75 of file blockMeshGen.C.

References controlDict_, createBlockMshDict(), fvSchemes_, fvSolution_, blockMeshParams::isPackMesh, params_, and runTime_.

Referenced by blockMeshGen().

75  {
76  // create dictionaries needed in memory
77  bool writeDicts;
78  if (params_->isPackMesh)
79  writeDicts = true;
80  else
81  writeDicts = false;
82  std::unique_ptr<getDicts> initFoam;
83  initFoam = std::unique_ptr<getDicts>(new getDicts());
84  controlDict_ = initFoam->createControlDict(writeDicts);
85  fvSchemes_ = initFoam->createFvSchemes(writeDicts);
86  fvSolution_ = initFoam->createFvSolution(writeDicts);
87 
88  // creates dictionaries in constant folder
89  createBlockMshDict(writeDicts);
90 
91  runTime_ =
92  std::unique_ptr<Foam::Time>(new Foam::Time(controlDict_.get(), ".", "."));
93 
94  //- 2d cartesian mesher cannot be run in parallel
95  Foam::argList::noParallel();
96 }
std::unique_ptr< Foam::Time > runTime_
Definition: blockMeshGen.H:101
std::unique_ptr< Foam::dictionary > fvSchemes_
Definition: blockMeshGen.H:105
blockMeshParams * params_
Definition of parameters pointer.
Definition: blockMeshGen.H:97
bool isPackMesh
Boolean for if the operation is packmesh.
std::unique_ptr< Foam::dictionary > fvSolution_
Definition: blockMeshGen.H:106
void createBlockMshDict(const bool &write)
Creates blockMeshDict from user arguments.
Definition: blockMeshGen.C:124
std::unique_ptr< Foam::dictionary > controlDict_
Definition: blockMeshGen.H:104

Member Data Documentation

◆ blockMshDict_

std::unique_ptr<Foam::dictionary> blockMeshGen::blockMshDict_
private

Definition at line 103 of file blockMeshGen.H.

Referenced by createBlockMshDict(), and createMeshFromSTL().

◆ controlDict_

std::unique_ptr<Foam::dictionary> blockMeshGen::controlDict_
private

Definition at line 104 of file blockMeshGen.H.

Referenced by initialize().

◆ dataSet

vtkSmartPointer<vtkDataSet> meshGen::dataSet
protectedinherited

Definition at line 57 of file meshGen.H.

◆ defaults

bool blockMeshGen::defaults
private

Definition at line 93 of file blockMeshGen.H.

Referenced by blockMeshGen().

◆ fmesh_

std::unique_ptr<Foam::fvMesh> blockMeshGen::fmesh_
private

Definition at line 102 of file blockMeshGen.H.

Referenced by createMeshFromSTL().

◆ fvSchemes_

std::unique_ptr<Foam::dictionary> blockMeshGen::fvSchemes_
private

Definition at line 105 of file blockMeshGen.H.

Referenced by initialize().

◆ fvSolution_

std::unique_ptr<Foam::dictionary> blockMeshGen::fvSolution_
private

Definition at line 106 of file blockMeshGen.H.

Referenced by initialize().

◆ gmData

◆ params_

blockMeshParams* blockMeshGen::params_
private

Definition at line 97 of file blockMeshGen.H.

Referenced by blockMeshGen(), createBlockMshDict(), createMeshFromSTL(), and initialize().

◆ runTime_

std::unique_ptr<Foam::Time> blockMeshGen::runTime_
private

Definition at line 101 of file blockMeshGen.H.

Referenced by createMeshFromSTL(), and initialize().


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