conformalVoronoiMeshIO.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2012-2015 OpenFOAM Foundation
6  \\/ M anipulation | Copyright (C) 2015 OpenCFD Ltd.
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "conformalVoronoiMesh.H"
27 #include "IOstreams.H"
28 #include "OFstream.H"
29 #include "pointMesh.H"
30 #include "pointFields.H"
31 #include "ListOps.H"
32 #include "polyMeshFilter.H"
33 #include "polyTopoChange.H"
34 #include "PrintTable.H"
35 #include "pointMesh.H"
36 #include "indexedVertexOps.H"
37 #include "DelaunayMeshTools.H"
38 #include "syncTools.H"
39 #include "faceSet.H"
40 #include "OBJstream.H"
41 
42 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
43 
45 (
46  const string& description
47 ) const
48 {
49  timeCheck(time(), description, foamyHexMeshControls().timeChecks());
50 }
51 
52 
54 (
55  const Time& runTime,
56  const string& description,
57  const bool check
58 )
59 {
60  if (check)
61  {
62  Info<< nl << "--- [ cpuTime "
63  << runTime.elapsedCpuTime() << " s, "
64  << "delta " << runTime.cpuTimeIncrement()<< " s";
65 
66  if (description != word::null)
67  {
68  Info<< ", " << description << " ";
69  }
70  else
71  {
72  Info<< " ";
73  }
74 
75  Info<< "] --- " << endl;
76 
77  memInfo m;
78 
79  if (m.valid())
80  {
81  PrintTable<word, label> memoryTable
82  (
83  "Memory Usage (kB): "
84  + description
85  );
86 
87  memoryTable.add("mSize", m.size());
88  memoryTable.add("mPeak", m.peak());
89  memoryTable.add("mRss", m.rss());
90 
91  Info<< incrIndent;
92  memoryTable.print(Info, true, true);
93  Info<< decrIndent;
94  }
95  }
96 }
97 
98 
99 void Foam::conformalVoronoiMesh::writeMesh(const fileName& instance)
100 {
102 
103  // Per cell the Delaunay vertex
104  labelList cellToDelaunayVertex;
105  // Per patch, per face the Delaunay vertex
106  labelListList patchToDelaunayVertex;
107 
108  labelList dualPatchStarts;
109 
110  {
112  labelList boundaryPts;
113  faceList faces;
114  labelList owner;
115  labelList neighbour;
117  PtrList<dictionary> patchDicts;
118  pointField cellCentres;
119  PackedBoolList boundaryFacesToRemove;
120 
121  calcDualMesh
122  (
123  points,
124  boundaryPts,
125  faces,
126  owner,
127  neighbour,
128  patchNames,
129  patchDicts,
130  cellCentres,
131  cellToDelaunayVertex,
132  patchToDelaunayVertex,
133  boundaryFacesToRemove
134  );
135 
136  Info<< nl << "Writing polyMesh to " << instance << endl;
137 
138  writeMesh
139  (
141  instance,
142  points,
143  boundaryPts,
144  faces,
145  owner,
146  neighbour,
147  patchNames,
148  patchDicts,
149  cellCentres,
150  boundaryFacesToRemove
151  );
152 
153  dualPatchStarts.setSize(patchDicts.size());
154 
155  forAll(dualPatchStarts, patchI)
156  {
157  dualPatchStarts[patchI] =
158  readLabel(patchDicts[patchI].lookup("startFace"));
159  }
160  }
161 
162  if (foamyHexMeshControls().writeCellShapeControlMesh())
163  {
164  cellShapeControls().shapeControlMesh().write();
165  }
166 
167  if (foamyHexMeshControls().writeBackgroundMeshDecomposition())
168  {
169  Info<< nl << "Writing " << "backgroundMeshDecomposition" << endl;
170 
171  // Have to explicitly update the mesh instance.
172  const_cast<fvMesh&>(decomposition_().mesh()).setInstance
173  (
174  time().timeName()
175  );
176 
177  decomposition_().mesh().write();
178  }
179 
180  if (foamyHexMeshControls().writeTetDualMesh())
181  {
182  label cellI = 0;
183  for
184  (
185  Finite_cells_iterator cit = finite_cells_begin();
186  cit != finite_cells_end();
187  ++cit
188  )
189  {
190  if
191  (
192  !cit->hasFarPoint()
193  && !is_infinite(cit)
194  )
195  {
196  cit->cellIndex() = cellI++;
197  }
198  }
199 
200  Info<< nl << "Writing " << "tetDualMesh" << endl;
201 
203  labelList cellMap;
204  autoPtr<polyMesh> tetMesh =
205  createMesh("tetDualMesh", vertexMap, cellMap);
206 
207  tetMesh().write();
208 
209 // // Determine map from Delaunay vertex to Dual mesh
210 // // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
211 //
212 // // From all Delaunay vertices to cell (positive index)
213 // // or patch face (negative index)
214 // labelList vertexToDualAddressing(number_of_vertices(), 0);
215 //
216 // forAll(cellToDelaunayVertex, cellI)
217 // {
218 // label vertI = cellToDelaunayVertex[cellI];
219 //
220 // if (vertexToDualAddressing[vertI] != 0)
221 // {
222 // FatalErrorInFunction
223 // << "Delaunay vertex " << vertI
224 // << " from cell " << cellI
225 // << " is already mapped to "
226 // << vertexToDualAddressing[vertI]
227 // << exit(FatalError);
228 // }
229 // vertexToDualAddressing[vertI] = cellI+1;
230 // }
231 //
232 // forAll(patchToDelaunayVertex, patchI)
233 // {
234 // const labelList& patchVertices = patchToDelaunayVertex[patchI];
235 //
236 // forAll(patchVertices, i)
237 // {
238 // label vertI = patchVertices[i];
239 //
240 // if (vertexToDualAddressing[vertI] > 0)
241 // {
242 // FatalErrorInFunction
243 // << "Delaunay vertex " << vertI
244 // << " from patch " << patchI
245 // << " local index " << i
246 // << " is already mapped to cell "
247 // << vertexToDualAddressing[vertI]-1
248 // << exit(FatalError);
249 // }
250 //
251 // // Vertex might be used by multiple faces. Which one to
252 // // use? For now last one wins.
253 // label dualFaceI = dualPatchStarts[patchI]+i;
254 // vertexToDualAddressing[vertI] = -dualFaceI-1;
255 // }
256 // }
257 //
258 //
259 // // Calculate tet mesh addressing
260 // // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
261 //
262 // pointField points;
263 // labelList boundaryPts(number_of_finite_cells(), -1);
264 // // From tet point back to Delaunay vertex index
265 // labelList pointToDelaunayVertex;
266 // faceList faces;
267 // labelList owner;
268 // labelList neighbour;
269 // wordList patchTypes;
270 // wordList patchNames;
271 // PtrList<dictionary> patchDicts;
272 // pointField cellCentres;
273 //
274 // calcTetMesh
275 // (
276 // points,
277 // pointToDelaunayVertex,
278 // faces,
279 // owner,
280 // neighbour,
281 // patchTypes,
282 // patchNames,
283 // patchDicts
284 // );
285 //
286 //
287 //
288 // // Calculate map from tet points to dual mesh cells/patch faces
289 // // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
290 //
291 // labelIOList pointDualAddressing
292 // (
293 // IOobject
294 // (
295 // "pointDualAddressing",
296 // instance,
297 // "tetDualMesh"/polyMesh::meshSubDir,
298 // runTime_,
299 // IOobject::NO_READ,
300 // IOobject::AUTO_WRITE,
301 // false
302 // ),
303 // UIndirectList<label>
304 // (
305 // vertexToDualAddressing,
306 // pointToDelaunayVertex
307 // )()
308 // );
309 //
310 // label pointI = findIndex(pointDualAddressing, -1);
311 // if (pointI != -1)
312 // {
313 // WarningInFunction
314 // << "Delaunay vertex " << pointI
315 // << " does not have a corresponding dual cell." << endl;
316 // }
317 //
318 // Info<< "Writing map from tetDualMesh points to Voronoi mesh to "
319 // << pointDualAddressing.objectPath() << endl;
320 // pointDualAddressing.write();
321 //
322 //
323 //
324 // // Write tet points corresponding to the Voronoi cell/face centre
325 // // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
326 // {
327 // // Read Voronoi mesh
328 // fvMesh mesh
329 // (
330 // IOobject
331 // (
332 // Foam::polyMesh::defaultRegion,
333 // instance,
334 // runTime_,
335 // IOobject::MUST_READ
336 // )
337 // );
338 // pointIOField dualPoints
339 // (
340 // IOobject
341 // (
342 // "dualPoints",
343 // instance,
344 // "tetDualMesh"/polyMesh::meshSubDir,
345 // runTime_,
346 // IOobject::NO_READ,
347 // IOobject::AUTO_WRITE,
348 // false
349 // ),
350 // points
351 // );
352 //
353 // forAll(pointDualAddressing, pointI)
354 // {
355 // label index = pointDualAddressing[pointI];
356 //
357 // if (index > 0)
358 // {
359 // label cellI = index-1;
360 // dualPoints[pointI] = mesh.cellCentres()[cellI];
361 // }
362 // else if (index < 0)
363 // {
364 // label faceI = -index-1;
365 // if (faceI >= mesh.nInternalFaces())
366 // {
367 // dualPoints[pointI] = mesh.faceCentres()[faceI];
368 // }
369 // }
370 // }
371 //
372 // Info<< "Writing tetDualMesh points mapped onto Voronoi mesh to "
373 // << dualPoints.objectPath() << endl
374 // << "Replace the polyMesh/points with these." << endl;
375 // dualPoints.write();
376 // }
377 //
378 //
379 // Info<< nl << "Writing tetDualMesh to " << instance << endl;
380 //
381 // PackedBoolList boundaryFacesToRemove;
382 // writeMesh
383 // (
384 // "tetDualMesh",
385 // instance,
386 // points,
387 // boundaryPts,
388 // faces,
389 // owner,
390 // neighbour,
391 // patchTypes,
392 // patchNames,
393 // patchDicts,
394 // cellCentres,
395 // boundaryFacesToRemove
396 // );
397  }
398 }
399 
400 
402 (
403  const IOobject& io,
404  const wordList& patchNames,
405  const PtrList<dictionary>& patchDicts
406 ) const
407 {
408  autoPtr<fvMesh> meshPtr
409  (
410  new fvMesh
411  (
412  io,
413  xferCopy(pointField()),
414  xferCopy(faceList()),
415  xferCopy(cellList())
416  )
417  );
418  fvMesh& mesh = meshPtr();
419 
421 
422  forAll(patches, patchI)
423  {
424  if
425  (
426  patchDicts.set(patchI)
427  && (
428  word(patchDicts[patchI].lookup("type"))
429  == processorPolyPatch::typeName
430  )
431  )
432  {
433  patches[patchI] = new processorPolyPatch
434  (
435  patchNames[patchI],
436  0, //patchSizes[p],
437  0, //patchStarts[p],
438  patchI,
439  mesh.boundaryMesh(),
440  readLabel(patchDicts[patchI].lookup("myProcNo")),
441  readLabel(patchDicts[patchI].lookup("neighbProcNo")),
443  );
444  }
445  else
446  {
447  patches[patchI] = polyPatch::New
448  (
449  patchDicts[patchI].lookup("type"),
450  patchNames[patchI],
451  0, //patchSizes[p],
452  0, //patchStarts[p],
453  patchI,
455  ).ptr();
456  }
457  }
458 
460 
461  return meshPtr;
462 }
463 
464 
466 (
467  const PtrList<dictionary>& patchDicts
468 ) const
469 {
470  // Check patch sizes
471  labelListList procPatchSizes
472  (
473  Pstream::nProcs(),
475  );
476 
477  forAll(patchDicts, patchI)
478  {
479  if
480  (
481  patchDicts.set(patchI)
482  && (
483  word(patchDicts[patchI].lookup("type"))
484  == processorPolyPatch::typeName
485  )
486  )
487  {
488  const label procNeighb =
489  readLabel(patchDicts[patchI].lookup("neighbProcNo"));
490 
491  procPatchSizes[Pstream::myProcNo()][procNeighb]
492  = readLabel(patchDicts[patchI].lookup("nFaces"));
493  }
494  }
495 
496  Pstream::gatherList(procPatchSizes);
497 
498  if (Pstream::master())
499  {
500  bool allMatch = true;
501 
502  forAll(procPatchSizes, procI)
503  {
504  const labelList& patchSizes = procPatchSizes[procI];
505 
506  forAll(patchSizes, patchI)
507  {
508  if (patchSizes[patchI] != procPatchSizes[patchI][procI])
509  {
510  allMatch = false;
511 
512  Info<< indent << "Patches " << procI << " and " << patchI
513  << " have different sizes: " << patchSizes[patchI]
514  << " and " << procPatchSizes[patchI][procI] << endl;
515  }
516  }
517  }
518 
519  if (allMatch)
520  {
521  Info<< indent << "All processor patches have matching numbers of "
522  << "faces" << endl;
523  }
524  }
525 }
526 
527 
529 (
531  labelList& boundaryPts,
532  faceList& faces,
533  const label nInternalFaces
534 ) const
535 {
536  Info<< incrIndent << indent << "Reordering points into internal/external"
537  << endl;
538 
539  labelList oldToNew(points.size(), label(0));
540 
541  // Find points that are internal
542  for (label fI = nInternalFaces; fI < faces.size(); ++fI)
543  {
544  const face& f = faces[fI];
545 
546  forAll(f, fpI)
547  {
548  oldToNew[f[fpI]] = 1;
549  }
550  }
551 
552  const label nInternalPoints = points.size() - sum(oldToNew);
553 
554  label countInternal = 0;
555  label countExternal = nInternalPoints;
556 
557  forAll(points, pI)
558  {
559  if (oldToNew[pI] == 0)
560  {
561  oldToNew[pI] = countInternal++;
562  }
563  else
564  {
565  oldToNew[pI] = countExternal++;
566  }
567  }
568 
569  Info<< indent
570  << "Number of internal points: " << countInternal << nl
571  << indent << "Number of external points: " << countExternal
572  << decrIndent << endl;
573 
574  inplaceReorder(oldToNew, points);
575  inplaceReorder(oldToNew, boundaryPts);
576 
577  forAll(faces, fI)
578  {
579  face& f = faces[fI];
580 
581  forAll(f, fpI)
582  {
583  f[fpI] = oldToNew[f[fpI]];
584  }
585  }
586 }
587 
588 
590 (
591  const word& meshName,
592  const fileName& instance,
593  const pointField& points,
594  faceList& faces,
595  const wordList& patchNames,
596  const PtrList<dictionary>& patchDicts
597 ) const
598 {
599  Info<< incrIndent << indent << "Reordering processor patches" << endl;
600 
601  Info<< incrIndent;
602  checkProcessorPatchesMatch(patchDicts);
603 
604  // Create dummy mesh with correct proc boundaries to do sorting
605  autoPtr<fvMesh> sortMeshPtr
606  (
607  createDummyMesh
608  (
609  IOobject
610  (
611  meshName,
612  instance,
613  runTime_,
616  false
617  ),
618  patchNames,
619  patchDicts
620  )
621  );
622  const fvMesh& sortMesh = sortMeshPtr();
623 
624  // Change the transform type on processors to coincident full match.
625 // forAll(sortMesh.boundaryMesh(), patchI)
626 // {
627 // const polyPatch& patch = sortMesh.boundaryMesh()[patchI];
628 //
629 // if (isA<processorPolyPatch>(patch))
630 // {
631 // const processorPolyPatch& cpPatch
632 // = refCast<const processorPolyPatch>(patch);
633 //
634 // processorPolyPatch& pPatch
635 // = const_cast<processorPolyPatch&>(cpPatch);
636 //
637 // pPatch.transform() = coupledPolyPatch::COINCIDENTFULLMATCH;
638 // }
639 // }
640 
641  // Rotation on new faces.
642  labelList rotation(faces.size(), label(0));
643  labelList faceMap(faces.size(), label(-1));
644 
645  PstreamBuffers pBufs(Pstream::nonBlocking);
646 
647  // Send ordering
648  forAll(sortMesh.boundaryMesh(), patchI)
649  {
650  const polyPatch& pp = sortMesh.boundaryMesh()[patchI];
651 
652  if (isA<processorPolyPatch>(pp))
653  {
654  refCast<const processorPolyPatch>(pp).initOrder
655  (
656  pBufs,
658  (
659  SubList<face>
660  (
661  faces,
662  readLabel(patchDicts[patchI].lookup("nFaces")),
663  readLabel(patchDicts[patchI].lookup("startFace"))
664  ),
665  points
666  )
667  );
668  }
669  }
670 
671  pBufs.finishedSends();
672 
673  Info<< incrIndent << indent << "Face ordering initialised..." << endl;
674 
675  // Receive and calculate ordering
676  bool anyChanged = false;
677 
678  forAll(sortMesh.boundaryMesh(), patchI)
679  {
680  const polyPatch& pp = sortMesh.boundaryMesh()[patchI];
681 
682  if (isA<processorPolyPatch>(pp))
683  {
684  const label nPatchFaces =
685  readLabel(patchDicts[patchI].lookup("nFaces"));
686  const label patchStartFace =
687  readLabel(patchDicts[patchI].lookup("startFace"));
688 
689  labelList patchFaceMap(nPatchFaces, label(-1));
690  labelList patchFaceRotation(nPatchFaces, label(0));
691 
692  bool changed = refCast<const processorPolyPatch>(pp).order
693  (
694  pBufs,
696  (
697  SubList<face>
698  (
699  faces,
700  nPatchFaces,
701  patchStartFace
702  ),
703  points
704  ),
705  patchFaceMap,
706  patchFaceRotation
707  );
708 
709  if (changed)
710  {
711  // Merge patch face reordering into mesh face reordering table
712  forAll(patchFaceRotation, patchFaceI)
713  {
714  rotation[patchFaceI + patchStartFace]
715  = patchFaceRotation[patchFaceI];
716  }
717 
718  forAll(patchFaceMap, patchFaceI)
719  {
720  if (patchFaceMap[patchFaceI] != patchFaceI)
721  {
722  faceMap[patchFaceI + patchStartFace]
723  = patchFaceMap[patchFaceI] + patchStartFace;
724  }
725  }
726 
727  anyChanged = true;
728  }
729  }
730  }
731 
732  Info<< incrIndent << indent << "Faces matched." << endl;
733 
734  reduce(anyChanged, orOp<bool>());
735 
736  if (anyChanged)
737  {
738  label nReorderedFaces = 0;
739 
740  forAll(faceMap, faceI)
741  {
742  if (faceMap[faceI] != -1)
743  {
744  nReorderedFaces++;
745  }
746  }
747 
748  if (nReorderedFaces > 0)
749  {
750  inplaceReorder(faceMap, faces);
751  }
752 
753  // Rotate faces (rotation is already in new face indices).
754  label nRotated = 0;
755 
756  forAll(rotation, faceI)
757  {
758  if (rotation[faceI] != 0)
759  {
760  faces[faceI] = rotateList(faces[faceI], rotation[faceI]);
761  nRotated++;
762  }
763  }
764 
765  Info<< indent << returnReduce(nReorderedFaces, sumOp<label>())
766  << " faces have been reordered" << nl
767  << indent << returnReduce(nRotated, sumOp<label>())
768  << " faces have been rotated"
769  << decrIndent << decrIndent
770  << decrIndent << decrIndent << endl;
771  }
772 }
773 
774 
776 (
777  const word& meshName,
778  const fileName& instance,
780  labelList& boundaryPts,
781  faceList& faces,
782  labelList& owner,
783  labelList& neighbour,
784  const wordList& patchNames,
785  const PtrList<dictionary>& patchDicts,
786  const pointField& cellCentres,
787  PackedBoolList& boundaryFacesToRemove
788 ) const
789 {
790  if (foamyHexMeshControls().objOutput())
791  {
793  (
794  time().path()/word(meshName + ".obj"),
795  points,
796  faces
797  );
798  }
799 
800  const label nInternalFaces = readLabel(patchDicts[0].lookup("startFace"));
801 
802  reorderPoints(points, boundaryPts, faces, nInternalFaces);
803 
804  if (Pstream::parRun())
805  {
806  reorderProcessorPatches
807  (
808  meshName,
809  instance,
810  points,
811  faces,
812  patchNames,
813  patchDicts
814  );
815  }
816 
817  Info<< incrIndent;
818  Info<< indent << "Constructing mesh" << endl;
819 
820  timeCheck("Before fvMesh construction");
821 
822  fvMesh mesh
823  (
824  IOobject
825  (
826  meshName,
827  instance,
828  runTime_,
831  ),
832  xferMove(points),
833  xferMove(faces),
834  xferMove(owner),
835  xferMove(neighbour)
836  );
837 
838  Info<< indent << "Adding patches to mesh" << endl;
839 
841 
842  label nValidPatches = 0;
843 
844  forAll(patches, p)
845  {
846  label totalPatchSize = readLabel(patchDicts[p].lookup("nFaces"));
847 
848  if
849  (
850  patchDicts.set(p)
851  && (
852  word(patchDicts[p].lookup("type"))
853  == processorPolyPatch::typeName
854  )
855  )
856  {
857  const_cast<dictionary&>(patchDicts[p]).set
858  (
859  "transform",
860  "noOrdering"
861  //"coincidentFullMatch"
862  );
863 
864  // Do not create empty processor patches
865  if (totalPatchSize > 0)
866  {
867  patches[nValidPatches] = new processorPolyPatch
868  (
869  patchNames[p],
870  patchDicts[p],
871  nValidPatches,
872  mesh.boundaryMesh(),
873  processorPolyPatch::typeName
874  );
875 
876  nValidPatches++;
877  }
878  }
879  else
880  {
881  // Check that the patch is not empty on every processor
882  reduce(totalPatchSize, sumOp<label>());
883 
884  if (totalPatchSize > 0)
885  {
886  patches[nValidPatches] = polyPatch::New
887  (
888  patchNames[p],
889  patchDicts[p],
890  nValidPatches,
892  ).ptr();
893 
894  nValidPatches++;
895  }
896  }
897  }
898 
899  patches.setSize(nValidPatches);
900 
902 
903 
904 
905  // Add zones to the mesh
906  addZones(mesh, cellCentres);
907 
908 
909 
910  Info<< indent << "Add pointZones" << endl;
911  {
912  label sz = mesh.pointZones().size();
913 
914  DynamicList<label> bPts(boundaryPts.size());
915 
916  forAll(dualMeshPointTypeNames_, typeI)
917  {
918  forAll(boundaryPts, ptI)
919  {
920  const label& bPtType = boundaryPts[ptI];
921 
922  if (bPtType == typeI)
923  {
924  bPts.append(ptI);
925  }
926  }
927 
928 // syncTools::syncPointList(mesh, bPts, maxEqOp<label>(), -1);
929 
930  Info<< incrIndent << indent
931  << "Adding " << bPts.size()
932  << " points of type " << dualMeshPointTypeNames_.words()[typeI]
933  << decrIndent << endl;
934 
935  mesh.pointZones().append
936  (
937  new pointZone
938  (
939  dualMeshPointTypeNames_.words()[typeI],
940  bPts,
941  sz + typeI,
942  mesh.pointZones()
943  )
944  );
945 
946  bPts.clear();
947  }
948  }
949 
950 
951 
952  // Add indirectPatchFaces to a face zone
953  Info<< indent << "Adding indirect patch faces set" << endl;
954 
956  (
957  mesh,
958  boundaryFacesToRemove,
959  orEqOp<unsigned int>()
960  );
961 
962  labelList addr(boundaryFacesToRemove.count());
963  label count = 0;
964 
965  forAll(boundaryFacesToRemove, faceI)
966  {
967  if (boundaryFacesToRemove[faceI])
968  {
969  addr[count++] = faceI;
970  }
971  }
972 
973  addr.setSize(count);
974 
975  faceSet indirectPatchFaces
976  (
977  mesh,
978  "indirectPatchFaces",
979  addr,
981  );
982 
983  indirectPatchFaces.sync(mesh);
984 
985 
986  Info<< decrIndent;
987 
988  timeCheck("Before fvMesh filtering");
989 
990  autoPtr<polyMeshFilter> meshFilter;
991 
992  label nInitialBadFaces = 0;
993 
994  if (foamyHexMeshControls().filterEdges())
995  {
996  Info<< nl << "Filtering edges on polyMesh" << nl << endl;
997 
998  meshFilter.reset(new polyMeshFilter(mesh, boundaryPts));
999 
1000  // Filter small edges only. This reduces the number of faces so that
1001  // the face filtering is sped up.
1002  nInitialBadFaces = meshFilter().filterEdges(0);
1003  {
1004  const autoPtr<fvMesh>& newMesh = meshFilter().filteredMesh();
1005 
1006  polyTopoChange meshMod(newMesh());
1007 
1008  autoPtr<mapPolyMesh> map = meshMod.changeMesh(mesh, false);
1009 
1010  polyMeshFilter::copySets(newMesh(), mesh);
1011  }
1012  }
1013 
1014  if (foamyHexMeshControls().filterFaces())
1015  {
1016  labelIOList boundaryPtsIO
1017  (
1018  IOobject
1019  (
1020  "pointPriority",
1021  instance,
1022  time(),
1025  ),
1027  );
1028 
1029  forAll(mesh.points(), ptI)
1030  {
1031  boundaryPtsIO[ptI] = mesh.pointZones().whichZone(ptI);
1032  }
1033 
1034 
1035  Info<< nl << "Filtering faces on polyMesh" << nl << endl;
1036 
1037  meshFilter.reset(new polyMeshFilter(mesh, boundaryPtsIO));
1038 
1039  meshFilter().filter(nInitialBadFaces);
1040  {
1041  const autoPtr<fvMesh>& newMesh = meshFilter().filteredMesh();
1042 
1043  polyTopoChange meshMod(newMesh());
1044 
1045  autoPtr<mapPolyMesh> map = meshMod.changeMesh(mesh, false);
1046 
1047  polyMeshFilter::copySets(newMesh(), mesh);
1048  }
1049  }
1050 
1051  timeCheck("After fvMesh filtering");
1052 
1053  mesh.setInstance(instance);
1054 
1055  if (!mesh.write())
1056  {
1058  << "Failed writing polyMesh."
1059  << exit(FatalError);
1060  }
1061  else
1062  {
1063  Info<< nl << "Written filtered mesh to "
1064  << mesh.polyMesh::instance() << nl
1065  << endl;
1066  }
1067 
1068  {
1069  pointScalarField boundaryPtsScalarField
1070  (
1071  IOobject
1072  (
1073  "boundaryPoints_collapsed",
1074  instance,
1075  time(),
1078  ),
1080  dimensionedScalar("min", dimless, scalar(labelMin))
1081  );
1082 
1083  labelIOList boundaryPtsIO
1084  (
1085  IOobject
1086  (
1087  "pointPriority",
1088  instance,
1089  time(),
1092  ),
1094  );
1095 
1096  forAll(mesh.points(), ptI)
1097  {
1098  boundaryPtsScalarField[ptI] = mesh.pointZones().whichZone(ptI);
1099  boundaryPtsIO[ptI] = mesh.pointZones().whichZone(ptI);
1100  }
1101 
1102  boundaryPtsScalarField.write();
1103  boundaryPtsIO.write();
1104  }
1105 
1106 // writeCellSizes(mesh);
1107 
1108 // writeCellAlignments(mesh);
1109 
1110 // writeCellCentres(mesh);
1111 
1112  findRemainingProtrusionSet(mesh);
1113 }
1114 
1115 
1117 (
1118  const fvMesh& mesh
1119 ) const
1120 {
1121  {
1122  timeCheck("Start writeCellSizes");
1123 
1124  Info<< nl << "Create targetCellSize volScalarField" << endl;
1125 
1126  volScalarField targetCellSize
1127  (
1128  IOobject
1129  (
1130  "targetCellSize",
1131  mesh.polyMesh::instance(),
1132  mesh,
1135  ),
1136  mesh,
1137  dimensionedScalar("cellSize", dimLength, 0),
1138  zeroGradientFvPatchScalarField::typeName
1139  );
1140 
1141  scalarField& cellSize = targetCellSize.internalField();
1142 
1143  const vectorField& C = mesh.cellCentres();
1144 
1145  forAll(cellSize, i)
1146  {
1147  cellSize[i] = cellShapeControls().cellSize(C[i]);
1148  }
1149 
1150  // Info<< nl << "Create targetCellVolume volScalarField" << endl;
1151 
1152  // volScalarField targetCellVolume
1153  // (
1154  // IOobject
1155  // (
1156  // "targetCellVolume",
1157  // mesh.polyMesh::instance(),
1158  // mesh,
1159  // IOobject::NO_READ,
1160  // IOobject::AUTO_WRITE
1161  // ),
1162  // mesh,
1163  // dimensionedScalar("cellVolume", dimLength, 0),
1164  // zeroGradientFvPatchScalarField::typeName
1165  // );
1166 
1167  // targetCellVolume.internalField() = pow3(cellSize);
1168 
1169  // Info<< nl << "Create actualCellVolume volScalarField" << endl;
1170 
1171  // volScalarField actualCellVolume
1172  // (
1173  // IOobject
1174  // (
1175  // "actualCellVolume",
1176  // mesh.polyMesh::instance(),
1177  // mesh,
1178  // IOobject::NO_READ,
1179  // IOobject::AUTO_WRITE
1180  // ),
1181  // mesh,
1182  // dimensionedScalar("cellVolume", dimVolume, 0),
1183  // zeroGradientFvPatchScalarField::typeName
1184  // );
1185 
1186  // actualCellVolume.internalField() = mesh.cellVolumes();
1187 
1188  // Info<< nl << "Create equivalentCellSize volScalarField" << endl;
1189 
1190  // volScalarField equivalentCellSize
1191  // (
1192  // IOobject
1193  // (
1194  // "equivalentCellSize",
1195  // mesh.polyMesh::instance(),
1196  // mesh,
1197  // IOobject::NO_READ,
1198  // IOobject::AUTO_WRITE
1199  // ),
1200  // mesh,
1201  // dimensionedScalar("cellSize", dimLength, 0),
1202  // zeroGradientFvPatchScalarField::typeName
1203  // );
1204 
1205  // equivalentCellSize.internalField() = pow
1206  // (
1207  // actualCellVolume.internalField(),
1208  // 1.0/3.0
1209  // );
1210 
1211  targetCellSize.correctBoundaryConditions();
1212  // targetCellVolume.correctBoundaryConditions();
1213  // actualCellVolume.correctBoundaryConditions();
1214  // equivalentCellSize.correctBoundaryConditions();
1215 
1216  targetCellSize.write();
1217  // targetCellVolume.write();
1218  // actualCellVolume.write();
1219  // equivalentCellSize.write();
1220  }
1221 
1222  // {
1223  // polyMesh tetMesh
1224  // (
1225  // IOobject
1226  // (
1227  // "tetDualMesh",
1228  // runTime_.constant(),
1229  // runTime_,
1230  // IOobject::MUST_READ
1231  // )
1232  // );
1233 
1234  // pointMesh ptMesh(tetMesh);
1235 
1236  // pointScalarField ptTargetCellSize
1237  // (
1238  // IOobject
1239  // (
1240  // "ptTargetCellSize",
1241  // runTime_.timeName(),
1242  // tetMesh,
1243  // IOobject::NO_READ,
1244  // IOobject::AUTO_WRITE
1245  // ),
1246  // ptMesh,
1247  // dimensionedScalar("ptTargetCellSize", dimLength, 0),
1248  // pointPatchVectorField::calculatedType()
1249  // );
1250 
1251  // scalarField& cellSize = ptTargetCellSize.internalField();
1252 
1253  // const vectorField& P = tetMesh.points();
1254 
1255  // forAll(cellSize, i)
1256  // {
1257  // cellSize[i] = cellShapeControls().cellSize(P[i]);
1258  // }
1259 
1260  // ptTargetCellSize.write();
1261  // }
1262 }
1263 
1264 
1266 (
1267  const fvMesh& mesh
1268 ) const
1269 {
1270 // Info<< nl << "Create cellAlignments volTensorField" << endl;
1271 //
1272 // volTensorField cellAlignments
1273 // (
1274 // IOobject
1275 // (
1276 // "cellAlignments",
1277 // mesh.polyMesh::instance(),
1278 // mesh,
1279 // IOobject::NO_READ,
1280 // IOobject::AUTO_WRITE
1281 // ),
1282 // mesh,
1283 // tensor::I,
1284 // zeroGradientFvPatchTensorField::typeName
1285 // );
1286 //
1287 // tensorField& cellAlignment = cellAlignments.internalField();
1288 //
1289 // const vectorField& C = mesh.cellCentres();
1290 //
1291 // vectorField xDir(cellAlignment.size());
1292 // vectorField yDir(cellAlignment.size());
1293 // vectorField zDir(cellAlignment.size());
1294 //
1295 // forAll(cellAlignment, i)
1296 // {
1297 // cellAlignment[i] = cellShapeControls().cellAlignment(C[i]);
1298 // xDir[i] = cellAlignment[i] & vector(1, 0, 0);
1299 // yDir[i] = cellAlignment[i] & vector(0, 1, 0);
1300 // zDir[i] = cellAlignment[i] & vector(0, 0, 1);
1301 // }
1302 //
1303 // OFstream xStr("xDir.obj");
1304 // OFstream yStr("yDir.obj");
1305 // OFstream zStr("zDir.obj");
1306 //
1307 // forAll(xDir, i)
1308 // {
1309 // meshTools::writeOBJ(xStr, C[i], C[i] + xDir[i]);
1310 // meshTools::writeOBJ(yStr, C[i], C[i] + yDir[i]);
1311 // meshTools::writeOBJ(zStr, C[i], C[i] + zDir[i]);
1312 // }
1313 //
1314 // cellAlignments.correctBoundaryConditions();
1315 //
1316 // cellAlignments.write();
1317 }
1318 
1319 
1321 (
1322  const fvMesh& mesh
1323 ) const
1324 {
1325  Info<< "Writing components of cellCentre positions to volScalarFields"
1326  << " ccx, ccy, ccz in " << runTime_.timeName() << endl;
1327 
1328  for (direction i=0; i<vector::nComponents; i++)
1329  {
1330  volScalarField cci
1331  (
1332  IOobject
1333  (
1334  "cc" + word(vector::componentNames[i]),
1335  runTime_.timeName(),
1336  mesh,
1339  ),
1340  mesh.C().component(i)
1341  );
1342 
1343  cci.write();
1344  }
1345 }
1346 
1347 
1349 (
1350  const polyMesh& mesh
1351 ) const
1352 {
1353  timeCheck("Start findRemainingProtrusionSet");
1354 
1355  const polyBoundaryMesh& patches = mesh.boundaryMesh();
1356 
1357  labelHashSet protrudingBoundaryPoints;
1358 
1359  forAll(patches, patchI)
1360  {
1361  const polyPatch& patch = patches[patchI];
1362 
1363  forAll(patch.localPoints(), pLPI)
1364  {
1365  label meshPtI = patch.meshPoints()[pLPI];
1366 
1367  const Foam::point& pt = patch.localPoints()[pLPI];
1368 
1369  if
1370  (
1371  geometryToConformTo_.wellOutside
1372  (
1373  pt,
1374  sqr(targetCellSize(pt))
1375  )
1376  )
1377  {
1378  protrudingBoundaryPoints.insert(meshPtI);
1379  }
1380  }
1381  }
1382 
1383  cellSet protrudingCells
1384  (
1385  mesh,
1386  "foamyHexMesh_remainingProtrusions",
1387  mesh.nCells()/1000
1388  );
1389 
1390  forAllConstIter(labelHashSet, protrudingBoundaryPoints, iter)
1391  {
1392  const label pointI = iter.key();
1393  const labelList& pCells = mesh.pointCells()[pointI];
1394 
1395  forAll(pCells, pCI)
1396  {
1397  protrudingCells.insert(pCells[pCI]);
1398  }
1399  }
1400 
1401  label protrudingCellsSize = protrudingCells.size();
1402 
1403  reduce(protrudingCellsSize, sumOp<label>());
1404 
1405  if (foamyHexMeshControls().objOutput() && protrudingCellsSize > 0)
1406  {
1407  Info<< nl << "Found " << protrudingCellsSize
1408  << " cells protruding from the surface, writing cellSet "
1409  << protrudingCells.name()
1410  << endl;
1411 
1412  protrudingCells.write();
1413  }
1414 
1415  return protrudingCells;
1416 }
1417 
1418 
1420 (
1421  const fileName& fName
1422 ) const
1423 {
1424  OBJstream os(fName);
1425 
1426  for
1427  (
1428  Delaunay::Finite_edges_iterator eit = finite_edges_begin();
1429  eit != finite_edges_end();
1430  ++eit
1431  )
1432  {
1433  Cell_handle c = eit->first;
1434  Vertex_handle vA = c->vertex(eit->second);
1435  Vertex_handle vB = c->vertex(eit->third);
1436 
1437  if (ptPairs_.isPointPair(vA, vB))
1438  {
1439  os.write
1440  (
1441  linePointRef(topoint(vA->point()), topoint(vB->point()))
1442  );
1443  }
1444  }
1445 }
1446 
1447 
1448 // ************************************************************************* //
Foam::pointField
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
Foam::ZoneMesh::clear
void clear()
Clear the zones.
Definition: ZoneMesh.C:408
Foam::polyMesh::points
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:979
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:48
Foam::GeometricField::component
tmp< GeometricField< cmptType, PatchField, GeoMesh > > component(const direction) const
Return a component of the field.
PrintTable.H
Foam::dimless
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
IOstreams.H
Useful combination of include files which define Sin, Sout and Serr and the use of IO streams general...
Foam::faceMap
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
Definition: blockMeshMergeFast.C:90
p
p
Definition: pEqn.H:62
Foam::labelList
List< label > labelList
A List of labels.
Definition: labelList.H:56
Foam::returnReduce
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Definition: PstreamReduceOps.H:86
Foam::conformalVoronoiMesh::writeCellCentres
void writeCellCentres(const fvMesh &mesh) const
Calculate and write the cell centres.
Foam::IOobject::AUTO_WRITE
@ AUTO_WRITE
Definition: IOobject.H:117
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
Foam::polyMesh::defaultRegion
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:306
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::fvMesh::addFvPatches
void addFvPatches(const List< polyPatch * > &, const bool validBoundary=true)
Add boundary patches. Constructor helper.
Definition: fvMesh.C:462
Foam::conformalVoronoiMesh::writeMesh
void writeMesh(const fileName &instance)
Prepare data and call writeMesh for polyMesh and.
Foam::UPstream::nProcs
static label nProcs(const label communicator=0)
Number of processes in parallel run.
Definition: UPstream.H:387
Foam::conformalVoronoiMesh::findRemainingProtrusionSet
labelHashSet findRemainingProtrusionSet(const polyMesh &mesh) const
Find the cellSet of the boundary cells which have points that.
Foam::conformalVoronoiMesh::timeCheck
static void timeCheck(const Time &runTime, const string &description=string::null, const bool check=true)
Write the elapsedCpuTime and memory usage, with an optional.
polyTopoChange.H
Foam::UPstream::parRun
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:377
Foam::VectorSpace< Vector< scalar >, scalar, 3 >::nComponents
@ nComponents
Number of components in this vector space.
Definition: VectorSpace.H:88
Foam::primitivePatch
PrimitivePatch< face, SubList, const pointField & > primitivePatch
Foam::primitivePatch.
Definition: primitivePatch.H:45
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:421
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::conformalVoronoiMesh::createDummyMesh
autoPtr< fvMesh > createDummyMesh(const IOobject &io, const wordList &patchNames, const PtrList< dictionary > &patchDicts) const
Create an empty fvMesh.
Foam::Vector< scalar >::componentNames
static const char * componentNames[]
Definition: Vector.H:79
Foam::HashSet< label, Hash< label > >
Foam::labelIOList
IOList< label > labelIOList
Label container classes.
Definition: labelIOList.H:42
Foam::incrIndent
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition: Ostream.H:228
syncTools.H
Foam::IOstream::print
virtual void print(Ostream &) const
Print description of IOstream to Ostream.
Definition: IOstream.C:126
Foam::topoint
pointFromPoint topoint(const Point &P)
Definition: pointConversion.H:70
Foam::IOobject::NO_WRITE
@ NO_WRITE
Definition: IOobject.H:118
Foam::MeshObject< polyMesh, UpdateableMeshObject, pointMesh >::New
static const pointMesh & New(const polyMesh &mesh)
Definition: MeshObject.C:44
C
volScalarField & C
Definition: readThermalProperties.H:104
patchDicts
PtrList< dictionary > patchDicts
Definition: readKivaGrid.H:537
forAllConstIter
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:39
OFstream.H
Foam::vectorField
Field< vector > vectorField
Specialisation of Field<T> for vector.
Definition: primitiveFieldsFwd.H:49
Foam::wordList
List< word > wordList
A List of words.
Definition: fileName.H:54
Foam::primitiveMesh::nCells
label nCells() const
Definition: primitiveMeshI.H:64
Foam::fvMesh::write
virtual bool write() const
Write mesh using IO settings from time.
Definition: fvMesh.C:873
Foam::conformalVoronoiMesh::writeCellSizes
void writeCellSizes(const fvMesh &mesh) const
Calculate and write a field of the target cell size,.
Foam::primitiveMesh::nPoints
label nPoints() const
Definition: primitiveMeshI.H:35
Foam::reduce
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
Definition: PstreamReduceOps.H:43
Foam::labelMin
static const label labelMin
Definition: label.H:61
Foam::UPstream::nonBlocking
@ nonBlocking
Definition: UPstream.H:68
Foam::inplaceReorder
void inplaceReorder(const labelUList &oldToNew, ListType &)
Inplace reorder the elements of a list.
Definition: ListOpsTemplates.C:102
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Foam::DelaunayMeshTools::writeObjMesh
void writeObjMesh(const fileName &fName, const pointField &points, const faceList &faces)
Write an OBJ mesh consisting of points and faces.
Foam::Ostream::write
virtual Ostream & write(const token &)=0
Write next token to stream.
Foam::IOobject::NO_READ
@ NO_READ
Definition: IOobject.H:111
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
Foam::polyMesh::pointZones
const pointZoneMesh & pointZones() const
Return point zone mesh.
Definition: polyMesh.H:457
faceSet.H
Foam::cellList
List< cell > cellList
list of cells
Definition: cellList.H:42
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:41
Foam::coupledPolyPatch::COINCIDENTFULLMATCH
@ COINCIDENTFULLMATCH
Definition: coupledPolyPatch.H:62
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
patchNames
wordList patchNames(nPatches)
Foam::ZoneMesh::whichZone
label whichZone(const label objectIndex) const
Given a global object index, return the zone it is in.
Definition: ZoneMesh.C:226
Foam::conformalVoronoiMesh::writePointPairs
void writePointPairs(const fileName &fName) const
meshPtr
static fvMesh * meshPtr
Definition: globalFoam.H:52
Foam::conformalVoronoiMesh::reorderPoints
void reorderPoints(pointField &points, labelList &boundaryPts, faceList &faces, const label nInternalFaces) const
Foam::FatalError
error FatalError
Foam::polyPatch::New
static autoPtr< polyPatch > New(const word &patchType, const word &name, const label size, const label start, const label index, const polyBoundaryMesh &bm)
Return a pointer to a new patch created on freestore from.
Definition: polyPatchNew.C:32
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
polyMeshFilter.H
Foam::rotateList
ListType rotateList(const ListType &list, const label n)
Rotate a list by n places. If n is positive rotate clockwise/right/down.
Definition: ListOpsTemplates.C:803
DelaunayMeshTools.H
Foam::decrIndent
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition: Ostream.H:235
Foam::indent
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:221
Foam::pointScalarField
GeometricField< scalar, pointPatchField, pointMesh > pointScalarField
Definition: pointFields.H:49
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Foam::UPstream::myProcNo
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:405
Foam::UPstream::master
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:399
Foam::polyMesh::setInstance
void setInstance(const fileName &)
Set the instance for mesh files.
Definition: polyMeshIO.C:32
Foam::autoPtr< Foam::fvMesh >
Foam::conformalVoronoiMesh::foamyHexMeshControls
const cvControls & foamyHexMeshControls() const
Return the foamyHexMeshControls object.
Definition: conformalVoronoiMeshI.H:578
Foam::labelListList
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:57
Foam::conformalVoronoiMesh::time
const Time & time() const
Return the Time object.
Definition: conformalVoronoiMeshI.H:530
Foam::primitiveMesh::cellCentres
const vectorField & cellCentres() const
Definition: primitiveMeshCellCentresAndVols.C:211
Foam::linePointRef
line< point, const point & > linePointRef
Line using referred points.
Definition: linePointRef.H:45
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
Foam::Pstream::gatherList
static void gatherList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Gather data but keep individual values separate.
Definition: gatherScatterList.C:49
Foam::DelaunayMesh< Delaunay >::labelTolabelPairHashTable
HashTable< label, labelPair, FixedList< label, 2 >::Hash<> > labelTolabelPairHashTable
Definition: DelaunayMesh.H:90
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:49
Foam::fvMesh::C
const volVectorField & C() const
Return cell centres as volVectorField.
Definition: fvMeshGeometry.C:369
Foam::xferMove
Xfer< T > xferMove(T &)
Construct by transferring the contents of the arg.
f
labelList f(nPoints)
Foam::faceList
List< face > faceList
Definition: faceListFwd.H:43
Foam::xferCopy
Xfer< T > xferCopy(const T &)
Construct by copying the contents of the arg.
Foam::Vector< scalar >
Foam::conformalVoronoiMesh::writeCellAlignments
void writeCellAlignments(const fvMesh &mesh) const
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::primitiveMesh::pointCells
const labelListList & pointCells() const
Definition: primitiveMeshPointCells.C:108
Foam::sum
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
Definition: DimensionedFieldFunctions.C:333
path
fileName path(UMean.rootPath()/UMean.caseName()/"graphs"/UMean.instance())
Foam::polyMeshFilter::copySets
static void copySets(const polyMesh &oldMesh, const polyMesh &newMesh)
Definition: polyMeshFilterTemplates.C:71
Foam::word::null
static const word null
An empty word.
Definition: word.H:77
Foam::readLabel
label readLabel(Istream &is)
Definition: label.H:64
Foam::direction
unsigned char direction
Definition: direction.H:43
meshName
static char meshName[]
Definition: globalFoam.H:7
Foam::constant::universal::c
const dimensionedScalar c
Speed of light in a vacuum.
patches
patches[0]
Definition: createSingleCellMesh.H:36
List
Definition: Test.C:19
ListOps.H
Various functions to operate on Lists.
Foam::syncTools::syncFaceList
static void syncFaceList(const polyMesh &mesh, UList< T > &l, const CombineOp &cop)
Synchronize values on all mesh faces.
Definition: syncTools.H:381
Foam::List::size
void size(const label)
Override size to be inconsistent with allocated storage.
Foam::labelHashSet
HashSet< label, Hash< label > > labelHashSet
A HashSet with label keys.
Definition: HashSet.H:210
Foam::DelaunayMeshTools::writeInternalDelaunayVertices
void writeInternalDelaunayVertices(const fileName &instance, const Triangulation &t)
Write the internal Delaunay vertices of the tessellation as a.
Foam::conformalVoronoiMesh::checkProcessorPatchesMatch
void checkProcessorPatchesMatch(const PtrList< dictionary > &patchDicts) const
pointFields.H
OBJstream.H
pointMesh.H
indexedVertexOps.H
Foam::conformalVoronoiMesh::reorderProcessorPatches
void reorderProcessorPatches(const word &meshName, const fileName &instance, const pointField &points, faceList &faces, const wordList &patchNames, const PtrList< dictionary > &patchDicts) const
Rotate the faces on processor patches if necessary.
conformalVoronoiMesh.H
writeMesh
void writeMesh(const string &msg, const meshRefinement &meshRefiner, const meshRefinement::debugType debugLevel, const meshRefinement::writeType writeLevel)
Definition: snappyHexMesh.C:580
lookup
stressControl lookup("compactNormalStress") >> compactNormalStress