Loading...
Searching...
No Matches
reconstructParMesh.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 | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2016 OpenFOAM Foundation
9 Copyright (C) 2016-2025 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27Application
28 reconstructParMesh
29
30Group
31 grpParallelUtilities
32
33Description
34 Reconstructs a mesh using geometric information only.
35
36 Writes point/face/cell procAddressing so afterwards reconstructPar can be
37 used to reconstruct fields.
38
39Usage
40 \b reconstructParMesh [OPTION]
41
42 Options:
43 - \par -fullMatch
44 Does geometric matching on all boundary faces. Assumes no point
45 ordering
46
47 - \par -procMatch
48 Assumes processor patches already in face order but not point order.
49 This is the pre v2106 default behaviour but might be removed if the new
50 topological method works well
51
52 - \par -mergeTol <tol>
53 Specifies non-default merge tolerance (fraction of mesh bounding box)
54 for above options
55
56 The default is to assume all processor boundaries are correctly ordered
57 (both faces and points) in which case no merge tolerance is needed.
58
59\*---------------------------------------------------------------------------*/
60
61#include "argList.H"
62#include "timeSelector.H"
63
64#include "IOobjectList.H"
65#include "labelIOList.H"
66#include "processorPolyPatch.H"
67#include "mapAddedPolyMesh.H"
68#include "polyMeshAdder.H"
69#include "faceCoupleInfo.H"
70#include "fvMeshAdder.H"
71#include "polyTopoChange.H"
72#include "topoSet.H"
73#include "regionProperties.H"
74#include "fvMeshTools.H"
75
76#include "faMeshReconstructor.H"
77#include "ignoreFaPatch.H"
78
79using namespace Foam;
80
81// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
82
83// Tolerance (as fraction of the bounding box). Needs to be fairly lax since
84// usually meshes get written with limited precision (6 digits)
85static const scalar defaultMergeTol = 1e-7;
86
87
88// Determine which faces are coupled. Uses geometric merge distance.
89// Looks either at all boundaryFaces (fullMatch) or only at the
90// procBoundaries for proci. Assumes that masterMesh contains already merged
91// all the processors < proci.
92autoPtr<faceCoupleInfo> determineCoupledFaces
93(
94 const bool fullMatch,
95 const label masterMeshProcStart,
96 const label masterMeshProcEnd,
97 const polyMesh& masterMesh,
98 const label meshToAddProcStart,
99 const label meshToAddProcEnd,
100 const polyMesh& meshToAdd,
101 const scalar mergeDist
102)
103{
104 if (fullMatch || masterMesh.nCells() == 0)
105 {
107 (
108 masterMesh,
109 meshToAdd,
110 mergeDist, // Absolute merging distance
111 true // Matching faces identical
112 );
113 }
114 else
115 {
116 // Pick up all patches on masterMesh ending in "toDDD" where DDD is
117 // the processor number proci.
118
119 const polyBoundaryMesh& masterPatches = masterMesh.boundaryMesh();
120
121
122 DynamicList<label> masterFaces(masterMesh.nBoundaryFaces());
123
124 for (const polyPatch& pp : masterPatches)
125 {
127 {
128 for
129 (
130 label proci=meshToAddProcStart;
131 proci<meshToAddProcEnd;
132 proci++
133 )
134 {
135 const string toProcString("to" + Foam::name(proci));
136 if (pp.name().ends_with(toProcString))
137 {
138 label meshFacei = pp.start();
139 forAll(pp, i)
140 {
141 masterFaces.append(meshFacei++);
142 }
143 break;
144 }
145 }
146
147 }
148 }
149
150
151 // Pick up all patches on meshToAdd ending in "procBoundaryDDDtoYYY"
152 // where DDD is the processor number proci and YYY is < proci.
153
154 const polyBoundaryMesh& addPatches = meshToAdd.boundaryMesh();
155
156 DynamicList<label> addFaces(meshToAdd.nBoundaryFaces());
157
158 for (const polyPatch& pp : addPatches)
159 {
161 {
162 bool isConnected = false;
163
164 for
165 (
166 label mergedProci=masterMeshProcStart;
167 !isConnected && (mergedProci < masterMeshProcEnd);
168 mergedProci++
169 )
170 {
171 for
172 (
173 label proci = meshToAddProcStart;
174 proci < meshToAddProcEnd;
175 proci++
176 )
177 {
178 const word fromProcString
179 (
180 processorPolyPatch::newName(proci, mergedProci)
181 );
182
183 if (pp.name() == fromProcString)
184 {
185 isConnected = true;
186 break;
187 }
188 }
189 }
190
191 if (isConnected)
192 {
193 label meshFacei = pp.start();
194 forAll(pp, i)
195 {
196 addFaces.append(meshFacei++);
197 }
198 }
199 }
200 }
201
203 (
204 masterMesh,
205 masterFaces,
206 meshToAdd,
207 addFaces,
208 mergeDist, // Absolute merging distance
209 true, // Matching faces identical?
210 false, // If perfect match are faces already ordered
211 // (e.g. processor patches)
212 false // are faces each on separate patch?
213 );
214 }
215}
216
217
218autoPtr<mapPolyMesh> mergeSharedPoints
219(
220 const scalar mergeDist,
221 polyMesh& mesh,
222 labelListList& pointProcAddressing
223)
224{
225 // Find out which sets of points get merged and create a map from
226 // mesh point to unique point.
227 Map<label> pointToMaster
228 (
230 (
231 mesh,
232 mergeDist
233 )
234 );
235
236 Info<< "mergeSharedPoints : detected " << pointToMaster.size()
237 << " points that are to be merged." << endl;
238
239 if (returnReduceAnd(pointToMaster.empty()))
240 {
241 return nullptr;
242 }
243
244 polyTopoChange meshMod(mesh);
245
246 fvMeshAdder::mergePoints(mesh, pointToMaster, meshMod);
247
248 // Change the mesh (no inflation). Note: parallel comms allowed.
249 autoPtr<mapPolyMesh> map = meshMod.changeMesh(mesh, false, true);
250
251 // Update fields. No inflation, parallel sync.
252 mesh.updateMesh(map());
253
254 // pointProcAddressing give indices into the master mesh so adapt them
255 // for changed point numbering.
256
257 // Adapt constructMaps for merged points.
258 forAll(pointProcAddressing, proci)
259 {
260 labelList& constructMap = pointProcAddressing[proci];
261
262 forAll(constructMap, i)
263 {
264 label oldPointi = constructMap[i];
265
266 // New label of point after changeMesh.
267 label newPointi = map().reversePointMap()[oldPointi];
268
269 if (newPointi < -1)
270 {
271 constructMap[i] = -newPointi-2;
272 }
273 else if (newPointi >= 0)
274 {
275 constructMap[i] = newPointi;
276 }
277 else
278 {
280 << "Problem. oldPointi:" << oldPointi
281 << " newPointi:" << newPointi << abort(FatalError);
282 }
283 }
284 }
285
286 return map;
287}
288
289
290boundBox procBounds
291(
292 const PtrList<Time>& databases,
293 const word& regionName
294)
295{
296 boundBox bb;
297
298 for (const Time& procDb : databases)
299 {
300 fileName pointsInstance
301 (
302 procDb.findInstance(polyMesh::meshDir(regionName), "points")
303 );
304
305 if (pointsInstance != procDb.timeName())
306 {
308 << "Your time was specified as " << procDb.timeName()
309 << " but there is no polyMesh/points in that time." << nl
310 << "(points file at " << pointsInstance << ')' << nl
311 << "Please rerun with the correct time specified"
312 << " (through the -constant, -time or -latestTime "
313 << "(at your option)."
314 << endl << exit(FatalError);
315 }
316
317 Info<< "Reading points from "
318 << procDb.caseName()
319 << " for time = " << procDb.timeName()
320 << nl << endl;
321
323 (
325 (
326 "points",
327 pointsInstance,
329 procDb,
333 )
334 );
335
336 bb.add(points);
337 }
338
339 return bb;
340}
341
342
343void writeDistribution
344(
345 Time& runTime,
346 const fvMesh& masterMesh,
347 const labelListList& cellProcAddressing
348)
349{
350 // Write the decomposition as labelList for use with 'manual'
351 // decomposition method.
352 labelIOList cellDecomposition
353 (
355 (
356 "cellDecomposition",
357 masterMesh.facesInstance(),
358 masterMesh,
362 ),
363 masterMesh.nCells()
364 );
365
366 forAll(cellProcAddressing, proci)
367 {
368 const labelList& pCells = cellProcAddressing[proci];
369 labelUIndList(cellDecomposition, pCells) = proci;
370 }
371
372 cellDecomposition.write();
373
374 Info<< nl << "Wrote decomposition to "
375 << cellDecomposition.objectRelPath()
376 << " for use in manual decomposition." << endl;
377
378 // Write as volScalarField for postprocessing.
379 // Change time to 0 if was 'constant'
380 {
381 const scalar oldTime = runTime.value();
382 const label oldIndex = runTime.timeIndex();
383 if (runTime.timeName() == runTime.constant() && oldIndex == 0)
384 {
385 runTime.setTime(0, oldIndex+1);
386 }
387
388 volScalarField cellDist
389 (
391 (
392 "cellDist",
393 runTime.timeName(),
394 masterMesh,
398 ),
399 masterMesh,
402 );
403
404 forAll(cellDecomposition, celli)
405 {
406 cellDist[celli] = cellDecomposition[celli];
407 }
408
409 cellDist.correctBoundaryConditions();
410 cellDist.write();
411
412 Info<< nl << "Wrote decomposition to "
413 << cellDist.objectRelPath()
414 << " (volScalarField) for visualization."
415 << endl;
416
417 // Restore time
418 runTime.setTime(oldTime, oldIndex);
419 }
420}
421
422
423void writeMesh
424(
425 const fvMesh& mesh,
426 const labelListList& cellProcAddressing
427)
428{
429 const fileName outputDir
430 (
431 mesh.time().path()
432 / mesh.time().timeName()
433 / mesh.regionName()
435 );
436
437 Info<< nl << "Writing merged mesh to "
438 << mesh.time().relativePath(outputDir) << nl << endl;
439
440 if (!mesh.write())
441 {
443 << "Failed writing polyMesh."
444 << exit(FatalError);
445 }
447}
448
449
450void writeMaps
451(
452 const label masterInternalFaces,
453 const labelUList& masterOwner,
454 const polyMesh& procMesh,
455 const labelUList& cellProcAddressing,
456 const labelUList& faceProcAddressing,
457 const labelUList& pointProcAddressing,
458 const labelUList& boundProcAddressing
459)
460{
461 const fileName outputDir
462 (
463 procMesh.time().caseName()
464 / procMesh.facesInstance()
465 / procMesh.regionName()
467 );
468
469 IOobject ioAddr
470 (
471 "procAddressing",
472 procMesh.facesInstance(),
474 procMesh,
478 );
479
480
481 Info<< "Writing addressing : " << outputDir << nl;
482
483 // From processor point to reconstructed mesh point
484
485 Info<< " pointProcAddressing" << endl;
486 ioAddr.rename("pointProcAddressing");
487 labelIOList::writeContents(ioAddr, pointProcAddressing);
488
489 // From processor face to reconstructed mesh face
490 // - add turning index to faceProcAddressing.
491 // See reconstructPar for meaning of turning index.
492
493 labelList faceProcAddr(faceProcAddressing);
494
495 forAll(faceProcAddr, procFacei)
496 {
497 const label masterFacei = faceProcAddr[procFacei];
498
499 if
500 (
501 !procMesh.isInternalFace(procFacei)
502 && masterFacei < masterInternalFaces
503 )
504 {
505 // proc face is now external but used to be internal face.
506 // Check if we have owner or neighbour.
507
508 label procOwn = procMesh.faceOwner()[procFacei];
509 label masterOwn = masterOwner[masterFacei];
510
511 if (cellProcAddressing[procOwn] == masterOwn)
512 {
513 // No turning. Offset by 1.
514 faceProcAddr[procFacei]++;
515 }
516 else
517 {
518 // Turned face.
519 faceProcAddr[procFacei] = -1 - faceProcAddr[procFacei];
520 }
521 }
522 else
523 {
524 // No turning. Offset by 1.
525 faceProcAddr[procFacei]++;
526 }
527 }
528
529 Info<< " faceProcAddressing" << endl;
530 ioAddr.rename("faceProcAddressing");
531 labelIOList::writeContents(ioAddr, faceProcAddr);
532 faceProcAddr.clear();
533
534 // From processor cell to reconstructed mesh cell
535 Info<< " cellProcAddressing" << endl;
536 ioAddr.rename("cellProcAddressing");
537 labelIOList::writeContents(ioAddr, cellProcAddressing);
538
539
540 // From processor patch to reconstructed mesh patch
541 Info<< " boundaryProcAddressing" << endl;
542 ioAddr.rename("boundaryProcAddressing");
543 labelIOList::writeContents(ioAddr, boundProcAddressing);
544
545 Info<< endl;
546}
547
548
549void printWarning()
550{
551 Info<<
552"Merge individual processor meshes back into one master mesh.\n"
553"Use if the original master mesh has been deleted or the processor meshes\n"
554"have been modified (topology change).\n"
555"This tool will write the resulting mesh to a new time step and construct\n"
556"xxxxProcAddressing files in the processor meshes so reconstructPar can be\n"
557"used to regenerate the fields on the master mesh.\n\n"
558"Not well tested & use at your own risk!\n\n";
559}
560
561
562// Used to determine the correspondence between the edges.
563// See faMeshReconstructor::calcAddressing
564void determineFaEdgeMapping
565(
566 const uindirectPrimitivePatch& onePatch,// reconstructed faMesh patch
567 const PtrList<faMesh>& procFaMeshes, // individual faMeshes
568 const labelListList& pointProcAddressing, // procPolyMesh to reconstructed
569
570 labelListList& faEdgeProcAddressing
571)
572{
573 // Determines ordering from faMesh to onePatch (= recsontructed faMesh)
574
575 // From two polyMesh points to patch-edge label
576 EdgeMap<label> pointsToOnePatchEdge(onePatch.nEdges());
577 {
578 const edgeList& edges = onePatch.edges();
579 const labelList& mp = onePatch.meshPoints();
580
581 forAll(edges, edgei)
582 {
583 const edge meshE(mp, edges[edgei]);
584 pointsToOnePatchEdge.insert(meshE, edgei);
585 }
586 }
587
588 // Now we can create the mapping from patch-edge on processor
589 // to patch-edge on onePatch
590
591 faEdgeProcAddressing.resize_nocopy(procFaMeshes.size());
592 forAll(procFaMeshes, proci)
593 {
594 const auto& procPatch = procFaMeshes[proci].patch();
595 const auto& edges = procPatch.edges();
596 const auto& mp = procPatch.meshPoints();
597 const auto& ppAddressing = pointProcAddressing[proci];
598
599 labelList& edgeProcAddr = faEdgeProcAddressing[proci];
600 edgeProcAddr.resize_nocopy(edges.size());
601
602 label edgei = 0;
603 for
604 (
605 ;
606 edgei < procPatch.nEdges(); //procPatch.nInternalEdges();
607 edgei++
608 )
609 {
610 const edge meshE(mp, edges[edgei]);
611 const edge onePatchE(ppAddressing, meshE);
612
613 edgeProcAddr[edgei] = pointsToOnePatchEdge[onePatchE];
614 }
615 }
616}
617
618
619void sortFaEdgeMapping
620(
621 const uindirectPrimitivePatch& onePatch,// reconstructed faMesh patch
622 const PtrList<faMesh>& procFaMeshes, // individual faMeshes
623 const labelListList& pointProcAddressing, // procPolyMesh to reconstructed
624
625 labelListList& faEdgeProcAddressing, // per proc the map to the master
626 labelListList& singlePatchEdgeLabels // per patch the map to the master
627)
628{
629 // From faMeshReconstructor.C - edge shuffling on patches
630
631 Map<label> remapGlobal;
632 remapGlobal.reserve(onePatch.nEdges());
633 for (label edgei = 0; edgei < onePatch.nInternalEdges(); ++edgei)
634 {
635 remapGlobal.insert(edgei, remapGlobal.size());
636 }
637
638
639 const faMesh& procMesh0 = procFaMeshes[0];
640 const label nNonProcessor = procMesh0.boundary().nNonProcessor();
641
642 singlePatchEdgeLabels.resize_nocopy(nNonProcessor);
643
644 forAll(singlePatchEdgeLabels, patchi)
645 {
646 if (isA<ignoreFaPatch>(procMesh0.boundary()[patchi]))
647 {
648 // These are not real edges
649 continue;
650 }
651
652 forAll(procFaMeshes, proci)
653 {
654 const faMesh& procMesh = procFaMeshes[proci];
655 const faPatch& fap = procMesh.boundary()[patchi];
656
657 labelList patchEdgeLabels(fap.edgeLabels());
658
659 // Renumber from local edges to global edges (natural order)
660 for (label& edgeId : patchEdgeLabels)
661 {
662 edgeId = faEdgeProcAddressing[proci][edgeId];
663 }
664
665 // Combine from all processors
666 singlePatchEdgeLabels[patchi].append(patchEdgeLabels);
667 }
668
669 // Sorted order will be the original non-decomposed order
670 Foam::sort(singlePatchEdgeLabels[patchi]);
671
672 for (const label sortedEdgei : singlePatchEdgeLabels[patchi])
673 {
674 remapGlobal.insert(sortedEdgei, remapGlobal.size());
675 }
676 }
677
678 {
679 // Use the map to rewrite the local faEdgeProcAddressing
680
681 labelListList newEdgeProcAddr(faEdgeProcAddressing);
682
683 forAll(procFaMeshes, proci)
684 {
685 const faMesh& procMesh = procFaMeshes[proci];
686
687 label edgei = procMesh.nInternalEdges();
688
689 for (const faPatch& fap : procMesh.boundary())
690 {
691 for (const label patchEdgei : fap.edgeLabels())
692 {
693 const label globalEdgei =
694 faEdgeProcAddressing[proci][patchEdgei];
695
696 const auto fnd = remapGlobal.cfind(globalEdgei);
697 if (fnd.good())
698 {
699 newEdgeProcAddr[proci][edgei] = fnd.val();
700 ++edgei;
701 }
702 else
703 {
705 << "Failed to find edge " << globalEdgei
706 << " this indicates a programming error" << nl
707 << exit(FatalError);
708 }
709 }
710 }
711 }
712 faEdgeProcAddressing = std::move(newEdgeProcAddr);
713 }
714}
715
716
717int main(int argc, char *argv[])
718{
720 (
721 "Reconstruct a mesh using geometric/topological information only"
722 );
723
724 // Enable -constant ... if someone really wants it
725 // Enable -withZero to prevent accidentally trashing the initial fields
726 timeSelector::addOptions(true, true); // constant(true), zero(true)
727
729
732 (
733 "addressing-only",
734 "Create procAddressing only without overwriting the mesh"
735 );
737 (
738 "mergeTol",
739 "scalar",
740 "The merge distance relative to the bounding box size (default 1e-7)"
741 );
743 (
744 "fullMatch",
745 "Do (slower) geometric matching on all boundary faces"
746 );
748 (
749 "procMatch",
750 "Do matching on processor faces only"
751 );
753 (
754 "cellDist",
755 "Write cell distribution as a labelList - for use with 'manual' "
756 "decomposition method or as a volScalarField for post-processing."
757 );
759 (
760 "no-finite-area",
761 "Suppress finiteArea mesh reconstruction",
762 true // Advanced option
763 );
764
765 #include "addAllRegionOptions.H"
766 #include "addAllFaRegionOptions.H"
767
768 // Prevent volume fields [with regionFaModels] from incidental
769 // triggering finite-area
771
772 // Prevent volume BCs from triggering finite-area
774
775 #include "setRootCase.H"
776 #include "createTime.H"
777
778 printWarning();
779
780 // ------------------------------------------------------------------------
781 // Configuration
782
783 const bool fullMatch = args.found("fullMatch");
784 const bool procMatch = args.found("procMatch");
785 const bool writeCellDist = args.found("cellDist");
786 bool doFiniteArea = !args.found("no-finite-area");
787 const bool writeAddrOnly = args.found("addressing-only");
788
789 const scalar mergeTol =
790 args.getOrDefault<scalar>("mergeTol", defaultMergeTol);
791
792 if (fullMatch)
793 {
794 Info<< "Use geometric matching on all boundary faces." << nl << endl;
795 }
796 else if (procMatch)
797 {
798 Info<< "Use geometric matching on correct procBoundaries only." << nl
799 << "This assumes a correct decomposition." << endl;
800 }
801 else
802 {
803 Info<< "Merge assuming correct, fully matched procBoundaries." << nl
804 << endl;
805 }
806
807 if (fullMatch || procMatch)
808 {
809 const scalar writeTol =
810 std::pow(scalar(10), -scalar(IOstream::defaultPrecision()));
811
812 Info<< "Merge tolerance : " << mergeTol << nl
813 << "Write tolerance : " << writeTol << endl;
814
815 if
816 (
817 runTime.writeFormat() == IOstreamOption::ASCII
818 && mergeTol < writeTol
819 )
820 {
822 << "Your current settings specify ASCII writing with "
823 << IOstream::defaultPrecision() << " digits precision." << endl
824 << "Your merging tolerance (" << mergeTol << ")"
825 << " is finer than this." << endl
826 << "Please change your writeFormat to binary"
827 << " or increase the writePrecision" << endl
828 << "or adjust the merge tolerance (-mergeTol)."
829 << exit(FatalError);
830 }
831 }
832
833 // Handle volume region selections
834 #include "getAllRegionOptions.H"
835
836 // Handle area region selections
837 #include "getAllFaRegionOptions.H"
838
839 if (!doFiniteArea)
840 {
841 areaRegionNames.clear(); // For consistency
842 }
843
844 // Determine the processor count
845 label nProcs{0};
846
847 if (regionNames.empty())
848 {
850 << "No regions specified or detected."
851 << exit(FatalError);
852 }
854 {
855 nProcs = fileHandler().nProcs(args.path());
856 }
857 else
858 {
859 nProcs = fileHandler().nProcs(args.path(), regionNames[0]);
860
861 if (regionNames.size() == 1)
862 {
863 Info<< "Using region: " << regionNames[0] << nl << endl;
864 }
865 }
866
867 if (!nProcs)
868 {
870 << "No processor* directories found"
871 << exit(FatalError);
872 }
873
874 // Read all time databases
875 PtrList<Time> databases(nProcs);
876
877 Info<< "Found " << nProcs << " processor directories" << endl;
878 forAll(databases, proci)
879 {
880 Info<< " Reading database "
881 << args.caseName()/("processor" + Foam::name(proci))
882 << endl;
883
884 databases.set
885 (
886 proci,
887 new Time
888 (
890 args.rootPath(),
891 args.caseName()/("processor" + Foam::name(proci)),
892 args.allowFunctionObjects(),
893 args.allowLibs()
894 )
895 );
896 }
897 Info<< endl;
898
899 // Use the times list from the master processor
900 // and select a subset based on the command-line options
902 (
903 databases[0].times(),
904 args
905 );
906
907 // Loop over all times
908 forAll(timeDirs, timei)
909 {
910 // Set time for global database
911 runTime.setTime(timeDirs[timei], timei);
912
913 // Set time for all databases
914 forAll(databases, proci)
915 {
916 databases[proci].setTime(timeDirs[timei], timei);
917 }
918
919 Info<< "Time = " << runTime.timeName() << endl;
920
921 // Check for any mesh changes
922 label nMeshChanged = 0;
923 boolList hasRegionMesh(regionNames.size(), false);
924 forAll(regionNames, regioni)
925 {
926 const word& regionName = regionNames[regioni];
927
928 IOobject facesIO
929 (
930 "faces",
931 databases[0].timeName(),
933 databases[0],
937 );
938
939 // Problem: faceCompactIOList recognises both 'faceList' and
940 // 'faceCompactList' so we should be lenient when doing
941 // typeHeaderOk
942
943 const bool ok = facesIO.typeHeaderOk<faceCompactIOList>(false);
944
945 if (ok)
946 {
947 hasRegionMesh[regioni] = true;
948 ++nMeshChanged;
949 }
950 }
951
952 // Check for any mesh changes
953 if (!nMeshChanged)
954 {
955 Info<< "No meshes" << nl << endl;
956 continue;
957 }
958
959 Info<< endl;
960
961 forAll(regionNames, regioni)
962 {
963 const word& regionName = regionNames[regioni];
965
966 if (!hasRegionMesh[regioni])
967 {
968 Info<< "region=" << regionName << " (no mesh)" << nl << endl;
969 continue;
970 }
971
972 if (regionNames.size() > 1)
973 {
974 Info<< "region=" << regionName << nl;
975 }
976
977
978 // Addressing from processor to reconstructed case
979 labelListList cellProcAddressing(nProcs);
980 labelListList faceProcAddressing(nProcs);
981 labelListList pointProcAddressing(nProcs);
982 labelListList boundProcAddressing(nProcs);
983
984
985 // Internal faces on the final reconstructed mesh
986 label masterInternalFaces;
987
988 // Owner addressing on the final reconstructed mesh
989 labelList masterOwner;
990
991 if (procMatch)
992 {
993 // Read point on individual processors to determine
994 // merge tolerance
995 // (otherwise single cell domains might give problems)
996
997 const boundBox bb = procBounds(databases, regionDir);
998 const scalar mergeDist = mergeTol*bb.mag();
999
1000 Info<< "Overall mesh bounding box : " << bb << nl
1001 << "Relative tolerance : " << mergeTol << nl
1002 << "Absolute matching distance : " << mergeDist << nl
1003 << endl;
1004
1005
1006 // Construct empty mesh.
1007 PtrList<fvMesh> masterMesh(nProcs);
1008
1009 for (label proci=0; proci<nProcs; proci++)
1010 {
1011 masterMesh.set
1012 (
1013 proci,
1014 new fvMesh
1015 (
1016 IOobject
1017 (
1018 regionName,
1019 runTime.timeName(),
1020 runTime,
1022 ),
1023 Zero
1024 )
1025 );
1026
1027 fvMesh meshToAdd
1028 (
1029 IOobject
1030 (
1031 regionName,
1032 databases[proci].timeName(),
1033 databases[proci]
1034 )
1035 );
1036
1037 // Initialize its addressing
1038 cellProcAddressing[proci] = identity(meshToAdd.nCells());
1039 faceProcAddressing[proci] = identity(meshToAdd.nFaces());
1040 pointProcAddressing[proci] = identity(meshToAdd.nPoints());
1041 boundProcAddressing[proci] =
1042 identity(meshToAdd.boundaryMesh().size());
1043
1044 // Find geometrically shared points/faces.
1045 autoPtr<faceCoupleInfo> couples = determineCoupledFaces
1046 (
1047 fullMatch,
1048 proci,
1049 proci,
1050 masterMesh[proci],
1051 proci,
1052 proci,
1053 meshToAdd,
1054 mergeDist
1055 );
1056
1057 // Add elements to mesh
1059 (
1060 masterMesh[proci],
1061 meshToAdd,
1062 couples()
1063 );
1064
1065 // Added processor
1066 renumber(map().addedCellMap(), cellProcAddressing[proci]);
1067 renumber(map().addedFaceMap(), faceProcAddressing[proci]);
1068 renumber(map().addedPointMap(), pointProcAddressing[proci]);
1069 renumber(map().addedPatchMap(), boundProcAddressing[proci]);
1070 }
1071 for (label step=2; step<nProcs*2; step*=2)
1072 {
1073 for (label proci=0; proci<nProcs; proci+=step)
1074 {
1075 label next = proci + step/2;
1076 if(next >= nProcs)
1077 {
1078 continue;
1079 }
1080
1081 Info<< "Merging mesh " << proci << " with "
1082 << next << endl;
1083
1084 // Find geometrically shared points/faces.
1085 autoPtr<faceCoupleInfo> couples = determineCoupledFaces
1086 (
1087 fullMatch,
1088 proci,
1089 next,
1090 masterMesh[proci],
1091 next,
1092 proci+step,
1093 masterMesh[next],
1094 mergeDist
1095 );
1096
1097 // Add elements to mesh
1099 (
1100 masterMesh[proci],
1101 masterMesh[next],
1102 couples()
1103 );
1104
1105 // Processors that were already in masterMesh
1106 for (label mergedI=proci; mergedI<next; mergedI++)
1107 {
1108 renumber
1109 (
1110 map().oldCellMap(),
1111 cellProcAddressing[mergedI]
1112 );
1113
1114 renumber
1115 (
1116 map().oldFaceMap(),
1117 faceProcAddressing[mergedI]
1118 );
1119
1120 renumber
1121 (
1122 map().oldPointMap(),
1123 pointProcAddressing[mergedI]
1124 );
1125
1126 // Note: boundary is special since can contain -1.
1127 renumber
1128 (
1129 map().oldPatchMap(),
1130 boundProcAddressing[mergedI]
1131 );
1132 }
1133
1134 // Added processor
1135 for
1136 (
1137 label addedI=next;
1138 addedI < Foam::min(proci+step, nProcs);
1139 addedI++
1140 )
1141 {
1142 renumber
1143 (
1144 map().addedCellMap(),
1145 cellProcAddressing[addedI]
1146 );
1147
1148 renumber
1149 (
1150 map().addedFaceMap(),
1151 faceProcAddressing[addedI]
1152 );
1153
1154 renumber
1155 (
1156 map().addedPointMap(),
1157 pointProcAddressing[addedI]
1158 );
1159
1160 renumber
1161 (
1162 map().addedPatchMap(),
1163 boundProcAddressing[addedI]
1164 );
1165 }
1166
1167 masterMesh.set(next, nullptr);
1168 }
1169 }
1170
1171 for (label proci=0; proci<nProcs; proci++)
1172 {
1173 Info<< "Reading mesh to add from "
1174 << databases[proci].caseName()
1175 << " for time = " << databases[proci].timeName()
1176 << nl << nl << endl;
1177 }
1178
1179 // See if any points on the mastermesh have become connected
1180 // because of connections through processor meshes.
1181 mergeSharedPoints(mergeDist,masterMesh[0],pointProcAddressing);
1182
1183 // Save some properties on the reconstructed mesh
1184 masterInternalFaces = masterMesh[0].nInternalFaces();
1185 masterOwner = masterMesh[0].faceOwner();
1186
1187
1188 if (writeAddrOnly)
1189 {
1190 Info<< nl
1191 << "Disabled writing of merged mesh (-addressing-only)"
1192 << nl << nl;
1193 }
1194 else
1195 {
1196 // Write reconstructed mesh
1197 writeMesh(masterMesh[0], cellProcAddressing);
1198 if (writeCellDist)
1199 {
1200 writeDistribution
1201 (
1202 runTime,
1203 masterMesh[0],
1204 cellProcAddressing
1205 );
1206 }
1207 }
1208 }
1209 else
1210 {
1211 // Load all meshes
1212 PtrList<fvMesh> fvMeshes(nProcs);
1213 for (label proci=0; proci<nProcs; proci++)
1214 {
1215 fvMeshes.set
1216 (
1217 proci,
1218 new fvMesh
1219 (
1220 IOobject
1221 (
1222 regionName,
1223 databases[proci].timeName(),
1224 databases[proci]
1225 )
1226 )
1227 );
1228 }
1229
1230 // Construct pointers to meshes
1231 UPtrList<polyMesh> meshes(fvMeshes.size());
1232 forAll(fvMeshes, proci)
1233 {
1234 meshes.set(proci, &fvMeshes[proci]);
1235 }
1236
1237 // Get pairs of patches to stitch. These pairs have to
1238 // - have ordered, opposite faces (so one to one correspondence)
1239 List<DynamicList<label>> localPatch;
1240 List<DynamicList<label>> remoteProc;
1241 List<DynamicList<label>> remotePatch;
1242 const label nGlobalPatches = polyMeshAdder::procPatchPairs
1243 (
1244 meshes,
1245 localPatch,
1246 remoteProc,
1247 remotePatch
1248 );
1249
1250 // Collect matching boundary faces on patches-to-stitch
1251 labelListList localBoundaryFace;
1252 labelListList remoteFaceProc;
1253 labelListList remoteBoundaryFace;
1255 (
1256 meshes,
1257 localPatch,
1258 remoteProc,
1259 remotePatch,
1260 localBoundaryFace,
1261 remoteFaceProc,
1262 remoteBoundaryFace
1263 );
1264
1265 // All matched faces assumed to have vertex0 matched
1266 labelListList remoteFaceStart(meshes.size());
1267 forAll(meshes, proci)
1268 {
1269 const labelList& procFaces = localBoundaryFace[proci];
1270 remoteFaceStart[proci].setSize(procFaces.size(), 0);
1271 }
1272
1273
1274
1275 labelListList patchMap(meshes.size());
1276 labelListList pointZoneMap(meshes.size());
1277 labelListList faceZoneMap(meshes.size());
1278 labelListList cellZoneMap(meshes.size());
1279 forAll(meshes, proci)
1280 {
1281 const polyMesh& mesh = meshes[proci];
1282 patchMap[proci] = identity(mesh.boundaryMesh().size());
1283
1284 // Remove excess patches
1285 patchMap[proci].setSize(nGlobalPatches);
1286
1287 pointZoneMap[proci] = identity(mesh.pointZones().size());
1288 faceZoneMap[proci] = identity(mesh.faceZones().size());
1289 cellZoneMap[proci] = identity(mesh.cellZones().size());
1290 }
1291
1292
1293 refPtr<fvMesh> masterMeshPtr;
1294 {
1295 // Do in-place addition on proc0.
1296
1297 const labelList oldFaceOwner(fvMeshes[0].faceOwner());
1298
1300 (
1301 0, // index of mesh to modify (== mesh_)
1302 fvMeshes,
1303 oldFaceOwner,
1304
1305 // Coupling info
1306 localBoundaryFace,
1307 remoteFaceProc,
1308 remoteBoundaryFace,
1309
1310 boundProcAddressing,
1311 cellProcAddressing,
1312 faceProcAddressing,
1313 pointProcAddressing
1314 );
1315
1316 // Remove zero-faces processor patches
1317 const polyBoundaryMesh& pbm = fvMeshes[0].boundaryMesh();
1318 labelList oldToNew(pbm.size(), -1);
1319 label newi = 0;
1320 // Non processor patches first
1321 forAll(pbm, patchi)
1322 {
1323 const auto& pp = pbm[patchi];
1324 if (!isA<processorPolyPatch>(pp) || pp.size())
1325 {
1326 oldToNew[patchi] = newi++;
1327 }
1328 }
1329 const label nNonProcPatches = newi;
1330
1331 // Move all deletable patches to the end
1332 forAll(oldToNew, patchi)
1333 {
1334 if (oldToNew[patchi] == -1)
1335 {
1336 oldToNew[patchi] = newi++;
1337 }
1338 }
1340 (
1341 fvMeshes[0],
1342 oldToNew,
1343 nNonProcPatches,
1344 false
1345 );
1346
1347 masterMeshPtr.cref(fvMeshes[0]);
1348 }
1349
1350
1351 const fvMesh& masterMesh = masterMeshPtr();
1352
1353 // Number of internal faces on the final reconstructed mesh
1354 masterInternalFaces = masterMesh.nInternalFaces();
1355
1356 // Owner addressing on the final reconstructed mesh
1357 masterOwner = masterMesh.faceOwner();
1358
1359
1360 // Write reconstructed mesh
1361 // Override:
1362 // - caseName
1363 // - processorCase flag
1364 // so the resulting mesh goes to the correct location (even with
1365 // collated). The better way of solving this is to construct
1366 // (zero) mesh on the undecomposed runTime.
1367
1368 if (writeAddrOnly)
1369 {
1370 Info<< nl
1371 << "Disabled writing of merged mesh (-addressing-only)"
1372 << nl << nl;
1373 }
1374 else
1375 {
1376 Time& masterTime = const_cast<Time&>(masterMesh.time());
1377
1378 const word oldCaseName = masterTime.caseName();
1379 masterTime.caseName() = runTime.caseName();
1380 const bool oldProcCase(masterTime.processorCase(false));
1381
1382 // Write reconstructed mesh
1383 writeMesh(masterMesh, cellProcAddressing);
1384 if (writeCellDist)
1385 {
1386 writeDistribution
1387 (
1388 runTime,
1389 masterMesh,
1390 cellProcAddressing
1391 );
1392 }
1393 masterTime.caseName() = oldCaseName;
1394 masterTime.processorCase(oldProcCase);
1395 }
1396 }
1397
1398
1399 // Write the addressing
1400
1401 Info<< "Reconstructing addressing from processor meshes"
1402 << " to the newly reconstructed mesh" << nl << endl;
1403
1404
1405 // Read finite-area
1406 // For each named area region that exists create a HashTable
1407 // entry that will contain the PtrList for all processors
1408
1409 PtrList<polyMesh> procMeshes(databases.size());
1410 HashTable<PtrList<faMesh>> procAreaRegionMeshes;
1411 procAreaRegionMeshes.reserve(areaRegionNames.size());
1412
1413 forAll(databases, proci)
1414 {
1415 Info<< "Processor " << proci << nl
1416 << "Read processor mesh: "
1417 << (databases[proci].caseName()/regionDir) << endl;
1418
1419 const polyMesh& procMesh = procMeshes.emplace
1420 (
1421 proci,
1422 IOobject
1423 (
1424 regionName,
1425 databases[proci].timeName(),
1426 databases[proci]
1427 )
1428 );
1429
1430 writeMaps
1431 (
1432 masterInternalFaces,
1433 masterOwner,
1434 procMesh,
1435 cellProcAddressing[proci],
1436 faceProcAddressing[proci],
1437 pointProcAddressing[proci],
1438 boundProcAddressing[proci]
1439 );
1440
1441
1442 if (doFiniteArea)
1443 {
1444 const word boundaryInst =
1445 procMesh.time().findInstance
1446 (
1447 procMesh.meshDir(),
1448 "boundary"
1449 );
1450
1451 for (const word& areaName : areaRegionNames)
1452 {
1453 IOobject io
1454 (
1455 "faBoundary",
1456 boundaryInst,
1457 faMesh::meshDir(procMesh, areaName),
1458 procMesh.time(),
1462 );
1463
1464 if (io.typeHeaderOk<faBoundaryMesh>(true))
1465 {
1466 // Always based on the volume decomposition!
1467 auto& procFaMeshes = procAreaRegionMeshes(areaName);
1468 procFaMeshes.resize(databases.size());
1469
1470 procFaMeshes.set
1471 (
1472 proci,
1473 new faMesh(areaName, procMesh)
1474 );
1475 }
1476 }
1477 }
1478 }
1479
1480
1481 // A finite-area mapping exists if procFaMeshes was filled
1482
1483 // Re-read reconstructed polyMesh. Note: could probably be avoided
1484 // by merging into loops above.
1485 refPtr<polyMesh> masterPolyMeshPtr;
1486
1487 if (!procAreaRegionMeshes.empty())
1488 {
1489 masterPolyMeshPtr.reset
1490 (
1491 new polyMesh
1492 (
1493 IOobject
1494 (
1495 regionName,
1496 runTime.timeName(),
1497 runTime
1498 ),
1499 true
1500 )
1501 );
1502 }
1503
1504 // Process any finite-area meshes
1505 for (const auto& iter : procAreaRegionMeshes.csorted())
1506 {
1507 const auto& areaName = iter.key();
1508 const auto& procFaMeshes = iter.val();
1509
1510 const polyMesh& masterMesh = masterPolyMeshPtr();
1511
1512 // Addressing from processor to reconstructed case
1513 labelListList faFaceProcAddressing(nProcs);
1514 labelListList faEdgeProcAddressing(nProcs);
1515 labelListList faPointProcAddressing(nProcs);
1516 labelListList faBoundProcAddressing(nProcs);
1517
1518
1519 // boundProcAddressing
1520 // ~~~~~~~~~~~~~~~~~~~
1521
1522 forAll(procFaMeshes, proci)
1523 {
1524 const auto& procMesh = procFaMeshes[proci];
1525 const auto& bm = procMesh.boundary();
1526
1527 faBoundProcAddressing[proci] = identity(bm.size());
1528 // Mark processor patches
1529 faBoundProcAddressing[proci].slice(bm.nNonProcessor()) = -1;
1530 }
1531
1532
1533 // faceProcAddressing
1534 // ~~~~~~~~~~~~~~~~~~
1535
1536 DynamicList<label> masterFaceLabels;
1537 label nPatchFaces = 0;
1538 forAll(procFaMeshes, proci)
1539 {
1540 const auto& procMesh = procFaMeshes[proci];
1541 const auto& procPolyFaces = procMesh.faceLabels();
1542 const auto& fpa = faceProcAddressing[proci];
1543
1544 labelList& faceAddr = faFaceProcAddressing[proci];
1545 faceAddr.resize_nocopy(procPolyFaces.size());
1546
1547 // Map to masterPolyMesh faces
1548 forAll(procPolyFaces, i)
1549 {
1550 const label facei = procPolyFaces[i];
1551 masterFaceLabels.append(fpa[facei]);
1552 faceAddr[i] = nPatchFaces++;
1553 }
1554 }
1555
1556
1557 // faMesh itself
1558 // ~~~~~~~~~~~~~
1559
1560 // Set up to read-if-present
1561 IOobject io(masterMesh);
1563
1564 // Construct without patches
1565 faMesh masterFaMesh
1566 (
1567 areaName,
1568 masterMesh,
1569 std::move(masterFaceLabels),
1570 io
1571 );
1572
1573 const auto& masterPatch = masterFaMesh.patch();
1574
1575
1576 // pointProcAddressing
1577 // ~~~~~~~~~~~~~~~~~~~
1578
1579 const auto& mpm = masterPatch.meshPointMap();
1580 forAll(procFaMeshes, proci)
1581 {
1582 const auto& procPatch = procFaMeshes[proci].patch();
1583 const auto& mp = procPatch.meshPoints();
1584
1585 auto& pointAddr = faPointProcAddressing[proci];
1586 pointAddr.resize_nocopy(mp.size());
1587
1588 forAll(mp, i)
1589 {
1590 pointAddr[i] = mpm[pointProcAddressing[proci][mp[i]]];
1591 }
1592 }
1593
1594
1595 // edgeProcAddressing
1596 // ~~~~~~~~~~~~~~~~~~
1597 // 1. straight mapping from proc faMesh back to reconstructed
1598 // faMesh
1599 determineFaEdgeMapping
1600 (
1601 masterPatch, // reconstructed faMesh patch
1602 procFaMeshes, // individual faMeshes
1603 pointProcAddressing, // procPolyMesh to reconstructed
1604
1605 faEdgeProcAddressing
1606 );
1607
1608
1609 // Per patch all edges on the masterPatch
1610 labelListList singlePatchEdgeLabels;
1611
1612 // 2. edgeProcAddressing : fix ordering on patches
1613 sortFaEdgeMapping
1614 (
1615 masterPatch, // reconstructed faMesh patch
1616 procFaMeshes, // individual faMeshes
1617 pointProcAddressing, // procPolyMesh to reconstructed
1618
1619 faEdgeProcAddressing, // per proc the map to the master
1620 singlePatchEdgeLabels // per patch the map to the master
1621 );
1622
1623
1624 const faMesh& procMesh0 = procFaMeshes[0];
1625
1626
1627 // Add faPatches
1628 // ~~~~~~~~~~~~~
1629
1630 faPatchList completePatches(singlePatchEdgeLabels.size());
1631 label nPatches = 0;
1632 forAll(completePatches, patchi)
1633 {
1634 const auto& patchEdgeLabels = singlePatchEdgeLabels[patchi];
1635
1636 const faPatch& fap = procMesh0.boundary()[patchi];
1637
1638 if (isA<ignoreFaPatch>(fap))
1639 {
1640 // These are not real edges
1641 continue;
1642 }
1643
1644 const label neiPolyPatchId = fap.ngbPolyPatchIndex();
1645
1646 completePatches.set
1647 (
1648 nPatches,
1649 fap.clone
1650 (
1651 masterFaMesh.boundary(),
1652 patchEdgeLabels,
1653 nPatches, // index
1654 neiPolyPatchId
1655 )
1656 );
1657
1658 ++nPatches;
1659 }
1660
1661 completePatches.resize(nPatches);
1662
1663 // Serial mesh - no parallel communication
1664
1665 const bool oldParRun = UPstream::parRun(false);
1666
1667 masterFaMesh.addFaPatches(completePatches);
1668
1669 UPstream::parRun(oldParRun); // Restore parallel state
1670
1671
1672 // Write mesh & individual addressing
1673 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1674
1676 (
1677 masterMesh.time().timeName(),
1678 masterFaMesh,
1679 masterFaMesh.faceLabels()
1680 );
1681
1682 forAll(procFaMeshes, proci)
1683 {
1684 const faMesh& procMesh = procFaMeshes[proci];
1685
1687 (
1688 IOobject
1689 (
1690 "procAddressing",
1691 masterMesh.time().timeName(),
1693 procMesh.thisDb(),
1697 ),
1698 faBoundProcAddressing[proci],
1699 faFaceProcAddressing[proci],
1700 faPointProcAddressing[proci],
1701 faEdgeProcAddressing[proci]
1702 );
1703 }
1704 }
1705 }
1706 }
1707
1708 Info<< "End\n" << endl;
1709
1710 return 0;
1711}
1712
1713
1714// ************************************************************************* //
Required Classes.
Required Classes.
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
const polyBoundaryMesh & pbm
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition DynamicList.H:68
void append(const T &val)
Copy append an element to the end of this list.
Map from edge (expressed as its endpoints) to value. Hashing (and ==) on an edge is symmetric.
Definition edgeHashes.H:59
A HashTable similar to std::unordered_map.
Definition HashTable.H:124
bool set(const Key &key, const T &obj)
Copy assign a new entry, overwriting existing entries.
Definition HashTableI.H:174
bool empty() const noexcept
True if the hash table is empty.
Definition HashTable.H:353
const_iterator cfind(const Key &key) const
Find and return an const_iterator set at the hashed entry.
Definition HashTableI.H:113
void reserve(label numEntries)
Reserve space for at least the specified number of elements (not the number of buckets) and regenerat...
Definition HashTable.C:729
bool insert(const Key &key, const T &obj)
Copy insert a new entry, not overwriting existing entries.
Definition HashTableI.H:152
label size() const noexcept
The number of elements in table.
Definition HashTable.H:358
UPtrList< const node_type > csorted() const
Const access to the hash-table contents in sorted order (sorted by keys).
Definition HashTable.C:181
bool emplace(const Key &key, Args &&... args)
Emplace insert a new entry, not overwriting existing entries.
Definition HashTableI.H:129
static void writeContents(const IOobject &io, const UList< label > &content)
@ NO_REGISTER
Do not request registration (bool: false).
@ NO_READ
Nothing to be read.
@ READ_IF_PRESENT
Reading is optional [identical to LAZY_READ].
@ MUST_READ
Reading required.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
@ ASCII
"ascii" (normal default)
static unsigned int defaultPrecision() noexcept
Return the default precision.
Definition IOstream.H:437
void resize_nocopy(const label len)
Adjust allocated size of list without necessarily.
Definition ListI.H:171
void append(const T &val)
Append an element at the end of the list.
Definition List.H:497
void setSize(label n)
Alias for resize().
Definition List.H:536
A HashTable to objects of type <T> with a label key.
Definition Map.H:54
label nEdges() const
Number of edges in patch.
const edgeList & edges() const
Return list of edges, address into LOCAL point list.
label nInternalEdges() const
Number of internal edges.
const labelList & meshPoints() const
Return labelList of mesh points in patch.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition PtrList.H:67
const T * set(const label i) const
Return const pointer to element (can be nullptr), or nullptr for out-of-range access (ie,...
Definition PtrList.H:171
void resize(const label newLen)
Adjust size of PtrList.
Definition PtrList.C:124
bool processorCase() const noexcept
True if this is a processor case.
Definition TimePathsI.H:52
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition Time.H:75
word findInstance(const fileName &directory, const word &name=word::null, IOobjectOption::readOption rOpt=IOobjectOption::MUST_READ, const word &stopInstance=word::null, const bool constant_fallback=true) const
Return time instance (location) of directory containing the file name (eg, used in reading mesh data)...
Definition Time.C:725
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
Definition Time.C:714
static word controlDictName
The default control dictionary name (normally "controlDict").
Definition Time.H:267
const fileName & caseName() const noexcept
The case name.
Definition TimePathsI.H:78
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
static bool parRun(const bool on) noexcept
Set as parallel run on/off.
Definition UPstream.H:1669
A list of pointers to objects of type <T>, without allocation/deallocation management of the pointers...
Definition UPtrList.H:101
label size() const noexcept
The number of entries in the list.
Definition UPtrListI.H:106
static void addVerboseOption(const string &usage="", bool advanced=false)
Enable a 'verbose' bool option, with usage information.
Definition argList.C:535
static void addBoolOption(const word &optName, const string &usage="", bool advanced=false)
Add a bool option to validOptions with usage information.
Definition argList.C:389
static void noParallel()
Remove the parallel options.
Definition argList.C:599
static void addOption(const word &optName, const string &param="", const string &usage="", bool advanced=false)
Add an option to validOptions with usage information.
Definition argList.C:400
static void addNote(const string &note)
Add extra notes for the usage information.
Definition argList.C:477
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
static autoPtr< T > New(Args &&... args)
Construct autoPtr with forwarding arguments.
Definition autoPtr.H:178
A bounding box defined in terms of min/max extrema points.
Definition boundBox.H:71
void add(const boundBox &bb)
Extend to include the second box.
Definition boundBoxI.H:323
scalar mag() const
The magnitude/length of the bounding box diagonal.
Definition boundBoxI.H:198
An edge is a list of two vertex labels. This can correspond to a directed graph edge or an edge on a ...
Definition edge.H:62
Finite area boundary mesh, which is a faPatch list with registered IO, a reference to the associated ...
label nNonProcessor() const
The number of patches before the first processor patch.
static void writeMesh(const word &timeName, const faMesh &fullMesh, const labelUList &singlePatchFaceLabels)
Write mesh information.
static void writeAddressing(const IOobject &io, const labelUList &faBoundaryProcAddr, const labelUList &faFaceProcAddr, const labelUList &faPointProcAddr, const labelUList &faEdgeProcAddr)
Write proc addressing.
Finite area mesh (used for 2-D non-Euclidian finite area method) defined using a patch of faces on a ...
Definition faMesh.H:140
virtual const objectRegistry & thisDb() const
Reference to the mesh database.
Definition faMesh.H:1209
const faBoundaryMesh & boundary() const noexcept
Return constant reference to boundary mesh.
Definition faMeshI.H:24
fileName meshDir() const
Return the local mesh directory (dbDir()/meshSubDir).
Definition faMesh.C:1022
label nInternalEdges() const noexcept
Number of internal faces.
Definition faMeshI.H:42
static word meshSubDir
The mesh sub-directory name (usually "faMesh").
Definition faMesh.H:750
Finite area patch class. Used for 2-D non-Euclidian finite area method.
Definition faPatch.H:76
virtual autoPtr< faPatch > clone(const faBoundaryMesh &bm) const
Construct and return a clone, resetting the boundary mesh.
Definition faPatch.H:257
const labelList & edgeLabels() const noexcept
Return the list of edges.
Definition faPatch.H:335
label ngbPolyPatchIndex() const noexcept
The neighbour polyPatch index.
Definition faPatch.H:359
A class for handling file names.
Definition fileName.H:75
static autoPtr< mapAddedPolyMesh > add(fvMesh &mesh0, const fvMesh &mesh1, const faceCoupleInfo &coupleInfo, const bool validBoundary=true, const bool fullyMapped=false)
Inplace add mesh to fvMesh. Maps all stored fields. Returns map.
Definition fvMeshAdder.C:67
static void reorderPatches(fvMesh &, const labelList &oldToNew, const label nPatches, const bool validBoundary)
Reorder and remove trailing patches.
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
const Time & time() const
Return the top-level database.
Definition fvMesh.H:360
static const word & zeroGradientType() noexcept
The type name for zeroGradient patch fields.
const Time & time() const noexcept
Return time registry.
A polyBoundaryMesh is a polyPatch list with registered IO, a reference to the associated polyMesh,...
static void mergePoints(const polyMesh &, const Map< label > &pointToMaster, polyTopoChange &meshMod)
Helper: Merge points.
static void patchFacePairs(const UPtrList< polyMesh > &meshes, const List< DynamicList< label > > &localPatch, const List< DynamicList< label > > &remoteMesh, const List< DynamicList< label > > &remotePatch, labelListList &localBoundaryFace, labelListList &remoteFaceMesh, labelListList &remoteBoundaryFace)
Helper: expand list of coupled patches into pairs of coupled faces.
static Map< label > findSharedPoints(const polyMesh &, const scalar mergeTol)
Find topologically and geometrically shared points.
static label procPatchPairs(const UPtrList< polyMesh > &meshes, List< DynamicList< label > > &localPatch, List< DynamicList< label > > &remoteMesh, List< DynamicList< label > > &remotePatch)
Helper: find pairs of processor patches. Return number of non-processor patches.
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition polyMesh.H:609
static const word & regionName(const word &region)
The mesh region name or word::null if polyMesh::defaultRegion.
Definition polyMesh.C:796
const fileName & facesInstance() const
Return the current instance directory for faces.
Definition polyMesh.C:844
static word defaultRegion
Return the default region name.
Definition polyMesh.H:406
fileName meshDir() const
Return the local mesh directory (dbDir()/meshSubDir).
Definition polyMesh.C:826
virtual const labelList & faceOwner() const
Return face owner.
Definition polyMesh.C:1101
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh").
Definition polyMesh.H:411
A patch is a list of labels that address the faces in the global face list.
Definition polyPatch.H:73
Direct mesh changes based on v1.3 polyTopoChange syntax.
bool isInternalFace(const label faceIndex) const noexcept
Return true if given face label is internal to the mesh.
label nBoundaryFaces() const noexcept
Number of boundary faces (== nFaces - nInternalFaces).
label nInternalFaces() const noexcept
Number of internal faces.
label nPoints() const noexcept
Number of mesh points.
label nCells() const noexcept
Number of mesh cells.
label nFaces() const noexcept
Number of mesh faces.
static word newName(const label myProcNo, const label neighbProcNo)
Return the name of a processorPolyPatch ("procBoundary..") constructed from the pair of processor IDs...
A class for managing references or pointers (no reference counting).
Definition refPtr.H:54
const T & cref() const
Return const reference to the object or to the contents of a (non-null) managed pointer.
Definition refPtrI.H:216
void reset(T *p=nullptr) noexcept
Delete managed pointer and set to new given pointer.
Definition refPtrI.H:314
static void addOptions(const bool constant=true, const bool withZero=false)
Add timeSelector options to argList::validOptions.
instantList select(const instantList &times) const
Select a list of Time values that are within the ranges.
static void removeFiles(const polyMesh &)
Helper: remove all sets files from mesh instance.
Definition topoSet.C:693
A class for handling words, derived from Foam::string.
Definition word.H:66
static const word null
An empty word.
Definition word.H:84
dynamicFvMesh & mesh
engineTime & runTime
Info<< "Creating field kinetic energy K\n"<< endl;volScalarField K("K", 0.5 *magSqr(U));if(U.nOldTimes()){ volVectorField *Uold=&U.oldTime();volScalarField *Kold=&K.oldTime(); *Kold==0.5 *magSqr(*Uold);while(Uold->nOldTimes()) { Uold=&Uold-> oldTime()
Foam::word regionName(args.getOrDefault< word >("region", Foam::polyMesh::defaultRegion))
Foam::PtrList< Foam::fvMesh > meshes(regionNames.size())
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
const auto & io
const word & regionDir
wordList areaRegionNames
wordList regionNames
const pointField & points
word timeName
Definition getTimeIndex.H:3
const dimensionedScalar mp
Proton mass.
bool allowFaModels() noexcept
The enable/disable state for regionFaModel (default: true).
Namespace for OpenFOAM.
List< edge > edgeList
List of edge.
Definition edgeList.H:32
const dimensionSet dimless
Dimensionless.
List< labelList > labelListList
List of labelList.
Definition labelList.H:38
List< label > labelList
A List of labels.
Definition List.H:62
GeometricField< scalar, fvPatchField, volMesh > volScalarField
UIndirectList< label > labelUIndList
UIndirectList of labels.
IOList< label > labelIOList
IO for a List of label.
Definition labelIOList.H:32
refPtr< fileOperation > fileHandler(std::nullptr_t)
Delete current file handler - forwards to fileOperation::handler().
CompactIOList< face > faceCompactIOList
Compact IO for a List of face.
Definition faceIOList.H:35
messageStream Info
Information stream (stdout output on master, null elsewhere).
IntListType renumber(const labelUList &oldToNew, const IntListType &input)
Renumber the values within a list.
vectorIOField pointIOField
pointIOField is a vectorIOField.
List< instant > instantList
List of instants.
Definition instantList.H:41
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
static void writeMaps(Ostream &os, const word &key, const labelListList &maps)
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
bool returnReduceAnd(const bool value, const int communicator=UPstream::worldComm)
Perform logical (and) MPI Allreduce on a copy. Uses UPstream::reduceAnd.
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
labelList identity(const label len, label start=0)
Return an identity map of the given length with (map[i] == i), works like std::iota() but returning a...
void sort(UList< T > &list)
Sort the list.
Definition UList.C:283
errorManip< error > abort(error &err)
Definition errorManip.H:139
List< bool > boolList
A List of bools.
Definition List.H:60
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
PrimitivePatch< UIndirectList< face >, const pointField & > uindirectPrimitivePatch
A PrimitivePatch with UIndirectList for the faces, const reference for the point field.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
PtrList< faPatch > faPatchList
Store lists of faPatch as a PtrList.
Definition faPatch.H:64
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
UList< label > labelUList
A UList of labels.
Definition UList.H:75
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
label newPointi
label nPatches
Foam::argList args(argc, argv)
volScalarField & e
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299