conformalVoronoiMeshIO.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2012-2020 OpenFOAM Foundation
6  \\/ M anipulation |
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  patchDicts[patchi].lookup<label>("startFace");
159  }
160  }
161 
163  {
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.localObjectPath() << 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.localObjectPath() << 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 
401 Foam::autoPtr<Foam::fvMesh> Foam::conformalVoronoiMesh::createDummyMesh
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  pointField(),
414  faceList(),
415  cellList()
416  )
417  );
418  fvMesh& mesh = meshPtr();
419 
420  List<polyPatch*> patches(patchDicts.size());
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  0, // patchSizes[p],
436  0, // patchStarts[p],
437  patchi,
438  mesh.boundaryMesh(),
439  patchDicts[patchi].lookup<label>("myProcNo"),
440  patchDicts[patchi].lookup<label>("neighbProcNo")
441  );
442  }
443  else
444  {
445  patches[patchi] = polyPatch::New
446  (
447  patchDicts[patchi].lookup("type"),
448  patchNames[patchi],
449  0, // patchSizes[p],
450  0, // patchStarts[p],
451  patchi,
452  mesh.boundaryMesh()
453  ).ptr();
454  }
455  }
456 
457  mesh.addFvPatches(patches);
458 
459  return meshPtr;
460 }
461 
462 
463 void Foam::conformalVoronoiMesh::checkProcessorPatchesMatch
464 (
465  const PtrList<dictionary>& patchDicts
466 ) const
467 {
468  // Check patch sizes
469  labelListList procPatchSizes
470  (
471  Pstream::nProcs(),
473  );
474 
475  forAll(patchDicts, patchi)
476  {
477  if
478  (
479  patchDicts.set(patchi)
480  && (
481  word(patchDicts[patchi].lookup("type"))
482  == processorPolyPatch::typeName
483  )
484  )
485  {
486  const label procNeighb =
487  patchDicts[patchi].lookup<label>("neighbProcNo");
488 
489  procPatchSizes[Pstream::myProcNo()][procNeighb]
490  = patchDicts[patchi].lookup<label>("nFaces");
491  }
492  }
493 
494  Pstream::gatherList(procPatchSizes);
495 
496  if (Pstream::master())
497  {
498  bool allMatch = true;
499 
500  forAll(procPatchSizes, proci)
501  {
502  const labelList& patchSizes = procPatchSizes[proci];
503 
504  forAll(patchSizes, patchi)
505  {
506  if (patchSizes[patchi] != procPatchSizes[patchi][proci])
507  {
508  allMatch = false;
509 
510  Info<< indent << "Patches " << proci << " and " << patchi
511  << " have different sizes: " << patchSizes[patchi]
512  << " and " << procPatchSizes[patchi][proci] << endl;
513  }
514  }
515  }
516 
517  if (allMatch)
518  {
519  Info<< indent << "All processor patches have matching numbers of "
520  << "faces" << endl;
521  }
522  }
523 }
524 
525 
526 void Foam::conformalVoronoiMesh::reorderPoints
527 (
528  pointField& points,
529  labelList& boundaryPts,
530  faceList& faces,
531  const label nInternalFaces
532 ) const
533 {
534  Info<< incrIndent << indent << "Reordering points into internal/external"
535  << endl;
536 
537  labelList oldToNew(points.size(), label(0));
538 
539  // Find points that are internal
540  for (label fI = nInternalFaces; fI < faces.size(); ++fI)
541  {
542  const face& f = faces[fI];
543 
544  forAll(f, fpI)
545  {
546  oldToNew[f[fpI]] = 1;
547  }
548  }
549 
550  const label nInternalPoints = points.size() - sum(oldToNew);
551 
552  label countInternal = 0;
553  label countExternal = nInternalPoints;
554 
555  forAll(points, pI)
556  {
557  if (oldToNew[pI] == 0)
558  {
559  oldToNew[pI] = countInternal++;
560  }
561  else
562  {
563  oldToNew[pI] = countExternal++;
564  }
565  }
566 
567  Info<< indent
568  << "Number of internal points: " << countInternal << nl
569  << indent << "Number of external points: " << countExternal
570  << decrIndent << endl;
571 
572  inplaceReorder(oldToNew, points);
573  inplaceReorder(oldToNew, boundaryPts);
574 
575  forAll(faces, fI)
576  {
577  face& f = faces[fI];
578 
579  forAll(f, fpI)
580  {
581  f[fpI] = oldToNew[f[fpI]];
582  }
583  }
584 }
585 
586 
587 void Foam::conformalVoronoiMesh::reorderProcessorPatches
588 (
589  const word& meshName,
590  const fileName& instance,
591  const pointField& points,
592  faceList& faces,
593  labelList& owner,
594  labelList& neighbour,
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_,
614  IOobject::NO_READ,
616  false
617  ),
618  patchNames,
619  patchDicts
620  )
621  );
622  const fvMesh& sortMesh = sortMeshPtr();
623 
624  // Rotation on new faces.
625  labelList rotation(faces.size(), label(0));
626  labelList faceMap(faces.size(), label(-1));
627 
628  PstreamBuffers pBufs(Pstream::commsTypes::nonBlocking);
629 
630  // Send ordering
631  forAll(sortMesh.boundaryMesh(), patchi)
632  {
633  const polyPatch& pp = sortMesh.boundaryMesh()[patchi];
634 
635  if (isA<processorPolyPatch>(pp))
636  {
637  refCast<const processorPolyPatch>(pp).initOrder
638  (
639  pBufs,
641  (
642  SubList<face>
643  (
644  faces,
645  patchDicts[patchi].lookup<label>("nFaces"),
646  patchDicts[patchi].lookup<label>("startFace")
647  ),
648  points
649  )
650  );
651  }
652  }
653 
654  pBufs.finishedSends();
655 
656  Info<< incrIndent << indent << "Face ordering initialised..." << endl;
657 
658  // Receive and calculate ordering
659  bool anyChanged = false;
660 
661  forAll(sortMesh.boundaryMesh(), patchi)
662  {
663  const polyPatch& pp = sortMesh.boundaryMesh()[patchi];
664 
665  if (isA<processorPolyPatch>(pp))
666  {
667  const label nPatchFaces =
668  patchDicts[patchi].lookup<label>("nFaces");
669  const label patchStartFace =
670  patchDicts[patchi].lookup<label>("startFace");
671 
672  labelList patchFaceMap(nPatchFaces, label(-1));
673  labelList patchFaceRotation(nPatchFaces, label(0));
674 
675  bool changed = refCast<const processorPolyPatch>(pp).order
676  (
677  pBufs,
679  (
680  SubList<face>
681  (
682  faces,
683  nPatchFaces,
684  patchStartFace
685  ),
686  points
687  ),
688  patchFaceMap,
689  patchFaceRotation
690  );
691 
692  if (changed)
693  {
694  // Merge patch face reordering into mesh face reordering table
695  forAll(patchFaceRotation, patchFacei)
696  {
697  rotation[patchFacei + patchStartFace]
698  = patchFaceRotation[patchFacei];
699  }
700 
701  forAll(patchFaceMap, patchFacei)
702  {
703  if (patchFaceMap[patchFacei] != patchFacei)
704  {
705  faceMap[patchFacei + patchStartFace]
706  = patchFaceMap[patchFacei] + patchStartFace;
707  }
708  }
709 
710  anyChanged = true;
711  }
712  }
713  }
714 
715  Info<< incrIndent << indent << "Faces matched." << endl;
716 
717  reduce(anyChanged, orOp<bool>());
718 
719  if (anyChanged)
720  {
721  label nReorderedFaces = 0;
722 
723  forAll(faceMap, facei)
724  {
725  if (faceMap[facei] != -1)
726  {
727  nReorderedFaces++;
728  }
729  }
730 
731  if (nReorderedFaces > 0)
732  {
733  inplaceReorder(faceMap, faces);
734  inplaceReorder(faceMap, owner);
735  inplaceReorder(faceMap, neighbour);
736  }
737 
738  // Rotate faces (rotation is already in new face indices).
739  label nRotated = 0;
740 
741  forAll(rotation, facei)
742  {
743  if (rotation[facei] != 0)
744  {
745  faces[facei] = rotateList(faces[facei], rotation[facei]);
746  nRotated++;
747  }
748  }
749 
750  Info<< indent << returnReduce(nReorderedFaces, sumOp<label>())
751  << " faces have been reordered" << nl
752  << indent << returnReduce(nRotated, sumOp<label>())
753  << " faces have been rotated"
754  << decrIndent << decrIndent
755  << decrIndent << decrIndent << endl;
756  }
757 }
758 
759 
761 (
762  const word& meshName,
763  const fileName& instance,
764  pointField& points,
765  labelList& boundaryPts,
766  faceList& faces,
767  labelList& owner,
768  labelList& neighbour,
769  const wordList& patchNames,
770  const PtrList<dictionary>& patchDicts,
771  const pointField& cellCentres,
772  PackedBoolList& boundaryFacesToRemove
773 ) const
774 {
775  if (foamyHexMeshControls().objOutput())
776  {
778  (
779  time().path()/word(meshName + ".obj"),
780  points,
781  faces
782  );
783  }
784 
785  const label nInternalFaces = patchDicts[0].lookup<label>("startFace");
786 
787  reorderPoints(points, boundaryPts, faces, nInternalFaces);
788 
789  if (Pstream::parRun())
790  {
791  reorderProcessorPatches
792  (
793  meshName,
794  instance,
795  points,
796  faces,
797  owner,
798  neighbour,
799  patchNames,
800  patchDicts
801  );
802  }
803 
804  Info<< incrIndent;
805  Info<< indent << "Constructing mesh" << endl;
806 
807  timeCheck("Before fvMesh construction");
808 
809  fvMesh mesh
810  (
811  IOobject
812  (
813  meshName,
814  instance,
815  runTime_,
816  IOobject::NO_READ,
818  ),
819  Foam::move(points),
820  Foam::move(faces),
821  Foam::move(owner),
822  Foam::move(neighbour)
823  );
824 
825  Info<< indent << "Adding patches to mesh" << endl;
826 
827  List<polyPatch*> patches(patchNames.size());
828 
829  label nValidPatches = 0;
830 
831  forAll(patches, p)
832  {
833  label totalPatchSize = patchDicts[p].lookup<label>("nFaces");
834 
835  if
836  (
837  patchDicts.set(p)
838  && (
839  word(patchDicts[p].lookup("type"))
840  == processorPolyPatch::typeName
841  )
842  )
843  {
844  const_cast<dictionary&>(patchDicts[p]).set
845  (
846  "transform",
847  "coincidentFullMatch"
848  );
849 
850  // Do not create empty processor patches
851  if (totalPatchSize > 0)
852  {
853  patches[nValidPatches] = new processorPolyPatch
854  (
855  patchNames[p],
856  patchDicts[p],
857  nValidPatches,
858  mesh.boundaryMesh(),
859  processorPolyPatch::typeName
860  );
861 
862  nValidPatches++;
863  }
864  }
865  else
866  {
867  // Check that the patch is not empty on every processor
868  reduce(totalPatchSize, sumOp<label>());
869 
870  if (totalPatchSize > 0)
871  {
872  patches[nValidPatches] = polyPatch::New
873  (
874  patchNames[p],
875  patchDicts[p],
876  nValidPatches,
877  mesh.boundaryMesh()
878  ).ptr();
879 
880  nValidPatches++;
881  }
882  }
883  }
884 
885  patches.setSize(nValidPatches);
886 
887  mesh.addFvPatches(patches);
888 
889 
890 
891  // Add zones to the mesh
892  addZones(mesh, cellCentres);
893 
894 
895 
896  Info<< indent << "Add pointZones" << endl;
897  {
898  label sz = mesh.pointZones().size();
899 
900  DynamicList<label> bPts(boundaryPts.size());
901 
903  {
904  forAll(boundaryPts, ptI)
905  {
906  const label& bPtType = boundaryPts[ptI];
907 
908  if (bPtType == typeI)
909  {
910  bPts.append(ptI);
911  }
912  }
913 
914 // syncTools::syncPointList(mesh, bPts, maxEqOp<label>(), -1);
915 
916  Info<< incrIndent << indent
917  << "Adding " << bPts.size()
918  << " points of type " << dualMeshPointTypeNames_.words()[typeI]
919  << decrIndent << endl;
920 
921  mesh.pointZones().append
922  (
923  new pointZone
924  (
926  bPts,
927  sz + typeI,
928  mesh.pointZones()
929  )
930  );
931 
932  bPts.clear();
933  }
934  }
935 
936 
937 
938  // Add indirectPatchFaces to a face zone
939  Info<< indent << "Adding indirect patch faces set" << endl;
940 
942  (
943  mesh,
944  boundaryFacesToRemove,
945  orEqOp<unsigned int>()
946  );
947 
948  labelList addr(boundaryFacesToRemove.count());
949  label count = 0;
950 
951  forAll(boundaryFacesToRemove, facei)
952  {
953  if (boundaryFacesToRemove[facei])
954  {
955  addr[count++] = facei;
956  }
957  }
958 
959  addr.setSize(count);
960 
961  faceSet indirectPatchFaces
962  (
963  mesh,
964  "indirectPatchFaces",
965  addr,
967  );
968 
969  indirectPatchFaces.sync(mesh);
970 
971 
972  Info<< decrIndent;
973 
974  timeCheck("Before fvMesh filtering");
975 
976  autoPtr<polyMeshFilter> meshFilter;
977 
978  label nInitialBadFaces = 0;
979 
980  if (foamyHexMeshControls().filterEdges())
981  {
982  Info<< nl << "Filtering edges on polyMesh" << nl << endl;
983 
984  meshFilter.reset(new polyMeshFilter(mesh, boundaryPts));
985 
986  // Filter small edges only. This reduces the number of faces so that
987  // the face filtering is sped up.
988  nInitialBadFaces = meshFilter().filterEdges(0);
989  {
990  const autoPtr<fvMesh>& newMesh = meshFilter().filteredMesh();
991 
992  polyTopoChange meshMod(newMesh());
993 
994  autoPtr<mapPolyMesh> map = meshMod.changeMesh(mesh, false);
995 
996  polyMeshFilter::copySets(newMesh(), mesh);
997  }
998  }
999 
1000  if (foamyHexMeshControls().filterFaces())
1001  {
1002  labelIOList boundaryPtsIO
1003  (
1004  IOobject
1005  (
1006  "pointPriority",
1007  instance,
1008  time(),
1009  IOobject::NO_READ,
1011  ),
1012  labelList(mesh.nPoints(), labelMin)
1013  );
1014 
1015  forAll(mesh.points(), ptI)
1016  {
1017  boundaryPtsIO[ptI] = mesh.pointZones().whichZone(ptI);
1018  }
1019 
1020 
1021  Info<< nl << "Filtering faces on polyMesh" << nl << endl;
1022 
1023  meshFilter.reset(new polyMeshFilter(mesh, boundaryPtsIO));
1024 
1025  meshFilter().filter(nInitialBadFaces);
1026  {
1027  const autoPtr<fvMesh>& newMesh = meshFilter().filteredMesh();
1028 
1029  polyTopoChange meshMod(newMesh());
1030 
1031  autoPtr<mapPolyMesh> map = meshMod.changeMesh(mesh, false);
1032 
1033  polyMeshFilter::copySets(newMesh(), mesh);
1034  }
1035  }
1036 
1037  timeCheck("After fvMesh filtering");
1038 
1039  mesh.setInstance(instance);
1040 
1041  if (!mesh.write())
1042  {
1044  << "Failed writing polyMesh."
1045  << exit(FatalError);
1046  }
1047  else
1048  {
1049  Info<< nl << "Written filtered mesh to "
1050  << mesh.polyMesh::instance() << nl
1051  << endl;
1052  }
1053 
1054  {
1055  pointScalarField boundaryPtsScalarField
1056  (
1057  IOobject
1058  (
1059  "boundaryPoints_collapsed",
1060  instance,
1061  time(),
1062  IOobject::NO_READ,
1064  ),
1065  pointMesh::New(mesh),
1066  scalar(labelMin)
1067  );
1068 
1069  labelIOList boundaryPtsIO
1070  (
1071  IOobject
1072  (
1073  "pointPriority",
1074  instance,
1075  time(),
1076  IOobject::NO_READ,
1078  ),
1079  labelList(mesh.nPoints(), labelMin)
1080  );
1081 
1082  forAll(mesh.points(), ptI)
1083  {
1084  boundaryPtsScalarField[ptI] = mesh.pointZones().whichZone(ptI);
1085  boundaryPtsIO[ptI] = mesh.pointZones().whichZone(ptI);
1086  }
1087 
1088  boundaryPtsScalarField.write();
1089  boundaryPtsIO.write();
1090  }
1091 
1092 // writeCellSizes(mesh);
1093 
1094 // writeCellAlignments(mesh);
1095 
1096 // writeCellCentres(mesh);
1097 
1099 }
1100 
1101 
1103 (
1104  const fvMesh& mesh
1105 ) const
1106 {
1107  {
1108  timeCheck("Start writeCellSizes");
1109 
1110  Info<< nl << "Create targetCellSize volScalarField" << endl;
1111 
1112  volScalarField targetCellSize
1113  (
1114  IOobject
1115  (
1116  "targetCellSize",
1117  mesh.polyMesh::instance(),
1118  mesh,
1121  ),
1122  mesh,
1124  zeroGradientFvPatchScalarField::typeName
1125  );
1126 
1127  scalarField& cellSize = targetCellSize.primitiveFieldRef();
1128 
1129  const vectorField& C = mesh.cellCentres();
1130 
1131  forAll(cellSize, i)
1132  {
1133  cellSize[i] = cellShapeControls().cellSize(C[i]);
1134  }
1135 
1136  // Info<< nl << "Create targetCellVolume volScalarField" << endl;
1137 
1138  // volScalarField targetCellVolume
1139  // (
1140  // IOobject
1141  // (
1142  // "targetCellVolume",
1143  // mesh.polyMesh::instance(),
1144  // mesh,
1145  // IOobject::NO_READ,
1146  // IOobject::AUTO_WRITE
1147  // ),
1148  // mesh,
1149  // dimensionedScalar(dimLength, 0),
1150  // zeroGradientFvPatchScalarField::typeName
1151  // );
1152 
1153  // targetCellVolume.primitiveFieldRef() = pow3(cellSize);
1154 
1155  // Info<< nl << "Create actualCellVolume volScalarField" << endl;
1156 
1157  // volScalarField actualCellVolume
1158  // (
1159  // IOobject
1160  // (
1161  // "actualCellVolume",
1162  // mesh.polyMesh::instance(),
1163  // mesh,
1164  // IOobject::NO_READ,
1165  // IOobject::AUTO_WRITE
1166  // ),
1167  // mesh,
1168  // dimensionedScalar(dimVolume, 0),
1169  // zeroGradientFvPatchScalarField::typeName
1170  // );
1171 
1172  // actualCellVolume.primitiveFieldRef() = mesh.cellVolumes();
1173 
1174  // Info<< nl << "Create equivalentCellSize volScalarField" << endl;
1175 
1176  // volScalarField equivalentCellSize
1177  // (
1178  // IOobject
1179  // (
1180  // "equivalentCellSize",
1181  // mesh.polyMesh::instance(),
1182  // mesh,
1183  // IOobject::NO_READ,
1184  // IOobject::AUTO_WRITE
1185  // ),
1186  // mesh,
1187  // dimensionedScalar(dimLength, 0),
1188  // zeroGradientFvPatchScalarField::typeName
1189  // );
1190 
1191  // equivalentCellSize.primitiveFieldRef() = pow
1192  // (
1193  // actualCellVolume.primitiveField(),
1194  // 1.0/3.0
1195  // );
1196 
1197  targetCellSize.correctBoundaryConditions();
1198  // targetCellVolume.correctBoundaryConditions();
1199  // actualCellVolume.correctBoundaryConditions();
1200  // equivalentCellSize.correctBoundaryConditions();
1201 
1202  targetCellSize.write();
1203  // targetCellVolume.write();
1204  // actualCellVolume.write();
1205  // equivalentCellSize.write();
1206  }
1207 
1208  // {
1209  // polyMesh tetMesh
1210  // (
1211  // IOobject
1212  // (
1213  // "tetDualMesh",
1214  // runTime_.constant(),
1215  // runTime_,
1216  // IOobject::MUST_READ
1217  // )
1218  // );
1219 
1220  // pointMesh ptMesh(tetMesh);
1221 
1222  // pointScalarField ptTargetCellSize
1223  // (
1224  // IOobject
1225  // (
1226  // "ptTargetCellSize",
1227  // runTime_.timeName(),
1228  // tetMesh,
1229  // IOobject::NO_READ,
1230  // IOobject::AUTO_WRITE
1231  // ),
1232  // ptMesh,
1233  // dimensionedScalar(dimLength, 0),
1234  // pointPatchVectorField::calculatedType()
1235  // );
1236 
1237  // scalarField& cellSize = ptTargetCellSize.primitiveFieldRef();
1238 
1239  // const vectorField& P = tetMesh.points();
1240 
1241  // forAll(cellSize, i)
1242  // {
1243  // cellSize[i] = cellShapeControls().cellSize(P[i]);
1244  // }
1245 
1246  // ptTargetCellSize.write();
1247  // }
1248 }
1249 
1250 
1252 (
1253  const fvMesh& mesh
1254 ) const
1255 {
1256 // Info<< nl << "Create cellAlignments volTensorField" << endl;
1257 //
1258 // volTensorField cellAlignments
1259 // (
1260 // IOobject
1261 // (
1262 // "cellAlignments",
1263 // mesh.polyMesh::instance(),
1264 // mesh,
1265 // IOobject::NO_READ,
1266 // IOobject::AUTO_WRITE
1267 // ),
1268 // mesh,
1269 // tensor::I,
1270 // zeroGradientFvPatchTensorField::typeName
1271 // );
1272 //
1273 // tensorField& cellAlignment = cellAlignments.primitiveFieldRef();
1274 //
1275 // const vectorField& C = mesh.cellCentres();
1276 //
1277 // vectorField xDir(cellAlignment.size());
1278 // vectorField yDir(cellAlignment.size());
1279 // vectorField zDir(cellAlignment.size());
1280 //
1281 // forAll(cellAlignment, i)
1282 // {
1283 // cellAlignment[i] = cellShapeControls().cellAlignment(C[i]);
1284 // xDir[i] = cellAlignment[i] & vector(1, 0, 0);
1285 // yDir[i] = cellAlignment[i] & vector(0, 1, 0);
1286 // zDir[i] = cellAlignment[i] & vector(0, 0, 1);
1287 // }
1288 //
1289 // OFstream xStr("xDir.obj");
1290 // OFstream yStr("yDir.obj");
1291 // OFstream zStr("zDir.obj");
1292 //
1293 // forAll(xDir, i)
1294 // {
1295 // meshTools::writeOBJ(xStr, C[i], C[i] + xDir[i]);
1296 // meshTools::writeOBJ(yStr, C[i], C[i] + yDir[i]);
1297 // meshTools::writeOBJ(zStr, C[i], C[i] + zDir[i]);
1298 // }
1299 //
1300 // cellAlignments.correctBoundaryConditions();
1301 //
1302 // cellAlignments.write();
1303 }
1304 
1305 
1307 (
1308  const fvMesh& mesh
1309 ) const
1310 {
1311  Info<< "Writing components of cellCentre positions to volScalarFields"
1312  << " ccx, ccy, ccz in " << runTime_.timeName() << endl;
1313 
1314  for (direction i=0; i<vector::nComponents; i++)
1315  {
1316  volScalarField cci
1317  (
1318  IOobject
1319  (
1320  "cc" + word(vector::componentNames[i]),
1321  runTime_.timeName(),
1322  mesh,
1325  ),
1326  mesh.C().component(i)
1327  );
1328 
1329  cci.write();
1330  }
1331 }
1332 
1333 
1335 (
1336  const polyMesh& mesh
1337 ) const
1338 {
1339  timeCheck("Start findRemainingProtrusionSet");
1340 
1341  const polyBoundaryMesh& patches = mesh.boundaryMesh();
1342 
1343  labelHashSet protrudingBoundaryPoints;
1344 
1345  forAll(patches, patchi)
1346  {
1347  const polyPatch& patch = patches[patchi];
1348 
1349  forAll(patch.localPoints(), pLPI)
1350  {
1351  label meshPtI = patch.meshPoints()[pLPI];
1352 
1353  const Foam::point& pt = patch.localPoints()[pLPI];
1354 
1355  if
1356  (
1357  geometryToConformTo_.wellOutside
1358  (
1359  pt,
1360  sqr(targetCellSize(pt))
1361  )
1362  )
1363  {
1364  protrudingBoundaryPoints.insert(meshPtI);
1365  }
1366  }
1367  }
1368 
1369  cellSet protrudingCells
1370  (
1371  mesh,
1372  "foamyHexMesh_remainingProtrusions",
1373  mesh.nCells()/1000
1374  );
1375 
1376  forAllConstIter(labelHashSet, protrudingBoundaryPoints, iter)
1377  {
1378  const label pointi = iter.key();
1379  const labelList& pCells = mesh.pointCells()[pointi];
1380 
1381  forAll(pCells, pCI)
1382  {
1383  protrudingCells.insert(pCells[pCI]);
1384  }
1385  }
1386 
1387  label protrudingCellsSize = protrudingCells.size();
1388 
1389  reduce(protrudingCellsSize, sumOp<label>());
1390 
1391  if (foamyHexMeshControls().objOutput() && protrudingCellsSize > 0)
1392  {
1393  Info<< nl << "Found " << protrudingCellsSize
1394  << " cells protruding from the surface, writing cellSet "
1395  << protrudingCells.name()
1396  << endl;
1397 
1398  protrudingCells.write();
1399  }
1400 
1401  return Foam::move(protrudingCells);
1402 }
1403 
1404 
1405 void Foam::conformalVoronoiMesh::writePointPairs
1406 (
1407  const fileName& fName
1408 ) const
1409 {
1410  OBJstream os(fName);
1411 
1412  for
1413  (
1414  Delaunay::Finite_edges_iterator eit = finite_edges_begin();
1415  eit != finite_edges_end();
1416  ++eit
1417  )
1418  {
1419  Cell_handle c = eit->first;
1420  Vertex_handle vA = c->vertex(eit->second);
1421  Vertex_handle vB = c->vertex(eit->third);
1422 
1423  if (ptPairs_.isPointPair(vA, vB))
1424  {
1425  os.write
1426  (
1427  linePointRef(topoint(vA->point()), topoint(vB->point()))
1428  );
1429  }
1430  }
1431 }
1432 
1433 
1434 // ************************************************************************* //
GeometricField< scalar, pointPatchField, pointMesh > pointScalarField
Definition: pointFields.H:49
labelHashSet findRemainingProtrusionSet(const polyMesh &mesh) const
Find the cellSet of the boundary cells which have points that.
static wordList words()
The set of names as a list of words.
Definition: NamedEnum.C:105
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:57
Delaunay::Cell_handle Cell_handle
void inplaceReorder(const labelUList &oldToNew, ListType &)
Inplace reorder the elements of a list.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
virtual Ostream & write(const char)=0
Write character.
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
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:221
Delaunay::Vertex_handle Vertex_handle
pointFromPoint topoint(const Point &P)
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
ListType rotateList(const ListType &list, const label n)
Rotate a list by n places. If n is positive rotate clockwise/right/down.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
static fvMesh * meshPtr
Definition: globalFoam.H:52
const cellShapeControl & cellShapeControls() const
Return the cellShapeControl object.
uint8_t direction
Definition: direction.H:45
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:429
List< face > faceList
Definition: faceListFwd.H:43
HashTable< label, labelPair, FixedList< label, 2 >::Hash<> > labelTolabelPairHashTable
Definition: DelaunayMesh.H:90
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:309
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:423
static void timeCheck(const Time &runTime, const string &description=string::null, const bool check=true)
Write the elapsedCpuTime and memory usage, with an optional.
void writeCellCentres(const fvMesh &mesh) const
Calculate and write the cell centres.
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:622
patches[0]
void writeMesh(const fileName &instance)
Prepare data and call writeMesh for polyMesh and.
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
Various functions to operate on Lists.
static const pointMesh & New(const polyMesh &mesh)
Definition: MeshObject.C:44
Useful combination of include files which define Sin, Sout and Serr and the use of IO streams general...
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
HashSet< label, Hash< label > > labelHashSet
A HashSet with label keys.
Definition: HashSet.H:211
void writeObjMesh(const fileName &fName, const pointField &points, const faceList &faces)
Write an OBJ mesh consisting of points and faces.
static const direction nComponents
Number of components in this vector space.
Definition: VectorSpace.H:99
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
stressControl lookup("compactNormalStress") >> compactNormalStress
void writeCellSizes(const fvMesh &mesh) const
Calculate and write a field of the target cell size,.
dynamicFvMesh & mesh
line< point, const point & > linePointRef
Line using referred points.
Definition: linePointRef.H:45
const pointField & points
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
wordList patchNames(nPatches)
static const word null
An empty word.
Definition: word.H:77
List< label > labelList
A List of labels.
Definition: labelList.H:56
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
static const char nl
Definition: Ostream.H:260
Delaunay ::Finite_cells_iterator Finite_cells_iterator
Switch writeCellShapeControlMesh() const
Write cellShapeControlMesh at output time.
Definition: cvControlsI.H:218
Field< bool > wellOutside(const pointField &samplePts, const scalarField &testDistSqr) const
Check if point is outside surfaces to conform to by at least.
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition: Ostream.H:235
labelList f(nPoints)
PrimitivePatch< SubList< face >, const pointField & > primitivePatch
Addressing for a faceList slice.
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
bool isPointPair(const Vertex_handle &vA, const Vertex_handle &vB) const
List< word > wordList
A List of words.
Definition: fileName.H:54
static void syncFaceList(const polyMesh &mesh, UList< T > &l, const CombineOp &cop)
Synchronize values on all mesh faces.
Definition: syncTools.H:381
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:399
static label nProcs(const label communicator=0)
Number of processes in parallel run.
Definition: UPstream.H:411
label patchi
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
messageStream Info
void writeCellAlignments(const fvMesh &mesh) const
static const NamedEnum< dualMeshPointType, 5 > dualMeshPointTypeNames_
cellShapeControlMesh & shapeControlMesh()
Field< vector > vectorField
Specialisation of Field<T> for vector.
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
const cvControls & foamyHexMeshControls() const
Return the foamyHexMeshControls object.
const Time & time() const
Return the Time object.
PtrList< dictionary > patchDicts
Definition: readKivaGrid.H:537
virtual bool write(const bool write=true) const
Write using setting from DB.
volScalarField & p
virtual void print(Ostream &) const
Print description of IOstream to Ostream.
Definition: IOstream.C:126
void writeInternalDelaunayVertices(const fileName &instance, const Triangulation &t)
Write the internal Delaunay vertices of the tessellation as a.
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition: Ostream.H:228
static void gatherList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Gather data but keep individual values separate.
List< cell > cellList
list of cells
Definition: cellList.H:42
static const label labelMin
Definition: label.H:61
IOList< label > labelIOList
Label container classes.
Definition: labelIOList.H:42
autoPtr< polyMesh > createMesh(const fileName &name, labelTolabelPairHashTable &vertexMap, labelList &cellMap, const bool writeDelaunayData=true) const
Create an fvMesh from the triangulation.
scalar cellSize(const point &pt) const
Return the cell size at the given location.
fileName path(UMean.rootPath()/UMean.caseName()/functionObjects::writeFile::outputPrefix/"graphs"/UMean.instance())