Loading...
Searching...
No Matches
redistributePar.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-2017 OpenFOAM Foundation
9 Copyright (C) 2015-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 redistributePar
29
30Group
31 grpParallelUtilities
32
33Description
34 Redistributes existing decomposed mesh and fields according to the current
35 settings in the decomposeParDict file.
36
37 Must be run on maximum number of source and destination processors.
38 Balances mesh and writes new mesh to new time directory.
39
40 Can optionally run in decompose/reconstruct mode to decompose/reconstruct
41 mesh and fields.
42
43Usage
44 \b redistributePar [OPTION]
45
46 Options:
47 - \par -decompose
48 Remove any existing \a processor subdirectories and decomposes the
49 mesh. Equivalent to running without processor subdirectories.
50
51 - \par -reconstruct
52 Reconstruct mesh and fields (like reconstructParMesh+reconstructPar).
53
54 - \par -newTimes
55 (in combination with -reconstruct) reconstruct only new times.
56
57 - \par -dry-run
58 (not in combination with -reconstruct) Test without actually
59 decomposing.
60
61 - \par -cellDist
62 not in combination with -reconstruct) Write the cell distribution
63 as a labelList, for use with 'manual'
64 decomposition method and as a volScalarField for visualization.
65
66 - \par -region <regionName>
67 Distribute named region.
68
69 - \par -allRegions
70 Distribute all regions in regionProperties. Does not check for
71 existence of processor*.
72
73\*---------------------------------------------------------------------------*/
74
75#include "argList.H"
76#include "sigFpe.H"
77#include "Time.H"
78#include "fvMesh.H"
79#include "fvMeshTools.H"
80#include "fvMeshDistribute.H"
81#include "fieldsDistributor.H"
82#include "decompositionMethod.H"
83#include "decompositionModel.H"
84#include "timeSelector.H"
85#include "PstreamReduceOps.H"
86#include "volFields.H"
87#include "surfaceFields.H"
89#include "IOobjectList.H"
90#include "globalIndex.H"
91#include "loadOrCreateMesh.H"
93#include "topoSet.H"
94#include "regionProperties.H"
95
98#include "hexRef8Data.H"
99#include "meshRefinement.H"
100#include "pointFields.H"
101
102#include "faMeshSubset.H"
103#include "faMeshTools.H"
104#include "faMeshDistributor.H"
105#include "faMeshesRegistry.H"
107
109
110using namespace Foam;
111
112// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
113
114// Use -verbose -verbose to switch on debug info. TBD.
115int debug(::Foam::debug::debugSwitch("redistributePar", 0));
116#define InfoOrPout (::debug ? Pout : Info.stream())
117
118
119// Allocate a new file handler on valid processors only
120// retaining the original IO ranks if possible
122getNewHandler(const boolUList& useProc, const bool verbose = true)
123{
125 (
126 fileOperation::New(fileHandler(), useProc, verbose)
127 );
128
129 if (::debug && handler)
130 {
131 Pout<< "Allocated " << handler().info()
132 << " myProcNo:" << UPstream::myProcNo(handler().comm())
133 << " ptr:" << Foam::name(handler.get()) << endl;
134 }
135
136 return handler;
137}
138
139
140// Allocate a new file handler on valid processors only
141// retaining the original IO ranks if possible
142void newHandler(const boolUList& useProc, refPtr<fileOperation>& handler)
143{
144 if (!handler)
145 {
146 handler = getNewHandler(useProc);
147 }
148}
149
150
151void createTimeDirs(const fileName& path)
152{
153 // Get current set of local processor's time directories. Uses
154 // fileHandler
155 instantList localTimeDirs(Time::findTimes(path, "constant"));
156
157 instantList masterTimeDirs;
158 if (Pstream::master())
159 {
160 masterTimeDirs = localTimeDirs;
161 }
162 Pstream::broadcast(masterTimeDirs);
163
164 // Sync any cached times (e.g. masterUncollatedFileOperation::times_)
165 // since only master would have done the findTimes
166 for (const instant& t : masterTimeDirs)
167 {
168 if (!localTimeDirs.contains(t))
169 {
170 const fileName timePath(path/t.name());
171
172 //Pout<< "Time:" << t << nl
173 // << " raw :" << timePath << nl
174 // << endl;
175 // Bypass fileHandler
176 Foam::mkDir(timePath);
177 }
178 }
179
180 // Just to make sure remove all state and re-scan
181 fileHandler().flush();
182 (void)Time::findTimes(path, "constant");
183}
184
185
186void copyUniform
187(
188 refPtr<fileOperation>& readHandler,
189 refPtr<fileOperation>& writeHandler,
190
191 const bool reconstruct,
192 const bool decompose,
193
194 const word& readTimeName,
195 const fileName& readCaseName,
196
197 const objectRegistry& readDb,
198 const objectRegistry& writeDb
199)
200{
201 // 3 modes: reconstruct, decompose, redistribute
202
203 // In reconstruct mode (separate reconstructed mesh):
204 // - read using readDb + readHandler
205 // - write using writeDb + writeHandler
206
207 // In decompose mode (since re-using processor0 mesh):
208 // - read using readDb + readCaseName + readHandler
209 // - write using writeDb + writeHandler
210
211 // In redistribute mode:
212 // - read using readDb + readHandler
213 // - write using writeDb + writeHandler
214
215 fileName readPath;
216
217 if (readHandler)
218 {
219 const label oldComm = UPstream::commWorld(readHandler().comm());
220
221 Time& readTime = const_cast<Time&>(readDb.time());
222 bool oldProcCase = readTime.processorCase();
223 string oldCaseName;
224 if (decompose)
225 {
226 //Pout<< "***Setting caseName to " << readCaseName
227 // << " to read undecomposed uniform" << endl;
228 oldCaseName = readTime.caseName();
229 readTime.caseName() = readCaseName;
230 oldProcCase = readTime.processorCase(false);
231 }
232
233 // Detect uniform/ at original database + time
234 readPath = readHandler().dirPath
235 (
236 false, // local directory
237 IOobject("uniform", readTimeName, readDb),
238 false // do not search in time
239 );
240
241
242 UPstream::commWorld(oldComm);
243
244 if (decompose)
245 {
246 // Reset caseName on master
247 //Pout<< "***Restoring caseName to " << oldCaseName << endl;
248 readTime.caseName() = oldCaseName;
249 readTime.processorCase(oldProcCase);
250 }
251 }
253
254 if (!readPath.empty())
255 {
256 InfoOrPout
257 << "Detected additional non-decomposed files in "
258 << readPath << endl;
259
260 // readPath: searching is the same for all file handlers. Typical:
261 // <case>/0.1/uniform (parent dir, decompose mode)
262 // <case>/processor1/0.1/uniform (redistribute/reconstruct mode)
263 // <case>/processors2/0.1/uniform ,,
264 // writePath:
265 // uncollated : <case>/0.1/uniform (reconstruct mode). Should only
266 // be done by master
267 // uncollated : <case>/processorXXX/0.1/uniform. Should be done by all.
268 // collated : <case>/processors2/0.1/uniform. Should be done by
269 // local master only.
270
271 const IOobject writeIO
272 (
273 "uniform",
274 writeDb.time().timeName(),
275 writeDb
276 );
277
278 // Switch to writeHandler
279 if (writeHandler)
280 {
281 auto oldHandler = fileOperation::fileHandler(writeHandler);
282
283 // Check: fileHandler.comm() is size 1 for uncollated
284 const label writeComm = fileHandler().comm();
285
286 if (reconstruct)
287 {
288 const bool oldParRun = UPstream::parRun(false);
289 const label oldNumProcs(fileHandler().nProcs());
290 const fileName writePath
291 (
292 fileHandler().objectPath
293 (
294 writeIO,
296 )
297 );
298 fileHandler().cp(readPath, writePath);
299 const_cast<fileOperation&>(fileHandler()).nProcs(oldNumProcs);
300 UPstream::parRun(oldParRun);
301 }
302 else
303 {
304 const fileName writePath
305 (
306 fileHandler().objectPath
307 (
308 writeIO,
310 )
311 );
312
313 if (::debug)
314 {
315 Pout<< " readPath :" << readPath << endl;
316 Pout<< " writePath :" << writePath << endl;
317 }
318
319 fileHandler().broadcastCopy
320 (
321 writeComm, // send to all in writeComm
322 UPstream::master(writeComm), // to use ioranks. Check!
323 readPath,
324 writePath
325 );
326 }
327 writeHandler = fileOperation::fileHandler(oldHandler);
328 }
329 }
330}
331
332
333void printMeshData(const polyMesh& mesh)
334{
335 // Collect all data on master
336
337 labelListList patchNeiProcNo(Pstream::nProcs());
338 labelListList patchSize(Pstream::nProcs());
339 const labelList& pPatches = mesh.globalData().processorPatches();
340 patchNeiProcNo[Pstream::myProcNo()].setSize(pPatches.size());
341 patchSize[Pstream::myProcNo()].setSize(pPatches.size());
342 forAll(pPatches, i)
343 {
345 (
346 mesh.boundaryMesh()[pPatches[i]]
347 );
348 patchNeiProcNo[Pstream::myProcNo()][i] = ppp.neighbProcNo();
349 patchSize[Pstream::myProcNo()][i] = ppp.size();
350 }
351 Pstream::gatherList(patchNeiProcNo);
352 Pstream::gatherList(patchSize);
353
354
355 // Print stats
356
357 const globalIndex globalCells(mesh.nCells());
358 const globalIndex globalBoundaryFaces(mesh.nBoundaryFaces());
359
360 label maxProcCells = 0;
361 label maxProcFaces = 0;
362 label totProcFaces = 0;
363 label maxProcPatches = 0;
364 label totProcPatches = 0;
365
366 for (const int proci : Pstream::allProcs())
367 {
368 const label nLocalCells = globalCells.localSize(proci);
369 const label nBndFaces = globalBoundaryFaces.localSize(proci);
370
371 InfoOrPout<< nl
372 << "Processor " << proci;
373
374 if (!nLocalCells)
375 {
376 InfoOrPout<< " (empty)" << endl;
377 continue;
378 }
379 else
380 {
381 InfoOrPout<< nl
382 << " Number of cells = " << nLocalCells << endl;
383 }
384
385 label nProcFaces = 0;
386 const labelList& nei = patchNeiProcNo[proci];
387
388 forAll(patchNeiProcNo[proci], i)
389 {
390 InfoOrPout
391 << " Number of faces shared with processor "
392 << patchNeiProcNo[proci][i] << " = "
393 << patchSize[proci][i] << nl;
394
395 nProcFaces += patchSize[proci][i];
396 }
397
398 {
399 InfoOrPout
400 << " Number of processor patches = " << nei.size() << nl
401 << " Number of processor faces = " << nProcFaces << nl
402 << " Number of boundary faces = "
403 << nBndFaces-nProcFaces << endl;
404 }
405
406 maxProcCells = Foam::max(maxProcCells, nLocalCells);
407 totProcFaces += nProcFaces;
408 totProcPatches += nei.size();
409 maxProcFaces = Foam::max(maxProcFaces, nProcFaces);
410 maxProcPatches = Foam::max(maxProcPatches, nei.size());
411 }
412
413 // Summary stats
414
415 InfoOrPout
416 << nl
417 << "Number of processor faces = " << (totProcFaces/2) << nl
418 << "Max number of cells = " << maxProcCells;
419
420 if (maxProcCells != globalCells.totalSize())
421 {
422 scalar avgValue = scalar(globalCells.totalSize())/Pstream::nProcs();
423
424 InfoOrPout
425 << " (" << 100.0*(maxProcCells-avgValue)/avgValue
426 << "% above average " << avgValue << ')';
427 }
428 InfoOrPout<< nl;
429
430 InfoOrPout<< "Max number of processor patches = " << maxProcPatches;
431 if (totProcPatches)
432 {
433 scalar avgValue = scalar(totProcPatches)/Pstream::nProcs();
434
435 InfoOrPout
436 << " (" << 100.0*(maxProcPatches-avgValue)/avgValue
437 << "% above average " << avgValue << ')';
438 }
439 InfoOrPout<< nl;
440
441 InfoOrPout<< "Max number of faces between processors = " << maxProcFaces;
442 if (totProcFaces)
443 {
444 scalar avgValue = scalar(totProcFaces)/Pstream::nProcs();
445
446 InfoOrPout
447 << " (" << 100.0*(maxProcFaces-avgValue)/avgValue
448 << "% above average " << avgValue << ')';
449 }
450 InfoOrPout<< nl << endl;
451}
452
453
454// Debugging: write volScalarField with decomposition for post processing.
455void writeDecomposition
456(
457 const word& name,
458 const fvMesh& mesh,
459 const labelUList& decomp
460)
461{
462 // Write the decomposition as labelList for use with 'manual'
463 // decomposition method.
465 (
466 // NB: mesh read from facesInstance
467 IOobject("cellDecomposition", mesh.facesInstance(), mesh),
468 decomp
469 );
470
471 InfoOrPout
472 << "Writing wanted cell distribution to volScalarField " << name
473 << " for postprocessing purposes." << nl << endl;
474
475 volScalarField procCells
476 (
478 (
479 name,
480 mesh.time().timeName(),
481 mesh,
485 ),
486 mesh,
487 dimensionedScalar("", dimless, -1),
489 );
490
491 forAll(procCells, celli)
492 {
493 procCells[celli] = decomp[celli];
494 }
495
496 procCells.correctBoundaryConditions();
497 procCells.write();
498}
499
500
501void determineDecomposition
502(
503 refPtr<fileOperation>& readHandler,
504 const Time& baseRunTime,
505 const fileName& decompDictFile, // optional location for decomposeParDict
506 const bool decompose, // decompose, i.e. read from undecomposed case
507 const fileName& proc0CaseName,
508 const fvMesh& mesh,
509 const bool writeCellDist,
510
511 label& nDestProcs,
512 labelList& decomp
513)
514{
515 // Switch to readHandler since decomposition method might do IO
516 // (e.g. read decomposeParDict)
517 auto oldHandler = fileOperation::fileHandler(readHandler);
518
519 // Read decomposeParDict (on all processors)
521 (
522 mesh,
523 decompDictFile
524 );
525
526 decompositionMethod& decomposer = method.decomposer();
527
528 if (!decomposer.parallelAware())
529 {
531 << "You have selected decomposition method \""
532 << decomposer.type() << "\n"
533 << " which does not synchronise decomposition across"
534 " processor patches.\n"
535 " You might want to select a decomposition method"
536 " that is aware of this. Continuing...." << endl;
537 }
538
539 Time& tm = const_cast<Time&>(mesh.time());
540
541 const bool oldProcCase = tm.processorCase();
542 if (decompose)
543 {
544 InfoOrPout
545 << "Setting caseName to " << baseRunTime.caseName()
546 << " to read decomposeParDict" << endl;
547 tm.caseName() = baseRunTime.caseName();
548 tm.processorCase(false);
549 }
550
551 scalarField cellWeights;
552 if (method.found("weightField"))
553 {
554 word weightName = method.get<word>("weightField");
555
556 volScalarField weights
557 (
559 (
560 weightName,
561 tm.timeName(),
562 mesh,
566 ),
567 mesh
568 );
569 cellWeights = weights.internalField();
570 }
571
572 nDestProcs = decomposer.nDomains();
573 decomp = decomposer.decompose(mesh, cellWeights);
574
575 readHandler = fileOperation::fileHandler(oldHandler);
576
577
578 if (decompose)
579 {
580 InfoOrPout
581 << "Restoring caseName" << endl;
582 tm.caseName() = proc0CaseName;
583 tm.processorCase(oldProcCase);
584 }
585
586 // Dump decomposition to volScalarField
587 if (writeCellDist)
588 {
589 const label oldNumProcs =
590 const_cast<fileOperation&>(fileHandler()).nProcs(nDestProcs);
591
592 // Note: on master make sure to write to processor0
593 if (decompose)
594 {
595 if (UPstream::master())
596 {
597 const bool oldParRun = UPstream::parRun(false);
598
599 InfoOrPout
600 << "Setting caseName to " << baseRunTime.caseName()
601 << " to write undecomposed cellDist" << endl;
602
603 tm.caseName() = baseRunTime.caseName();
604 tm.processorCase(false);
605 writeDecomposition("cellDist", mesh, decomp);
606 InfoOrPout<< "Restoring caseName" << endl;
607 tm.caseName() = proc0CaseName;
608 tm.processorCase(oldProcCase);
609
610 UPstream::parRun(oldParRun);
611 }
612 }
613 else
614 {
615 writeDecomposition("cellDist", mesh, decomp);
616 }
617
618 const_cast<fileOperation&>(fileHandler()).nProcs(oldNumProcs);
619 }
620}
621
622
623// Variant of GeometricField::correctBoundaryConditions that only
624// evaluates selected patch fields
625template<class CoupledPatchType, class GeoField>
626void correctCoupledBoundaryConditions(fvMesh& mesh)
627{
628 for (GeoField& fld : mesh.sorted<GeoField>())
629 {
630 fld.boundaryFieldRef().template evaluateCoupled<CoupledPatchType>();
631 }
632}
633
634
635// Inplace redistribute mesh and any fields
636autoPtr<mapDistributePolyMesh> redistributeAndWrite
637(
638 refPtr<fileOperation>& readHandler,
639 refPtr<fileOperation>& writeHandler,
640 const Time& baseRunTime,
641 const fileName& proc0CaseName,
642
643 // Controls
644 const bool doReadFields,
645 const bool decompose, // decompose, i.e. read from undecomposed case
646 const bool reconstruct,
647 const bool overwrite,
648
649 // Decomposition information
650 const label nDestProcs,
651 const labelList& decomp,
652
653 // Mesh information
654 const boolList& volMeshOnProc,
655 const fileName& volMeshInstance,
656 fvMesh& mesh
657)
658{
659 Time& runTime = const_cast<Time&>(mesh.time());
660 const bool oldProcCase = runTime.processorCase();
661
663 //Pout<< "Before distribution:" << endl;
664 //printMeshData(mesh);
665
666
667 // Storage of fields
668
669 PtrList<volScalarField> volScalarFields;
670 PtrList<volVectorField> volVectorFields;
671 PtrList<volSphericalTensorField> volSphereTensorFields;
672 PtrList<volSymmTensorField> volSymmTensorFields;
673 PtrList<volTensorField> volTensorFields;
674
675 PtrList<surfaceScalarField> surfScalarFields;
676 PtrList<surfaceVectorField> surfVectorFields;
677 PtrList<surfaceSphericalTensorField> surfSphereTensorFields;
678 PtrList<surfaceSymmTensorField> surfSymmTensorFields;
679 PtrList<surfaceTensorField> surfTensorFields;
680
686
687 PtrList<pointScalarField> pointScalarFields;
688 PtrList<pointVectorField> pointVectorFields;
689 PtrList<pointTensorField> pointTensorFields;
690 PtrList<pointSphericalTensorField> pointSphTensorFields;
691 PtrList<pointSymmTensorField> pointSymmTensorFields;
692
693 // Self-contained pointMesh for reading pointFields
694 const pointMesh oldPointMesh(mesh);
695
696 // Track how many (if any) pointFields are read/mapped
697 label nPointFields = 0;
698
699 refPtr<fileOperation> noWriteHandler;
700
701 parPointFieldDistributor pointDistributor
702 (
703 oldPointMesh, // source mesh
704 false, // savePoints=false (ie, delay until later)
705 //false // Do not write
706 noWriteHandler // Do not write
707 );
708
709
710 if (doReadFields)
711 {
712 // Create 0 sized mesh to do all the generation of zero sized
713 // fields on processors that have zero sized meshes. Note that this is
714 // only necessary on master but since polyMesh construction with
715 // Pstream::parRun does parallel comms we have to do it on all
716 // processors
717 autoPtr<fvMeshSubset> subsetterPtr;
718
719 // Missing a volume mesh somewhere?
720 if (volMeshOnProc.found(false))
721 {
722 // A zero-sized mesh with boundaries.
723 // This is used to create zero-sized fields.
724 subsetterPtr.reset(new fvMeshSubset(mesh, zero{}));
725 subsetterPtr().subMesh().init(true);
726 subsetterPtr().subMesh().globalData();
727 subsetterPtr().subMesh().tetBasePtIs();
728 subsetterPtr().subMesh().geometricD();
729 }
730
731
732 // Get original objects (before incrementing time!)
733 if (decompose)
734 {
735 InfoOrPout
736 << "Setting caseName to " << baseRunTime.caseName()
737 << " to read IOobjects" << endl;
738 runTime.caseName() = baseRunTime.caseName();
739 runTime.processorCase(false);
740 }
741
742 //IOobjectList objects(mesh, runTime.timeName());
743 // Swap to reading fileHandler and read IOobjects
744 IOobjectList objects;
745 if (readHandler)
746 {
747 auto oldHandler = fileOperation::fileHandler(readHandler);
748 const label oldComm = UPstream::commWorld(fileHandler().comm());
749
750 objects = IOobjectList(mesh, runTime.timeName());
751 readHandler = fileOperation::fileHandler(oldHandler);
752 UPstream::commWorld(oldComm);
753 }
754
755 if (decompose)
756 {
757 InfoOrPout
758 << "Restoring caseName" << endl;
759 runTime.caseName() = proc0CaseName;
760 runTime.processorCase(oldProcCase);
761 }
762
763 InfoOrPout
764 << "From time " << runTime.timeName()
765 << " mesh:" << mesh.objectRegistry::objectRelPath()
766 << " have objects:" << objects.names() << endl;
767
768 // We don't want to map the decomposition (mapping already tested when
769 // mapping the cell centre field)
770 (void)objects.erase("cellDist");
771
772
773 if (decompose)
774 {
775 runTime.caseName() = baseRunTime.caseName();
776 runTime.processorCase(false);
777 }
778
779 // Field reading
780
781 #undef doFieldReading
782 #define doFieldReading(Storage) \
783 { \
784 fieldsDistributor::readFields \
785 ( \
786 volMeshOnProc, readHandler, mesh, subsetterPtr, \
787 objects, Storage \
788 ); \
789 }
790
791 // volField
792 doFieldReading(volScalarFields);
793 doFieldReading(volVectorFields);
794 doFieldReading(volSphereTensorFields);
795 doFieldReading(volSymmTensorFields);
796 doFieldReading(volTensorFields);
797
798 // surfaceField
799 doFieldReading(surfScalarFields);
800 doFieldReading(surfVectorFields);
801 doFieldReading(surfSphereTensorFields);
802 doFieldReading(surfSymmTensorFields);
803 doFieldReading(surfTensorFields);
804
805 // Dimensioned internal fields
806 doFieldReading(dimScalarFields);
807 doFieldReading(dimVectorFields);
808 doFieldReading(dimSphereTensorFields);
809 doFieldReading(dimSymmTensorFields);
810 doFieldReading(dimTensorFields);
811
812 // pointFields
813 nPointFields = 0;
814
815 #undef doFieldReading
816 #define doFieldReading(Storage) \
817 { \
818 fieldsDistributor::readFields \
819 ( \
820 volMeshOnProc, readHandler, oldPointMesh, \
821 subsetterPtr, objects, Storage, \
822 true /* (deregister field) */ \
823 ); \
824 nPointFields += Storage.size(); \
825 }
826
827 doFieldReading(pointScalarFields);
828 doFieldReading(pointVectorFields);
829 doFieldReading(pointSphTensorFields);
830 doFieldReading(pointSymmTensorFields);
831 doFieldReading(pointTensorFields);
832 #undef doFieldReading
833
834
835 // Done reading
836
837 if (decompose)
838 {
839 runTime.caseName() = proc0CaseName;
840 runTime.processorCase(oldProcCase);
841 }
842 }
843
844 // Save pointMesh information before any topology changes occur!
845 if (nPointFields)
846 {
847 pointDistributor.saveMeshPoints();
848 }
849
850
851 // Read refinement data
852 autoPtr<hexRef8Data> refDataPtr;
853 {
854 // Read refinement data
855 if (decompose)
856 {
857 runTime.caseName() = baseRunTime.caseName();
858 runTime.processorCase(false);
859 }
861 (
862 "dummy",
863 volMeshInstance, //mesh.facesInstance(),
865 mesh,
869 );
870
871 if (readHandler)
872 {
873 auto oldHandler = fileOperation::fileHandler(readHandler);
874 const label oldComm = UPstream::commWorld(fileHandler().comm());
875
876 // Read
877 refDataPtr.reset(new hexRef8Data(io));
878
879 UPstream::commWorld(oldComm);
880 readHandler = fileOperation::fileHandler(oldHandler);
881 }
882 else
883 {
885 refDataPtr.reset(new hexRef8Data(io));
886 }
887
888 if (decompose)
889 {
890 runTime.caseName() = proc0CaseName;
891 runTime.processorCase(oldProcCase);
892 }
893
894 // Make sure all processors have valid data (since only some will
895 // read)
896 refDataPtr().sync(io);
897
898 // Now we've read refinement data we can remove it
900 }
901
902
903 // If faMeshesRegistry exists, it is also owned by the polyMesh and will
904 // be destroyed by clearGeom() in fvMeshDistribute::distribute()
905 //
906 // Rescue faMeshesRegistry from destruction by temporarily moving
907 // it to be locally owned.
908 std::unique_ptr<faMeshesRegistry> faMeshesRegistry_saved
909 (
911 );
912
913 // Mesh distribution engine
914 fvMeshDistribute distributor(mesh);
915
916 // Do all the distribution of mesh and fields
917 autoPtr<mapDistributePolyMesh> distMap = distributor.distribute(decomp);
918
919 // Restore ownership onto the polyMesh
920 faMeshesRegistry::Store(std::move(faMeshesRegistry_saved));
921
922
923 // Print some statistics
924 InfoOrPout<< "After distribution:" << endl;
925 printMeshData(mesh);
926
927 // Get other side of processor boundaries
928 do
929 {
930 #undef doCorrectCoupled
931 #define doCorrectCoupled(FieldType) \
932 correctCoupledBoundaryConditions<processorFvPatch, FieldType>(mesh);
933
934 doCorrectCoupled(volScalarField);
935 doCorrectCoupled(volVectorField);
936 doCorrectCoupled(volSphericalTensorField);
937 doCorrectCoupled(volSymmTensorField);
938 doCorrectCoupled(volTensorField);
939 #undef doCorrectCoupled
940 }
941 while (false);
942
943 // No update surface fields
944
945
946 // Map pointFields
947 if (nPointFields)
948 {
949 // Construct new pointMesh from distributed mesh
950 const pointMesh& newPointMesh = pointMesh::New(mesh);
951
952 pointDistributor.resetTarget(newPointMesh, distMap());
953
954 pointDistributor.distributeAndStore(pointScalarFields);
955 pointDistributor.distributeAndStore(pointVectorFields);
956 pointDistributor.distributeAndStore(pointSphTensorFields);
957 pointDistributor.distributeAndStore(pointSymmTensorFields);
958 pointDistributor.distributeAndStore(pointTensorFields);
959 }
960
961
962 // More precision (for points data)
964
965
966 if (!overwrite)
967 {
968 ++runTime;
969 mesh.setInstance(runTime.timeName());
970 }
971 else
972 {
973 mesh.setInstance(volMeshInstance);
974 }
975
976
977 // Register mapDistributePolyMesh for automatic writing...
979 (
981 (
982 "procAddressing",
983 mesh.facesInstance(),
985 mesh.thisDb(),
989 ),
990 distMap()
991 );
992
993
994 if (reconstruct)
995 {
996 auto oldHandler = fileOperation::fileHandler(writeHandler);
997
998 if (UPstream::master())
999 {
1000 InfoOrPout
1001 << "Setting caseName to " << baseRunTime.caseName()
1002 << " to write reconstructed mesh (and fields)." << endl;
1003 runTime.caseName() = baseRunTime.caseName();
1004 const bool oldProcCase(runTime.processorCase(false));
1005 const label oldNumProcs
1006 (
1007 const_cast<fileOperation&>(fileHandler()).nProcs(nDestProcs)
1008 );
1009 const bool oldParRun = UPstream::parRun(false);
1010
1011 mesh.write();
1013
1014 // Now we've written all. Reset caseName on master
1015 InfoOrPout<< "Restoring caseName" << endl;
1016 UPstream::parRun(oldParRun);
1017 const_cast<fileOperation&>(fileHandler()).nProcs(oldNumProcs);
1018 runTime.caseName() = proc0CaseName;
1019 runTime.processorCase(oldProcCase);
1020 }
1021
1022 writeHandler = fileOperation::fileHandler(oldHandler);
1023 }
1024 else
1025 {
1026 auto oldHandler = fileOperation::fileHandler(writeHandler);
1027
1028 const label oldNumProcs
1029 (
1030 const_cast<fileOperation&>(fileHandler()).nProcs(nDestProcs)
1031 );
1032 mesh.write();
1033 const_cast<fileOperation&>(fileHandler()).nProcs(oldNumProcs);
1034
1035 writeHandler = fileOperation::fileHandler(oldHandler);
1036
1038 }
1039 InfoOrPout
1040 << "Written redistributed mesh to "
1041 << mesh.facesInstance() << nl << endl;
1042
1043
1044 if (decompose)
1045 {
1046 // Decompose (1 -> N)
1047 // so {boundary,cell,face,point}ProcAddressing have meaning
1049 (
1050 mesh,
1051 distMap(),
1052 decompose,
1053 mesh.facesInstance(), //oldFacesInstance,
1054 writeHandler // to write *ProcAddressing
1055 );
1056 }
1057 else if (reconstruct)
1058 {
1059 // Reconstruct (N -> 1)
1060 // so {boundary,cell,face,point}ProcAddressing have meaning. Make sure
1061 // to write these to meshes containing the source meshes (i.e. using
1062 // the read handler)
1064 (
1065 mesh,
1066 distMap(),
1067 decompose,
1068 volMeshInstance, //oldFacesInstance,
1069 readHandler //writeHandler
1070 );
1071 }
1072 else
1073 {
1074 // Redistribute (N -> M)
1075 // {boundary,cell,face,point}ProcAddressing would be incorrect
1076 // - can either remove or redistribute previous
1078 }
1079
1080
1081 // Refinement data
1082 if (refDataPtr)
1083 {
1084 auto& refData = refDataPtr();
1085
1086 // Set instance
1087 IOobject io
1088 (
1089 "dummy",
1090 mesh.facesInstance(),
1092 mesh,
1096 );
1097 refData.sync(io);
1098
1099
1100 // Distribute
1101 refData.distribute(distMap());
1102
1103
1104 auto oldHandler = fileOperation::fileHandler(writeHandler);
1105
1106 if (reconstruct)
1107 {
1108 if (UPstream::master())
1109 {
1110 const bool oldParRun = UPstream::parRun(false);
1111 const label oldNumProcs
1112 (
1113 const_cast<fileOperation&>(fileHandler()).nProcs(nDestProcs)
1114 );
1115
1116 InfoOrPout
1117 << "Setting caseName to " << baseRunTime.caseName()
1118 << " to write reconstructed refinement data." << endl;
1119 runTime.caseName() = baseRunTime.caseName();
1120 const bool oldProcCase(runTime.processorCase(false));
1121
1122 refData.write();
1123
1124 // Now we've written all. Reset caseName on master
1125 InfoOrPout<< "Restoring caseName" << endl;
1126 runTime.caseName() = proc0CaseName;
1127 runTime.processorCase(oldProcCase);
1128
1129 const_cast<fileOperation&>(fileHandler()).nProcs(oldNumProcs);
1130 UPstream::parRun(oldParRun);
1131 }
1132 }
1133 else
1134 {
1135 const label oldNumProcs
1136 (
1137 const_cast<fileOperation&>(fileHandler()).nProcs(nDestProcs)
1138 );
1139 refData.write();
1140 const_cast<fileOperation&>(fileHandler()).nProcs(oldNumProcs);
1141 }
1142
1143 writeHandler = fileOperation::fileHandler(oldHandler);
1144 }
1145
1147 //{
1148 // // Read sets
1149 // if (decompose)
1150 // {
1151 // runTime.caseName() = baseRunTime.caseName();
1152 // runTime.processorCase(false);
1153 // }
1154 // IOobjectList objects(mesh, mesh.facesInstance(), "polyMesh/sets");
1155 //
1156 // PtrList<cellSet> cellSets;
1157 // ReadFields(objects, cellSets);
1158 //
1159 // if (decompose)
1160 // {
1161 // runTime.caseName() = proc0CaseName;
1162 // runTime.processorCase(oldProcCase);
1163 // }
1164 //
1165 // forAll(cellSets, i)
1166 // {
1167 // cellSets[i].distribute(distMap());
1168 // }
1169 //
1170 // if (reconstruct)
1171 // {
1172 // if (Pstream::master())
1173 // {
1174 // Pout<< "Setting caseName to " << baseRunTime.caseName()
1175 // << " to write reconstructed refinement data." << endl;
1176 // runTime.caseName() = baseRunTime.caseName();
1177 // const bool oldProcCase(runTime.processorCase(false));
1178 //
1179 // forAll(cellSets, i)
1180 // {
1181 // cellSets[i].distribute(distMap());
1182 // }
1183 //
1184 // // Now we've written all. Reset caseName on master
1185 // Pout<< "Restoring caseName" << endl;
1186 // runTime.caseName() = proc0CaseName;
1187 // runTime.processorCase(oldProcCase);
1188 // }
1189 // }
1190 // else
1191 // {
1192 // forAll(cellSets, i)
1193 // {
1194 // cellSets[i].distribute(distMap());
1195 // }
1196 // }
1197 //}
1198
1199
1200 return distMap;
1201}
1202
1203
1204/*---------------------------------------------------------------------------*\
1205 main
1206\*---------------------------------------------------------------------------*/
1207
1208int main(int argc, char *argv[])
1209{
1211 (
1212 "Redistribute decomposed mesh and fields according"
1213 " to the decomposeParDict settings.\n"
1214 "Optionally run in decompose/reconstruct mode"
1215 );
1216
1217 argList::noFunctionObjects(); // Never use function objects
1218
1219 // enable -constant ... if someone really wants it
1220 // enable -zeroTime to prevent accidentally trashing the initial fields
1221 timeSelector::addOptions(true, true);
1222
1223 #include "addAllRegionOptions.H"
1224
1225 #include "addOverwriteOption.H"
1226 argList::addBoolOption("decompose", "Decompose case");
1227 argList::addBoolOption("reconstruct", "Reconstruct case");
1229 (
1230 "Additional verbosity. (Can be used multiple times for debugging)"
1231 );
1233 (
1234 "Test without writing the decomposition. "
1235 "Changes -cellDist to only write volScalarField."
1236 );
1237 argList::addVerboseOption("Additional verbosity");
1239 (
1240 "cellDist",
1241 "Write cell distribution as a labelList - for use with 'manual' "
1242 "decomposition method or as a volScalarField for post-processing."
1243 );
1245 (
1246 "newTimes",
1247 "Only reconstruct new times (i.e. that do not exist already)"
1248 );
1250 (
1251 "Additional verbosity. (Can be used multiple times)"
1252 );
1254 (
1255 "no-finite-area",
1256 "Suppress finiteArea mesh/field handling",
1257 true // Advanced option
1258 );
1259
1260 // Prevent volume BCs from triggering finite-area
1262
1263 //- Disable caching of times/processor dirs etc. Cause massive parallel
1264 // problems when e.g decomposing.
1266
1267
1268 // Handle arguments
1269 // ~~~~~~~~~~~~~~~~
1270 // (replacement for setRootCase that does not abort)
1271
1272 argList args(argc, argv);
1273
1274 // ------------------------------------------------------------------------
1275
1276 // As much as possible avoid synchronised operation. To be looked at more
1277 // closely for the three scenarios:
1278 // - decompose - reads on master (and from parent directory) and sends
1279 // dictionary to slaves
1280 // - distribute - reads on potentially a different number of processors
1281 // than it writes to
1282 // - reconstruct - reads parallel, write on master only and to parent
1283 // directory
1284
1285 // Running data-distributed?
1286 // (processors cannot see all other processors' files)
1287 const bool hasDistributedFiles(fileHandler().distributed());
1288
1289
1290 // Set up loose processorsXXX directory matching (for collated) so e.g.
1291 // when checking for -latestTime we don't miss anything. Once we know
1292 // the time, actual number of processors etc we switch back to strict
1293 // matching.
1295
1296 // Need this line since we don't include "setRootCase.H"
1297 #include "foamDlOpenLibs.H"
1298
1299 const bool reconstruct = args.found("reconstruct");
1300 const bool writeCellDist = args.found("cellDist");
1301 const bool dryrun = args.dryRun();
1302 const bool newTimes = args.found("newTimes");
1303 const int optVerbose = args.verbose();
1304
1305 const bool doFiniteArea = !args.found("no-finite-area");
1306 bool decompose = args.found("decompose");
1307 bool overwrite = args.found("overwrite");
1308
1309 // Field restrictions...
1310 const wordRes selectedFields;
1311 const wordRes selectedLagrangianFields;
1312
1313
1314 if (!UPstream::parRun())
1315 {
1317 << ": This utility can only be run parallel"
1318 << exit(FatalError);
1319 }
1320
1321 if (decompose && reconstruct)
1322 {
1324 << "Cannot specify both -decompose and -reconstruct"
1325 << exit(FatalError);
1326 }
1327
1328 if (optVerbose)
1329 {
1330 // Report on output
1333
1334 if (optVerbose > 1)
1335 {
1336 Info<< "Additional debugging enabled" << nl << endl;
1337 ::debug = 1;
1338 }
1339 }
1340
1341 // Disable NaN setting and floating point error trapping. This is to avoid
1342 // any issues inside the field redistribution inside fvMeshDistribute
1343 // which temporarily moves processor faces into existing patches. These
1344 // will now not have correct values. After all bits have been assembled
1345 // the processor fields will be restored but until then there might
1346 // be uninitialised values which might upset some patch field constructors.
1347 // Workaround by disabling floating point error trapping. TBD: have
1348 // actual field redistribution instead of subsetting inside
1349 // fvMeshDistribute.
1350 Foam::sigFpe::unset(true);
1351
1352
1353 // File handlers (read/write)
1354 // ==========================
1355
1356 // Read handler on processors with a volMesh
1357 refPtr<fileOperation> volMeshReadHandler;
1358
1359 // Read handler on processors with an areaMesh
1360 refPtr<fileOperation> areaMeshReadHandler;
1361
1362 // Handler for master-only operation (read/writing from/to undecomposed)
1363 // - only the 'real' master, not io-rank masters
1364 refPtr<fileOperation> masterOnlyHandler;
1365
1367 {
1368 const bool oldParRun = UPstream::parRun(false);
1369
1370 masterOnlyHandler = fileOperation::NewUncollated();
1371
1372 UPstream::parRun(oldParRun);
1373 }
1374
1375 if (decompose)
1376 {
1377 InfoOrPout<< "Decomposing case (like decomposePar)"
1378 << nl << endl;
1379 }
1380 else if (reconstruct)
1381 {
1382 InfoOrPout<< "Reconstructing case (like reconstructParMesh)"
1383 << nl << endl;
1384 }
1385
1386 if (decompose || reconstruct)
1387 {
1388 // The UPstream::nProcs is either the source or destination procs
1390 InfoOrPout<< "Switching to exact matching for "
1392 << " processor directories"
1393 << nl << endl;
1394 }
1395 else
1396 {
1397 // Redistribute mode. Accept any processorsXXX naming since we don't
1398 // know yet what the source/target number of processors is
1400 InfoOrPout<< "Switching to matching any "
1401 << fileOperation::processorsBaseDir << " directory" << nl << endl;
1402 }
1403
1404
1405 if ((decompose || reconstruct) && !overwrite)
1406 {
1407 overwrite = true;
1408
1409 Warning
1410 << nl << " "
1411 << "Added implicit -overwrite for decompose/reconstruct modes"
1412 << nl << endl;
1413 }
1414
1415 if
1416 (
1417 fileHandler().ioRanks().contains(UPstream::myProcNo())
1418 && !Foam::isDir(args.rootPath())
1419 )
1420 {
1421 //FatalErrorInFunction
1423 << ": cannot open root directory " << args.rootPath()
1424 << endl;
1425 //<< exit(FatalError);
1426 }
1427
1428
1429 if (hasDistributedFiles)
1430 {
1431 InfoOrPout<< "Detected multiple roots i.e. non-nfs running"
1432 << nl << endl;
1433
1434 fileHandler().mkDir(args.globalPath());
1435 }
1436
1437
1438 // Check if we have processor directories. Ideally would like to
1439 // use fileHandler().dirPath here but we don't have runTime yet and
1440 // want to delay constructing runTime until we've synced all time
1441 // directories...
1442 const fileName procDir(fileHandler().filePath(args.path()));
1443 if (Foam::isDir(procDir))
1444 {
1445 if (decompose)
1446 {
1447 InfoOrPout<< "Removing existing processor directory:"
1448 << args.relativePath(procDir) << endl;
1449 fileHandler().rmDir(procDir);
1450 }
1451 }
1452 else
1453 {
1454 // Directory does not exist. If this happens on master -> decompose mode
1455 if (UPstream::master() && !reconstruct && !decompose)
1456 {
1457 decompose = true;
1458 InfoOrPout
1459 << "No processor directories; switching on decompose mode"
1460 << nl << endl;
1461 }
1462 }
1463 // If master changed to decompose mode make sure all nodes know about it
1464 Pstream::broadcast(decompose);
1465 if (decompose)
1466 {
1467 // The UPstream::nProcs is either the source or destination procs
1469 InfoOrPout
1470 << "Switching to exact matching for "
1472 << " processor directories"
1473 << nl << endl;
1474 }
1475
1476
1477
1478 // If running distributed we have problem of new processors not finding
1479 // a system/controlDict. However if we switch on the master-only reading
1480 // the problem becomes that the time directories are differing sizes and
1481 // e.g. latestTime will pick up a different time (which causes createTime.H
1482 // to abort). So for now make sure to have master times on all
1483 // processors
1484 if (!reconstruct)
1485 {
1486 InfoOrPout<< "Creating time directories on all processors"
1487 << nl << endl;
1488 createTimeDirs(args.path());
1489 }
1490
1491 // Construct time
1492 // ~~~~~~~~~~~~~~
1493
1494 // Replace #include "createTime.H" with our own version
1495 // that has MUST_READ instead of READ_MODIFIED
1496
1497 Info<< "Create time\n" << endl;
1499 (
1501 args,
1502 false, // no enableFunctionObjects
1503 true, // permit enableLibs
1504 IOobjectOption::MUST_READ // Without watching
1505 );
1506
1507
1508 refPtr<fileOperation> writeHandler;
1509
1510
1511 runTime.functionObjects().off(); // Extra safety?
1512 // Make sure that no runTime checking is done since fileOperation::addWatch
1513 // etc. does broadcast over world, even if constructed only on a subset
1514 // of procs
1515 runTime.runTimeModifiable(false);
1516
1517 // Save local processor0 casename
1518 const fileName proc0CaseName = runTime.caseName();
1519 const bool oldProcCase = runTime.processorCase();
1520
1521
1522 // Construct undecomposed Time
1523 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~
1524 // This will read the same controlDict but might have a different
1525 // set of times so enforce same times
1526
1527 if (hasDistributedFiles)
1528 {
1529 InfoOrPout<< "Creating time directories for undecomposed Time"
1530 << " on all processors" << nl << endl;
1531 createTimeDirs(args.globalPath());
1532 }
1533
1534
1535 InfoOrPout<< "Create undecomposed database" << nl << endl;
1536 Time baseRunTime
1537 (
1538 runTime.controlDict(),
1539 runTime.rootPath(),
1540 runTime.globalCaseName(),
1541 runTime.system(),
1542 runTime.constant(),
1543 false, // No function objects
1544 false, // No extra controlDict libs
1545 IOobjectOption::MUST_READ // Without watching
1546 );
1547 // Make sure that no runTime checking is done since fileOperation::addWatch
1548 // etc. does broadcast over world, even if constructed only on a subset
1549 // of procs
1550 baseRunTime.runTimeModifiable(false);
1551
1552
1553 wordHashSet masterTimeDirSet;
1554 if (newTimes)
1555 {
1556 instantList baseTimeDirs(baseRunTime.times());
1557 for (const instant& t : baseTimeDirs)
1558 {
1559 masterTimeDirSet.insert(t.name());
1560 }
1561 }
1562
1563
1564 // Allow override of decomposeParDict location
1565 const fileName decompDictFile =
1566 args.getOrDefault<fileName>("decomposeParDict", "");
1567
1568 // Get region names
1569 #include "getAllRegionOptions.H"
1570
1571 if (regionNames.size() == 1 && regionNames[0] != polyMesh::defaultRegion)
1572 {
1573 InfoOrPout<< "Using region: " << regionNames[0] << nl << endl;
1574 }
1575
1576
1577 // Demand driven lagrangian mapper
1578 autoPtr<parLagrangianDistributor> lagrangianDistributorPtr;
1579
1580 if (reconstruct)
1581 {
1582 // use the times list from the master processor
1583 // and select a subset based on the command-line options
1584 instantList timeDirs = timeSelector::select(runTime.times(), args);
1585 Pstream::broadcast(timeDirs);
1586
1587 if (timeDirs.empty())
1588 {
1590 << "No times selected"
1591 << exit(FatalError);
1592 }
1593
1594
1595 // Pass1 : reconstruct mesh and addressing
1596 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1597
1598
1599 InfoOrPout
1600 << nl
1601 << "Reconstructing mesh and addressing" << nl << endl;
1602
1603 forAll(regionNames, regioni)
1604 {
1605 const word& regionName = regionNames[regioni];
1606
1607 const fileName volMeshSubDir
1608 (
1610 );
1611 const fileName areaMeshSubDir
1612 (
1613 // Assume single-region area mesh
1615 );
1616
1617 InfoOrPout
1618 << nl
1619 << "Reconstructing mesh:"
1621
1622 bool areaMeshDetected = false;
1623
1624 // Loop over all times
1625 forAll(timeDirs, timeI)
1626 {
1627 // Set time for global database
1628 runTime.setTime(timeDirs[timeI], timeI);
1629 baseRunTime.setTime(timeDirs[timeI], timeI);
1630
1631 InfoOrPout<< "Time = " << runTime.timeName() << endl << endl;
1632
1633 // Where meshes are
1634 fileName volMeshInstance;
1635 fileName areaMeshInstance;
1636
1637 volMeshInstance = runTime.findInstance
1638 (
1639 volMeshSubDir,
1640 "faces",
1642 );
1643
1644 if (doFiniteArea)
1645 {
1646 areaMeshInstance = runTime.findInstance
1647 (
1648 areaMeshSubDir,
1649 "faceLabels",
1651 );
1652 }
1653
1655 (
1657 volMeshInstance,
1658 areaMeshInstance
1659 );
1660
1661
1662 // Check processors have meshes
1663 // - check for 'faces' file (polyMesh)
1664 // - check for 'faceLabels' file (faMesh)
1665 boolList volMeshOnProc;
1666 boolList areaMeshOnProc(UPstream::nProcs(), false);
1667
1668 volMeshOnProc = haveMeshFile
1669 (
1670 runTime,
1671 volMeshInstance/volMeshSubDir,
1672 "faces"
1673 );
1674
1675 // Create handler for reading
1676 newHandler(volMeshOnProc, volMeshReadHandler);
1677
1678 if (doFiniteArea)
1679 {
1680 areaMeshOnProc = haveMeshFile
1681 (
1682 runTime,
1683 areaMeshInstance/areaMeshSubDir,
1684 "faceLabels"
1685 );
1686 areaMeshDetected = areaMeshOnProc.found(true);
1687
1688 if (areaMeshOnProc == volMeshOnProc)
1689 {
1690 if (volMeshReadHandler)
1691 {
1692 // Use same reader for faMesh as for fvMesh
1693 areaMeshReadHandler.ref(volMeshReadHandler.ref());
1694 }
1695 }
1696 else
1697 {
1698 newHandler(areaMeshOnProc, areaMeshReadHandler);
1699 }
1700 }
1701
1702
1703 // Addressing back to reconstructed mesh as xxxProcAddressing.
1704 // - all processors have consistent faceProcAddressing
1705 // - processors without a mesh don't need faceProcAddressing
1706
1707
1708 // Note: filePath searches up on processors that don't have
1709 // processor if instance = constant so explicitly check
1710 // found filename.
1711 bool haveVolAddressing = false;
1712 if (volMeshOnProc[Pstream::myProcNo()])
1713 {
1714 // Read faces (just to know their size)
1715 faceCompactIOList faces
1716 (
1717 IOobject
1718 (
1719 "faces",
1720 volMeshInstance,
1721 volMeshSubDir,
1722 runTime,
1724 )
1725 );
1726
1727 // Check faceProcAddressing
1728 labelIOList faceProcAddressing
1729 (
1730 IOobject
1731 (
1732 "faceProcAddressing",
1733 volMeshInstance,
1734 volMeshSubDir,
1735 runTime,
1737 )
1738 );
1739
1740 haveVolAddressing =
1741 (
1742 faceProcAddressing.headerOk()
1743 && faceProcAddressing.size() == faces.size()
1744 );
1745 }
1746 else
1747 {
1748 // Have no mesh. Don't need addressing
1749 haveVolAddressing = true;
1750 }
1751
1752 bool haveAreaAddressing = false;
1753 if (areaMeshOnProc[Pstream::myProcNo()])
1754 {
1755 // Read faces (just to know their size)
1757 (
1758 IOobject
1759 (
1760 "faceLabels",
1761 areaMeshInstance,
1762 areaMeshSubDir,
1763 runTime,
1765 )
1766 );
1767
1768 // Check faceProcAddressing
1769 labelIOList faceProcAddressing
1770 (
1771 IOobject
1772 (
1773 "faceProcAddressing",
1774 areaMeshInstance,
1775 areaMeshSubDir,
1776 runTime,
1778 )
1779 );
1780
1781 haveAreaAddressing =
1782 (
1783 faceProcAddressing.headerOk()
1784 && faceProcAddressing.size() == faceLabels.size()
1785 );
1786 }
1787 else if (areaMeshDetected)
1788 {
1789 // Have no mesh. Don't need addressing
1790 haveAreaAddressing = true;
1791 }
1792
1793
1794 // Additionally check for master faces being readable. Could
1795 // do even more checks, e.g. global number of cells same
1796 // as cellProcAddressing
1797
1798 bool volMeshHaveUndecomposed = false;
1799 bool areaMeshHaveUndecomposed = false;
1800
1801 if (Pstream::master())
1802 {
1803 InfoOrPout
1804 << "Checking " << baseRunTime.caseName()
1805 << " for undecomposed volume and area meshes..."
1806 << endl;
1807
1808 const bool oldParRun = Pstream::parRun(false);
1809 const label oldNumProcs(fileHandler().nProcs());
1810
1811 // Volume
1812 {
1813 faceCompactIOList facesIO
1814 (
1815 IOobject
1816 (
1817 "faces",
1818 volMeshInstance,
1819 volMeshSubDir,
1820 baseRunTime,
1822 ),
1823 label(0)
1824 );
1825 volMeshHaveUndecomposed = facesIO.headerOk();
1826 }
1827
1828 // Area
1829 if (doFiniteArea)
1830 {
1831 labelIOList labelsIO
1832 (
1833 IOobject
1834 (
1835 "faceLabels",
1836 areaMeshInstance,
1837 areaMeshSubDir,
1838 baseRunTime,
1840 )
1841 );
1842 areaMeshHaveUndecomposed = labelsIO.headerOk();
1843 }
1844
1845 const_cast<fileOperation&>
1846 (
1847 fileHandler()
1848 ).nProcs(oldNumProcs);
1849 Pstream::parRun(oldParRun); // Restore parallel state
1850 }
1851
1853 (
1855 volMeshHaveUndecomposed,
1856 areaMeshHaveUndecomposed
1857 );
1858
1859 // Report
1860 {
1861 InfoOrPout
1862 << " volume mesh ["
1863 << volMeshHaveUndecomposed << "] : "
1864 << volMeshInstance << nl
1865 << " area mesh ["
1866 << areaMeshHaveUndecomposed << "] : "
1867 << areaMeshInstance << nl
1868 << endl;
1869 }
1870
1871
1872 if
1873 (
1874 !volMeshHaveUndecomposed
1875 || !returnReduceAnd(haveVolAddressing)
1876 )
1877 {
1878 InfoOrPout
1879 << "No undecomposed mesh. Creating from: "
1880 << volMeshInstance << endl;
1881
1882 if (areaMeshHaveUndecomposed)
1883 {
1884 areaMeshHaveUndecomposed = false;
1885 InfoOrPout
1886 << "Also ignore any undecomposed area mesh"
1887 << endl;
1888 }
1889
1891 (
1892 IOobject
1893 (
1894 regionName,
1895 volMeshInstance,
1896 runTime,
1898 ),
1899 volMeshReadHandler
1900 );
1901 fvMeshTools::setBasicGeometry(volMeshPtr());
1902 fvMesh& mesh = volMeshPtr();
1903
1904
1905 InfoOrPout<< nl << "Reconstructing mesh" << nl << endl;
1906
1907 // Reconstruct (1 processor)
1908 const label nDestProcs(1);
1909 const labelList finalDecomp(mesh.nCells(), Zero);
1910
1911 redistributeAndWrite
1912 (
1913 volMeshReadHandler,
1914 masterOnlyHandler, //writeHandler,
1915 baseRunTime,
1916 proc0CaseName,
1917
1918 // Controls
1919 false, // do not read fields
1920 false, // do not read undecomposed case on proc0
1921 true, // write redistributed files to proc0
1922 overwrite,
1923
1924 // Decomposition information
1925 nDestProcs,
1926 finalDecomp,
1927
1928 // For finite-volume
1929 volMeshOnProc,
1930 volMeshInstance,
1931 mesh
1932 );
1933 }
1934
1935
1936 // Similarly for finiteArea
1937 // - may or may not have undecomposed mesh
1938 // - may or may not have decomposed meshes
1939
1940 if
1941 (
1942 areaMeshOnProc.found(true) // ie, areaMeshDetected
1943 &&
1944 (
1945 !areaMeshHaveUndecomposed
1946 || !returnReduceAnd(haveAreaAddressing)
1947 )
1948 )
1949 {
1950 InfoOrPout
1951 << "Loading area mesh from "
1952 << areaMeshInstance << endl;
1953
1954 InfoOrPout<< " getting volume mesh support" << endl;
1955
1957 (
1958 IOobject
1959 (
1960 regionName,
1961 baseRunTime.timeName(),
1962 baseRunTime,
1964 ),
1965 true // read on master only
1966 );
1967 fvMeshTools::setBasicGeometry(baseMeshPtr());
1968
1970 (
1971 IOobject
1972 (
1973 regionName,
1974 baseMeshPtr().facesInstance(),
1975 runTime,
1977 ),
1978 volMeshReadHandler
1979 );
1980 fvMesh& mesh = volMeshPtr();
1981
1982 // Read volume proc addressing back to base mesh
1984 (
1986 );
1987
1988
1989 //autoPtr<faMesh> areaMeshPtr =
1990 // faMeshTools::loadOrCreateMesh
1991 //(
1992 // IOobject
1993 // (
1994 // regionName,
1995 // areaMeshInstance,
1996 // runTime,
1997 // IOobjectOption::MUST_READ
1998 // ),
1999 // mesh, // <- The referenced polyMesh (from above)
2000 // decompose
2001 //);
2003 (
2004 IOobject
2005 (
2006 regionName,
2007 areaMeshInstance,
2008 runTime,
2010 ),
2011 mesh, // <- The referenced polyMesh (from above)
2012 areaMeshReadHandler
2013 );
2014 faMesh& areaMesh = areaMeshPtr();
2015
2018
2019 autoPtr<faMesh> areaBaseMeshPtr;
2020
2021 // Reconstruct using polyMesh distribute map
2022 mapDistributePolyMesh faDistMap
2023 (
2025 (
2026 areaMesh,
2027 distMap(), // The polyMesh distMap
2028 baseMeshPtr(), // Target polyMesh
2029 areaBaseMeshPtr
2030 )
2031 );
2032
2033 faMeshTools::forceDemandDriven(areaBaseMeshPtr());
2034 faMeshTools::unregisterMesh(areaBaseMeshPtr());
2035
2036
2037 if (Pstream::master())
2038 {
2039 InfoOrPout
2040 << "Setting caseName to " << baseRunTime.caseName()
2041 << " to write reconstructed area mesh." << endl;
2042 runTime.caseName() = baseRunTime.caseName();
2043 const bool oldProcCase(runTime.processorCase(false));
2044 const bool oldParRun(Pstream::parRun(false));
2045 const label oldNumProcs(fileHandler().nProcs());
2046
2047 areaBaseMeshPtr().write();
2048
2049 // Now we've written all. Reset caseName on master
2050 InfoOrPout<< "Restoring caseName" << endl;
2051 const_cast<fileOperation&>
2052 (
2053 fileHandler()
2054 ).nProcs(oldNumProcs);
2055 Pstream::parRun(oldParRun);
2056 runTime.caseName() = proc0CaseName;
2057 runTime.processorCase(oldProcCase);
2058 }
2059
2060 // Update for the reconstructed procAddressing
2062 (
2063 areaBaseMeshPtr(), // Reconstruct location
2064 faDistMap,
2065 false, // decompose=false
2066 writeHandler, //writeHandler,
2067 areaMeshPtr.get() // procMesh
2068 );
2069 }
2070 }
2071
2072 // Make sure all is finished writing until re-reading in pass2
2073 // below
2074 fileHandler().flush();
2075
2076
2077 // Pass2 : read mesh and addressing and reconstruct fields
2078 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
2079
2080 InfoOrPout
2081 << nl
2082 << "Reconstructing fields" << nl << endl;
2083
2084 runTime.setTime(timeDirs[0], 0);
2085 baseRunTime.setTime(timeDirs[0], 0);
2086 InfoOrPout<< "Time = " << runTime.timeName() << endl << endl;
2087
2088
2089 // Read undecomposed mesh on master and 'empty' mesh
2090 // (zero faces, point, cells but valid patches and zones) on slaves.
2091 // This is a bit of tricky code and hidden inside fvMeshTools for
2092 // now.
2093 InfoOrPout<< "Reading undecomposed mesh (on master)" << endl;
2094 //autoPtr<fvMesh> baseMeshPtr = fvMeshTools::newMesh
2095 //(
2096 // IOobject
2097 // (
2098 // regionName,
2099 // baseRunTime.timeName(),
2100 // baseRunTime,
2101 // IOobjectOption::MUST_READ
2102 // ),
2103 // true // read on master only
2104 //);
2105 fileName facesInstance;
2106 fileName pointsInstance;
2108 (
2109 IOobject
2110 (
2111 regionName,
2112 baseRunTime.timeName(),
2113 baseRunTime,
2115 ),
2116 facesInstance,
2117 pointsInstance
2118 );
2119
2121 (
2122 IOobject
2123 (
2124 regionName,
2125 facesInstance, //baseRunTime.timeName(),
2126 baseRunTime,
2128 ),
2129 masterOnlyHandler // read on master only
2130 );
2131
2132 if (::debug)
2133 {
2134 Pout<< "Undecomposed mesh :"
2135 << " instance:" << baseMeshPtr().facesInstance()
2136 << " nCells:" << baseMeshPtr().nCells()
2137 << " nFaces:" << baseMeshPtr().nFaces()
2138 << " nPoints:" << baseMeshPtr().nPoints()
2139 << endl;
2140 }
2141
2142 InfoOrPout<< "Reading local, decomposed mesh" << endl;
2144 (
2145 IOobject
2146 (
2147 regionName,
2148 baseMeshPtr().facesInstance(),
2149 runTime,
2151 ),
2152 volMeshReadHandler // read on fvMesh processors
2153 );
2154 fvMesh& mesh = volMeshPtr();
2155
2156
2157 // Similarly for finiteArea
2158 autoPtr<faMesh> areaBaseMeshPtr;
2159 autoPtr<faMesh> areaMeshPtr;
2160 autoPtr<faMeshDistributor> faDistributor;
2161 mapDistributePolyMesh areaDistMap;
2162
2163 if (areaMeshDetected)
2164 {
2165 //areaBaseMeshPtr = faMeshTools::newMesh
2166 //(
2167 // IOobject
2168 // (
2169 // regionName,
2170 // baseRunTime.timeName(),
2171 // baseRunTime,
2172 // IOobjectOption::MUST_READ
2173 // ),
2174 // baseMeshPtr(),
2175 // true // read on master only
2176 //);
2177 areaBaseMeshPtr = faMeshTools::loadOrCreateMesh
2178 (
2179 IOobject
2180 (
2181 regionName,
2182 baseMeshPtr().facesInstance(),
2183 baseRunTime,
2185 ),
2186 baseMeshPtr(),
2187 masterOnlyHandler
2188 );
2189
2190 //areaMeshPtr = faMeshTools::loadOrCreateMesh
2191 //(
2192 // IOobject
2193 // (
2194 // regionName,
2195 // areaBaseMeshPtr().facesInstance(),
2196 // runTime,
2197 // IOobjectOption::MUST_READ
2198 // ),
2199 // mesh,
2200 // decompose
2201 //);
2202 areaMeshPtr = faMeshTools::loadOrCreateMesh
2203 (
2204 IOobject
2205 (
2206 regionName,
2207 areaBaseMeshPtr().facesInstance(),
2208 runTime,
2210 ),
2211 mesh,
2212 areaMeshReadHandler
2213 );
2214
2215 areaDistMap =
2217 (
2218 areaMeshPtr(),
2219 areaBaseMeshPtr
2220 );
2221
2222 faMeshTools::forceDemandDriven(areaMeshPtr());
2223
2224 // Create an appropriate field distributor
2225 faDistributor.reset
2226 (
2228 (
2229 areaMeshPtr(), // source
2230 areaBaseMeshPtr(), // target
2231 areaDistMap,
2232 masterOnlyHandler // only write on master
2233 )
2234 );
2235 // Report some messages. Tbd.
2237 }
2238
2239
2240 // Read addressing back to base mesh
2242 distMap = fvMeshTools::readProcAddressing(mesh, baseMeshPtr);
2243
2244 // Construct field mapper
2245 auto fvDistributorPtr =
2247 (
2248 mesh, // source
2249 baseMeshPtr(), // target
2250 distMap(),
2251 masterOnlyHandler // Write on master only
2252 );
2253
2254 // Construct point field mapper
2255 const auto& basePointMesh = pointMesh::New(baseMeshPtr());
2256 const auto& procPointMesh = pointMesh::New(mesh);
2257
2258 auto pointFieldDistributorPtr =
2260 (
2261 procPointMesh, // source
2262 basePointMesh, // target
2263 distMap(),
2264 false, // delay
2265 //UPstream::master() // Write reconstructed on master
2266 masterOnlyHandler // Write on master only
2267 );
2268
2269
2270 // Since we start from Times[0] and not runTime.timeName() we
2271 // might overlook point motion in the first timestep
2272 // (since mesh.readUpdate() below will not be triggered). Instead
2273 // detect points by hand
2274 if (mesh.pointsInstance() != mesh.facesInstance())
2275 {
2276 InfoOrPout
2277 << " Detected initial mesh motion;"
2278 << " reconstructing points" << nl
2279 << endl;
2280 fvDistributorPtr().reconstructPoints();
2281 }
2282
2283
2284 // Loop over all times
2285 forAll(timeDirs, timeI)
2286 {
2287 if (newTimes && masterTimeDirSet.found(timeDirs[timeI].name()))
2288 {
2289 InfoOrPout
2290 << "Skipping time " << timeDirs[timeI].name()
2291 << nl << endl;
2292 continue;
2293 }
2294
2295 // Set time for global database
2296 runTime.setTime(timeDirs[timeI], timeI);
2297 baseRunTime.setTime(timeDirs[timeI], timeI);
2298
2299 InfoOrPout<< "Time = " << runTime.timeName() << endl << endl;
2300
2301
2302 // Check if any new meshes need to be read.
2303 polyMesh::readUpdateState procStat = mesh.readUpdate();
2304
2305 if (procStat == polyMesh::POINTS_MOVED)
2306 {
2307 InfoOrPout
2308 << " Detected mesh motion; reconstructing points"
2309 << nl << endl;
2310 fvDistributorPtr().reconstructPoints();
2311 }
2312 else if
2313 (
2314 procStat == polyMesh::TOPO_CHANGE
2315 || procStat == polyMesh::TOPO_PATCH_CHANGE
2316 )
2317 {
2318 InfoOrPout
2319 << " Detected topology change;"
2320 << " reconstructing addressing" << nl << endl;
2321
2322 if (baseMeshPtr)
2323 {
2324 // Cannot do a baseMesh::readUpdate() since not all
2325 // processors will have mesh files. So instead just
2326 // recreate baseMesh
2327 baseMeshPtr.clear();
2328 //baseMeshPtr = fvMeshTools::newMesh
2329 //(
2330 // IOobject
2331 // (
2332 // regionName,
2333 // baseRunTime.timeName(),
2334 // baseRunTime,
2335 // IOobjectOption::MUST_READ
2336 // ),
2337 // true // read on master only
2338 //);
2339 baseMeshPtr = fvMeshTools::loadOrCreateMesh
2340 (
2341 IOobject
2342 (
2343 regionName,
2344 baseRunTime.timeName(),
2345 baseRunTime,
2347 ),
2348 masterOnlyHandler // read on master only
2349 );
2350 if (::debug)
2351 {
2352 Pout<< "Undecomposed mesh :"
2353 << " nCells:" << baseMeshPtr().nCells()
2354 << " nFaces:" << baseMeshPtr().nFaces()
2355 << " nPoints:" << baseMeshPtr().nPoints()
2356 << endl;
2357 }
2358 }
2359
2360 // Re-read procAddressing
2361 distMap =
2363
2364 // Reset field mappers
2365
2366 fvDistributorPtr.reset
2367 (
2369 (
2370 mesh, // source
2371 baseMeshPtr(), // target
2372 distMap(),
2373 masterOnlyHandler // Write on master only
2374 )
2375 );
2376
2377 // Construct point field mapper
2378 const auto& basePointMesh = pointMesh::New(baseMeshPtr());
2379 const auto& procPointMesh = pointMesh::New(mesh);
2380
2381 pointFieldDistributorPtr.reset
2382 (
2384 (
2385 procPointMesh, // source
2386 basePointMesh, // target
2387 distMap(),
2388 false, // delay until later
2389 masterOnlyHandler // Write on master only
2390 )
2391 );
2392
2393 lagrangianDistributorPtr.reset();
2394
2395 if (areaMeshPtr)
2396 {
2397 InfoOrPout
2398 << " Discarding finite-area addressing"
2399 << " (TODO)" << nl << endl;
2400
2401 areaBaseMeshPtr.reset();
2402 areaMeshPtr.reset();
2403 faDistributor.reset();
2404 areaDistMap.clear();
2405 }
2406 }
2407
2408
2409 // Get list of objects
2410 IOobjectList objects(mesh, runTime.timeName());
2411
2412 // Mesh fields (vol, surface, volInternal)
2413 fvDistributorPtr()
2414 .distributeAllFields(objects, selectedFields);
2415
2416 // pointfields
2417 // - distribute and write (verbose)
2418 pointFieldDistributorPtr()
2419 .distributeAllFields(objects, selectedFields);
2420
2421
2422 // Clouds (note: might not be present on all processors)
2424 (
2425 lagrangianDistributorPtr,
2426 baseMeshPtr(),
2427 mesh,
2428 distMap(),
2429 selectedLagrangianFields
2430 //masterOnlyHandler
2431 );
2432
2433 if (faDistributor)
2434 {
2435 faDistributor()
2436 .distributeAllFields(objects, selectedFields);
2437 }
2438
2439
2440 // If there are any "uniform" directories copy them from
2441 // the master processor
2442
2443 copyUniform
2444 (
2445 volMeshReadHandler, //masterOnlyHandler, // readHandler
2446 masterOnlyHandler, // writeHandler
2447
2448 true, // reconstruct
2449 false, // decompose
2450
2451 mesh.time().timeName(),
2452 word::null, // optional caseName for reading
2453 mesh,
2454 baseMeshPtr()
2455 );
2456 // Non-region specific. Note: should do outside region loop
2457 // but would then have to replicate the whole time loop ...
2458 copyUniform
2459 (
2460 volMeshReadHandler, //masterOnlyHandler, // readHandler,
2461 masterOnlyHandler, // writeHandler
2462
2463 true, // reconstruct
2464 false, // decompose
2465
2466 mesh.time().timeName(),
2467 word::null, // optional caseName for reading
2468 mesh.time(), // runTime
2469 baseMeshPtr().time() // baseRunTime
2470 );
2471 }
2472 }
2473 }
2474 else
2475 {
2476 // decompose or redistribution mode.
2477 // decompose : master : read from parent dir
2478 // slave : dummy mesh
2479 // redistribute : all read mesh or dummy mesh
2480
2481 Time& readRunTime =
2482 (
2483 (decompose)
2484 ? baseRunTime
2485 : runTime
2486 );
2487
2488 // Time coming from processor0 (or undecomposed if no processor0)
2489 scalar masterTime = timeSelector::selectIfPresent
2490 (
2491 readRunTime,
2492 args
2493 )[0].value();
2494 Pstream::broadcast(masterTime);
2495 InfoOrPout
2496 << "Setting time to that of master or undecomposed case : "
2497 << masterTime << endl;
2498 runTime.setTime(masterTime, 0);
2499 baseRunTime.setTime(masterTime, 0);
2500
2501
2502
2503
2504 // Save old time name (since might be incremented)
2505 const word oldTimeName(runTime.timeName());
2506
2507 forAll(regionNames, regioni)
2508 {
2509 const word& regionName = regionNames[regioni];
2510
2511 const fileName volMeshSubDir
2512 (
2514 );
2515 const fileName areaMeshSubDir
2516 (
2517 // Assume single-region area mesh
2519 );
2520
2521 InfoOrPout
2522 << nl << nl
2523 << (decompose ? "Decomposing" : "Redistributing")
2524 << " mesh:" << polyMesh::regionName(regionName) << nl << endl;
2525
2526
2527 // Get time instance directory
2528 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~
2529 // At this point we should be able to read at least a mesh on
2530 // processor0. Note the changing of the processor0 casename to
2531 // enforce it to read/write from the undecomposed case
2532
2533 fileName volMeshMasterInstance;
2534 fileName areaMeshMasterInstance;
2535
2536 // Assume to be true
2537 bool volMeshHaveUndecomposed = true;
2538 bool areaMeshHaveUndecomposed = doFiniteArea;
2539
2540 if (Pstream::master())
2541 {
2542 if (decompose)
2543 {
2544 InfoOrPout
2545 << "Checking undecomposed mesh in case: "
2546 << baseRunTime.caseName() << endl;
2547 runTime.caseName() = baseRunTime.caseName();
2548 runTime.processorCase(false);
2549 }
2550
2551 const bool oldParRun = Pstream::parRun(false);
2552 const label oldNumProcs(fileHandler().nProcs());
2553 volMeshMasterInstance = readRunTime.findInstance
2554 (
2555 volMeshSubDir,
2556 "faces",
2558 );
2559
2560 if (doFiniteArea)
2561 {
2562 areaMeshMasterInstance = readRunTime.findInstance
2563 (
2564 areaMeshSubDir,
2565 "faceLabels",
2567 );
2568
2569 // Note: findInstance returns "constant" even if not found,
2570 // so recheck now for a false positive.
2571
2572 if ("constant" == areaMeshMasterInstance)
2573 {
2574 const boolList areaMeshOnProc
2575 (
2577 (
2578 readRunTime,
2579 areaMeshMasterInstance/areaMeshSubDir,
2580 "faceLabels",
2581 false // verbose=false
2582 )
2583 );
2584
2585 if (areaMeshOnProc.empty() || !areaMeshOnProc[0])
2586 {
2587 areaMeshHaveUndecomposed = false;
2588 }
2589 }
2590 }
2591
2592 const_cast<fileOperation&>(fileHandler()).nProcs(oldNumProcs);
2593 Pstream::parRun(oldParRun); // Restore parallel state
2594
2595 if (decompose)
2596 {
2597 InfoOrPout
2598 << " volume mesh ["
2599 << volMeshHaveUndecomposed << "] : "
2600 << volMeshMasterInstance << nl
2601 << " area mesh ["
2602 << areaMeshHaveUndecomposed << "] : "
2603 << areaMeshMasterInstance << nl
2604 << nl << nl;
2605
2606 // Restoring caseName
2607 InfoOrPout<< "Restoring caseName" << endl;
2608 runTime.caseName() = proc0CaseName;
2609 runTime.processorCase(oldProcCase);
2610 }
2611 }
2612
2614 (
2616 volMeshHaveUndecomposed,
2617 areaMeshHaveUndecomposed,
2618 volMeshMasterInstance,
2619 areaMeshMasterInstance
2620 );
2621
2622 // Check processors have meshes
2623 // - check for 'faces' file (polyMesh)
2624 // - check for 'faceLabels' file (faMesh)
2625 boolList volMeshOnProc;
2626 boolList areaMeshOnProc;
2627
2628 if (decompose)
2629 {
2630 // Already determined above that master can read 'faces' file.
2631 // This avoids doing all the casename setting/restoring again.
2632 volMeshOnProc.setSize(UPstream::nProcs(), false);
2633 volMeshOnProc[UPstream::masterNo()] = volMeshHaveUndecomposed;
2634 }
2635 else
2636 {
2637 // All check if can read 'faces' file
2638 volMeshOnProc = haveMeshFile
2639 (
2640 runTime,
2641 volMeshMasterInstance/volMeshSubDir,
2642 "faces"
2643 );
2644 }
2645
2646 // Create handler for reading
2647 newHandler(volMeshOnProc, volMeshReadHandler);
2648
2649
2650 // Now we've determined which processors are reading switch back
2651 // to exact matching of 'processorsXXX' directory names.
2652 // - this determines the processorsXXX fileNames
2653 // - the XXX comes from the number of read processors
2654 // - also adapt the masterOnlyReader (used in copyUniform)
2655
2657 (
2658 findIndices(volMeshOnProc, true).size()
2659 );
2660
2661
2662 if (doFiniteArea)
2663 {
2664 if (decompose)
2665 {
2666 // Already determined above that master can read
2667 // 'faceLabels' file.
2668 areaMeshOnProc.setSize(UPstream::nProcs(), false);
2669 areaMeshOnProc[UPstream::masterNo()] =
2670 areaMeshHaveUndecomposed;
2671 }
2672 else
2673 {
2674 areaMeshOnProc = haveMeshFile
2675 (
2676 runTime,
2677 areaMeshMasterInstance/areaMeshSubDir,
2678 "faceLabels"
2679 );
2680 }
2681
2682 // Create handler for reading
2683 if (areaMeshOnProc == volMeshOnProc)
2684 {
2685 if (volMeshReadHandler)
2686 {
2687 // Use same reader for faMesh as for fvMesh
2688 areaMeshReadHandler.ref(volMeshReadHandler.ref());
2689 }
2690 }
2691 else
2692 {
2693 newHandler(areaMeshOnProc, areaMeshReadHandler);
2694 }
2695 }
2696
2697
2698 // Prior to loadOrCreateMesh, note which meshes already exist
2699 // for the current file handler.
2700 // - where mesh would be written if it didn't exist already.
2701 fileNameList volMeshDir(Pstream::nProcs());
2702 {
2703 volMeshDir[Pstream::myProcNo()] =
2704 (
2705 fileHandler().objectPath
2706 (
2707 IOobject
2708 (
2709 "faces",
2710 volMeshMasterInstance/volMeshSubDir,
2711 runTime
2712 ),
2714 ).path()
2715 );
2716
2717 Pstream::allGatherList(volMeshDir);
2718
2719 if (optVerbose && Pstream::master())
2720 {
2721 Info<< "Per processor faces dirs:" << nl
2722 << '(' << nl;
2723
2724 for (const int proci : Pstream::allProcs())
2725 {
2726 Info<< " "
2727 << runTime.relativePath(volMeshDir[proci]);
2728
2729 if (!volMeshOnProc[proci])
2730 {
2731 Info<< " [missing]";
2732 }
2733 Info<< nl;
2734 }
2735 Info<< ')' << nl << endl;
2736 }
2737 }
2738
2739 fileNameList areaMeshDir(Pstream::nProcs());
2740 if (doFiniteArea)
2741 {
2742 areaMeshDir[Pstream::myProcNo()] =
2743 (
2744 fileHandler().objectPath
2745 (
2746 IOobject
2747 (
2748 "faceLabels",
2749 areaMeshMasterInstance/areaMeshSubDir,
2750 runTime
2751 ),
2753 ).path()
2754 );
2755
2756 Pstream::allGatherList(areaMeshDir);
2757
2758 if (optVerbose && Pstream::master())
2759 {
2760 Info<< "Per processor faceLabels dirs:" << nl
2761 << '(' << nl;
2762
2763 for (const int proci : Pstream::allProcs())
2764 {
2765 Info<< " "
2766 << runTime.relativePath(areaMeshDir[proci]);
2767
2768 if (!areaMeshOnProc[proci])
2769 {
2770 Info<< " [missing]";
2771 }
2772 Info<< nl;
2773 }
2774 Info<< ')' << nl << endl;
2775 }
2776 }
2777
2778
2779 // Load mesh (or create dummy one)
2780 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
2781
2782 if (decompose)
2783 {
2784 InfoOrPout
2785 << "Setting caseName to " << baseRunTime.caseName()
2786 << " to read undecomposed mesh" << endl;
2787 runTime.caseName() = baseRunTime.caseName();
2788 runTime.processorCase(false);
2789 }
2790
2792 (
2793 IOobject
2794 (
2795 regionName,
2796 volMeshMasterInstance,
2797 runTime,
2799 ),
2800 volMeshReadHandler
2801 );
2802 fvMesh& mesh = volMeshPtr();
2803
2804
2805 // Area mesh
2806
2807 autoPtr<faMesh> areaMeshPtr;
2808
2809 // Decomposing: must have an undecomposed mesh
2810 // Redistributing: have any proc mesh
2811 if
2812 (
2813 doFiniteArea
2814 &&
2815 (
2816 decompose
2817 ? areaMeshHaveUndecomposed
2818 : areaMeshOnProc.found(true)
2819 )
2820 )
2821 {
2822 areaMeshPtr = faMeshTools::loadOrCreateMesh
2823 (
2824 IOobject
2825 (
2826 regionName,
2827 areaMeshMasterInstance,
2828 runTime,
2830 ),
2831 mesh, // <- The referenced polyMesh (from above)
2832 areaMeshReadHandler
2833 );
2834
2835 faMeshTools::forceDemandDriven(*areaMeshPtr);
2836 faMeshTools::unregisterMesh(*areaMeshPtr);
2837 }
2838
2839
2840 if (decompose)
2841 {
2842 InfoOrPout<< "Restoring caseName" << endl;
2843 runTime.caseName() = proc0CaseName;
2844 runTime.processorCase(oldProcCase);
2845 }
2846
2847 const label nOldCells = mesh.nCells();
2848
2849
2850 // Determine decomposition
2851 // ~~~~~~~~~~~~~~~~~~~~~~~
2852
2853 label nDestProcs;
2854 labelList finalDecomp;
2855 determineDecomposition
2856 (
2857 volMeshReadHandler, // how to read decomposeParDict
2858 baseRunTime,
2859 decompDictFile,
2860 decompose,
2861 proc0CaseName,
2862 mesh,
2863 writeCellDist,
2864
2865 nDestProcs,
2866 finalDecomp
2867 );
2868
2869 if (dryrun)
2870 {
2871 continue;
2872 }
2873
2874 if (!writeHandler && nDestProcs < fileHandler().nProcs())
2875 {
2876 boolList isWriteProc(UPstream::nProcs(), false);
2877 isWriteProc.slice(0, nDestProcs) = true;
2878 InfoOrPout
2879 << " dest procs ["
2880 << isWriteProc << "]" << nl
2881 << endl;
2882 newHandler(isWriteProc, writeHandler);
2883 }
2884
2885 // Area fields first. Read and deregister
2887 if (areaMeshPtr)
2888 {
2889 areaFields.read
2890 (
2891 baseRunTime,
2892 proc0CaseName,
2893 decompose,
2894
2895 areaMeshOnProc,
2896 areaMeshReadHandler,
2897 areaMeshMasterInstance,
2898 (*areaMeshPtr)
2899 );
2900 }
2901
2902
2903 // Detect lagrangian fields
2904 if (decompose)
2905 {
2906 InfoOrPout
2907 << "Setting caseName to " << baseRunTime.caseName()
2908 << " to read lagrangian" << endl;
2909 if (UPstream::master())
2910 {
2911 // Change case name but only on the master - this will
2912 // hopefully cause the slaves to not read.
2913 runTime.caseName() = baseRunTime.caseName();
2914 }
2915 else
2916 {
2917 // Explicitly make sure that casename is not recognised as
2918 // a processor case since that has special handling for
2919 // caching processor directories etc.
2920 runTime.caseName() = "#invalid-name#";
2921 }
2922 runTime.processorCase(false);
2923 }
2924
2925 // Read lagrangian fields and store on cloud (objectRegistry)
2927 (
2929 (
2930 mesh,
2931 selectedLagrangianFields
2932 )
2933 );
2934
2935 if (decompose)
2936 {
2937 InfoOrPout<< "Restoring caseName" << endl;
2938 runTime.caseName() = proc0CaseName;
2939 runTime.processorCase(oldProcCase);
2940 }
2941
2942
2943 // Load fields, do all distribution (mesh and fields)
2944 // - but not lagrangian fields; these are done later
2945 autoPtr<mapDistributePolyMesh> distMap = redistributeAndWrite
2946 (
2947 volMeshReadHandler, // readHandler
2948 writeHandler, // writeHandler,
2949 baseRunTime,
2950 proc0CaseName,
2951
2952 // Controls
2953 true, // read fields
2954 decompose, // decompose, i.e. read from undecomposed case
2955 false, // no reconstruction
2956 overwrite,
2957
2958 // Decomposition information
2959 nDestProcs,
2960 finalDecomp,
2961
2962 // For finite volume
2963 volMeshOnProc,
2964 volMeshMasterInstance,
2965 mesh
2966 );
2967
2968 // Redistribute any clouds
2970 (
2971 lagrangianDistributorPtr,
2972 mesh,
2973 nOldCells,
2974 distMap(),
2975 clouds
2976 );
2977
2978
2979 // Redistribute area fields
2980
2981 mapDistributePolyMesh faDistMap;
2982 autoPtr<faMesh> areaProcMeshPtr;
2983
2984 if (areaMeshPtr)
2985 {
2987 (
2988 areaMeshPtr(),
2989 distMap(),
2990 areaProcMeshPtr
2991 );
2992
2993 // Force recreation of everything that might vaguely
2994 // be used by patches:
2995
2996 faMeshTools::forceDemandDriven(areaProcMeshPtr());
2997
2998
2999 if (reconstruct)
3000 {
3001 if (Pstream::master())
3002 {
3003 InfoOrPout
3004 << "Setting caseName to " << baseRunTime.caseName()
3005 << " to write reconstructed mesh (and fields)."
3006 << endl;
3007 runTime.caseName() = baseRunTime.caseName();
3008 const bool oldProcCase(runTime.processorCase(false));
3009 const bool oldParRun = UPstream::parRun(false);
3010 const label oldNumProcs(fileHandler().nProcs());
3011
3012 areaProcMeshPtr->write();
3013
3014 // Now we've written all. Reset caseName on master
3015 InfoOrPout<< "Restoring caseName" << endl;
3016 const_cast<fileOperation&>
3017 (
3018 fileHandler()
3019 ).nProcs(oldNumProcs);
3020 UPstream::parRun(oldParRun);
3021 runTime.caseName() = proc0CaseName;
3022 runTime.processorCase(oldProcCase);
3023 }
3024 }
3025 else
3026 {
3027 auto oldHandler = fileOperation::fileHandler(writeHandler);
3028
3030 (
3031 IOobject
3032 (
3033 "procAddressing",
3034 areaProcMeshPtr->facesInstance(),
3036 areaProcMeshPtr->thisDb()
3037 ),
3038 faDistMap
3039 );
3040
3041 areaProcMeshPtr->write();
3042
3043 writeHandler = fileOperation::fileHandler(oldHandler);
3044
3045 if (decompose)
3046 {
3048 (
3049 areaProcMeshPtr(),
3050 faDistMap,
3051 decompose,
3052 writeHandler
3053 );
3054 }
3055 }
3056
3057 InfoOrPout
3058 << "Written redistributed mesh to "
3059 << areaProcMeshPtr->facesInstance() << nl << endl;
3060
3061 faMeshDistributor distributor
3062 (
3063 areaMeshPtr(), // source
3064 areaProcMeshPtr(), // target
3065 faDistMap,
3066 writeHandler
3067 );
3068
3069 areaFields.redistributeAndWrite(distributor, true);
3070 }
3071
3072
3073 // Get reference to standard write handler
3074 refPtr<fileOperation> defaultHandler;
3075 if (writeHandler)
3076 {
3077 defaultHandler.ref(writeHandler.ref());
3078 }
3079 else
3080 {
3081 defaultHandler.ref(const_cast<fileOperation&>(fileHandler()));
3082 }
3083
3084
3085 copyUniform
3086 (
3087 volMeshReadHandler, // read handler
3088 defaultHandler, //TBD: should be all IOranks
3089
3090 reconstruct, // reconstruct
3091 decompose, // decompose
3092
3093 oldTimeName,
3094 (decompose ? baseRunTime.caseName() : proc0CaseName),
3095 mesh, // read location is mesh (but oldTimeName)
3096 mesh // write location is mesh
3097 );
3098 }
3099
3100
3101 // Get reference to standard write handler
3102 refPtr<fileOperation> defaultHandler;
3103 if (writeHandler)
3104 {
3105 defaultHandler.ref(writeHandler.ref());
3106 }
3107 else
3108 {
3109 defaultHandler.ref(const_cast<fileOperation&>(fileHandler()));
3110 }
3111
3112 copyUniform
3113 (
3114 volMeshReadHandler, // read handler
3115 defaultHandler, //TBD: should be all IOranks
3116
3117 reconstruct, // reconstruct (=false)
3118 decompose, // decompose
3119
3120 oldTimeName, // provided read time
3121 (decompose ? baseRunTime.caseName() : proc0CaseName),
3122 readRunTime,
3123 runTime // writing location
3124 );
3125 }
3126
3127
3128 InfoOrPout<< "End\n" << endl;
3129
3130 return 0;
3131}
3132
3133
3134// ************************************************************************* //
Inter-processor communication reduction functions.
Required Classes.
Info<< nl;Info<< "Write faMesh in vtk format:"<< nl;{ vtk::uindirectPatchWriter writer(aMesh.patch(), fileName(aMesh.time().globalPath()/vtkBaseFileName));writer.writeGeometry();globalIndex procAddr(aMesh.nFaces());labelList cellIDs;if(UPstream::master()) { cellIDs.resize(procAddr.totalSize());for(const labelRange &range :procAddr.ranges()) { auto slice=cellIDs.slice(range);slice=identity(range);} } writer.beginCellData(4);writer.writeProcIDs();writer.write("cellID", cellIDs);writer.write("area", aMesh.S().field());writer.write("normal", aMesh.faceAreaNormals());writer.beginPointData(1);writer.write("normal", aMesh.pointAreaNormals());Info<< " "<< writer.output().name()<< nl;}{ vtk::lineWriter writer(aMesh.points(), aMesh.edges(), fileName(aMesh.time().globalPath()/(vtkBaseFileName+"-edges")));writer.writeGeometry();writer.beginCellData(4);writer.writeProcIDs();{ Field< scalar > fld(faMeshTools::flattenEdgeField(aMesh.magLe(), true))
labelList faceLabels(nFaceLabels)
bool erase(iterator &iter)
Erase entry specified by given iterator and delete the allocated pointer.
bool insert(const Key &key)
Insert a new entry, not overwriting existing entries.
Definition HashSet.H:229
bool found(const Key &key) const
Same as contains().
Definition HashTable.H:1370
static void writeContents(const IOobject &io, const UList< T > &content)
Write contents. The IOobject is never registered.
Definition IOList.C:194
A IOmapDistributePolyMesh wrapper for using referenced external data.
static void writeContents(const IOobject &io, const mapDistributePolyMesh &map)
Write contents. The IOobject is never registered.
List of IOobjects with searching and retrieving facilities. Implemented as a HashTable,...
wordList names() const
The unsorted names of the IOobjects.
@ NO_REGISTER
Do not request registration (bool: false).
@ REGISTER
Request registration (bool: true).
@ NO_READ
Nothing to be read.
@ READ_IF_PRESENT
Reading is optional [identical to LAZY_READ].
@ MUST_READ
Reading required.
@ LAZY_READ
Reading is optional [identical to READ_IF_PRESENT].
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
@ AUTO_WRITE
Automatically write from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
static unsigned int minPrecision(unsigned int prec) noexcept
Set the minimum default precision.
Definition IOstream.H:459
void setSize(label n)
Alias for resize().
Definition List.H:536
static bool Store(std::unique_ptr< faMeshesRegistry > &&ptr)
static FOAM_NO_DANGLING_REFERENCE const pointMesh & New(const polyMesh &mesh, Args &&... args)
static std::unique_ptr< faMeshesRegistry > Release(const word &objName, const polyMesh &mesh, bool checkout=true)
static void allGatherList(UList< T > &values, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Gather data, but keep individual values separate. Uses MPI_Allgather or manual communication.
static void gatherList(UList< T > &values, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Gather data, but keep individual values separate.
static void broadcast(Type &value, const int communicator=UPstream::worldComm)
Broadcast content (contiguous or non-contiguous) to all communicator ranks. Does nothing in non-paral...
static void broadcasts(const int communicator, Type &value, Args &&... values)
Broadcast multiple items to all communicator ranks. Does nothing in non-parallel.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition PtrList.H:67
instantList times() const
Search the case for valid time directories.
Definition TimePaths.C:230
static instantList findTimes(const fileName &directory, const word &constantDirName="constant")
Search a given directory for valid time directories.
Definition TimePaths.C:109
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
Switch runTimeModifiable() const noexcept
Supports re-reading.
Definition Time.H:585
virtual void setTime(const Time &t)
Reset the time and time-index to those of the given time.
Definition Time.C:905
const fileName & caseName() const noexcept
The case name.
Definition TimePathsI.H:78
bool found(const T &val, label pos=0) const
Same as contains().
Definition UList.H:983
bool empty() const noexcept
True if List is empty (ie, size() is zero).
Definition UList.H:701
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
static label commWorld() noexcept
Communicator for all ranks (respecting any local worlds).
Definition UPstream.H:1101
static int myProcNo(const label communicator=worldComm)
Rank of this process in the communicator (starting from masterNo()). Negative if the process is not a...
Definition UPstream.H:1706
static bool parRun(const bool on) noexcept
Set as parallel run on/off.
Definition UPstream.H:1669
static constexpr int masterNo() noexcept
Relative rank for the master process - is always 0.
Definition UPstream.H:1691
static bool master(const label communicator=worldComm)
True if process corresponds to the master rank in the communicator.
Definition UPstream.H:1714
static label nProcs(const label communicator=worldComm)
Number of ranks in parallel run (for given communicator). It is 1 for serial run.
Definition UPstream.H:1697
static label worldComm
Communicator for all ranks. May differ from commGlobal() if local worlds are in use.
Definition UPstream.H:1069
static rangeType allProcs(const label communicator=worldComm)
Range of process indices for all processes.
Definition UPstream.H:1857
Finite area area (element) fields.
Mesh data needed to do the Finite Area discretisation.
Definition areaFaMesh.H:50
Extract command arguments and options from the supplied argc and argv parameters.
Definition argList.H:119
static void noFunctionObjects(bool addWithOption=false)
Remove '-noFunctionObjects' option and ignore any occurrences.
Definition argList.C:562
static void addVerboseOption(const string &usage="", bool advanced=false)
Enable a 'verbose' bool option, with usage information.
Definition argList.C:535
static void addDryRunOption(const string &usage, bool advanced=false)
Enable a 'dry-run' bool option, with usage information.
Definition argList.C:519
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 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
void clear() noexcept
Same as reset(nullptr).
Definition autoPtr.H:255
T * get() noexcept
Return pointer to managed object without nullptr checking.
Definition autoPtr.H:216
void reset(T *p=nullptr) noexcept
Delete managed object and set to new given pointer.
Definition autoPtrI.H:37
static autoPtr< T > New(Args &&... args)
Construct autoPtr with forwarding arguments.
Definition autoPtr.H:178
Abstract base class for domain decomposition.
static label nDomains(const dictionary &decompDict, const word &regionName="")
Return region-specific or top-level numberOfSubdomains entry.
virtual labelList decompose(const pointField &points, const scalarField &pointWeights=scalarField::null()) const
Return the wanted processor number for every coordinate, using uniform or specified point weights.
virtual bool parallelAware() const =0
Is method parallel aware?
MeshObject wrapper of decompositionMethod.
static const decompositionModel & New(const polyMesh &mesh, const fileName &decompDictFile="", const dictionary *fallback=nullptr)
Read and register on mesh, optionally with alternative decomposeParDict path/name or with fallback co...
decompositionMethod & decomposer() const
Return demand-driven decomposition method.
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T. FatalIOError if not found, or if the number of tokens is incorrect.
bool found(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find an entry (const access) with the given keyword.
Holds a reference to the original mesh (the baseMesh) and optionally to a subset of that mesh (the su...
static int verbose_
Output verbosity when writing.
static mapDistributePolyMesh distribute(const faMesh &oldMesh, const mapDistributePolyMesh &distMap, const polyMesh &tgtPolyMesh, autoPtr< faMesh > &newMeshPtr)
Distribute mesh according to the given (volume) mesh distribution.
static void forceDemandDriven(faMesh &mesh)
Force creation of everything that might vaguely be used by patches.
static mapDistributePolyMesh readProcAddressing(const faMesh &mesh, const faMesh *baseMeshPtr)
Read decompose/reconstruct addressing.
static void writeProcAddressing(const faMesh &mesh, const mapDistributePolyMesh &faDistMap, const bool decompose, refPtr< fileOperation > &writeHandler, const faMesh *procMesh=nullptr)
Write decompose/reconstruct addressing.
static void unregisterMesh(const faMesh &mesh)
Unregister the faMesh from its associated polyMesh to prevent triggering on polyMesh changes etc.
static autoPtr< faMesh > loadOrCreateMesh(const word &areaName, const IOobject &io, const polyMesh &pMesh, const bool decompose, const bool verbose=false)
Finite area mesh (used for 2-D non-Euclidian finite area method) defined using a patch of faces on a ...
Definition faMesh.H:140
fileName meshDir() const
Return the local mesh directory (dbDir()/meshSubDir).
Definition faMesh.C:1022
static word meshSubDir
The mesh sub-directory name (usually "faMesh").
Definition faMesh.H:750
A class for handling file names.
Definition fileName.H:75
An encapsulation of filesystem-related operations.
static const fileOperation & fileHandler()
Return the current file handler. Will create the default file handler if necessary.
static autoPtr< fileOperation > NewUncollated()
The commonly used uncollatedFileOperation.
static int cacheLevel() noexcept
Return cache level.
static word processorsBaseDir
Return the processors directory name (usually "processors").
static autoPtr< fileOperation > New(const word &handlerType, bool verbose=false)
Select fileHandler-type. Uses defaultFileHandler if the handlerType is empty.
static int nProcsFilter() noexcept
Return collated 'processorsDDD' filtering.
Sends/receives parts of mesh+fvfields to neighbouring processors. Used in load balancing.
Holds a reference to the original mesh (the baseMesh) and optionally to a subset of that mesh (the su...
static autoPtr< fvMesh > newMesh(const IOobject &io, const bool masterOnlyReading, const bool verbose=false)
Read mesh or create dummy mesh (0 cells, >0 patches).
static autoPtr< mapDistributePolyMesh > readProcAddressing(const fvMesh &procMesh, const fvMesh *baseMeshPtr)
Read procAddressing components (reconstructing).
static void setBasicGeometry(fvMesh &mesh)
Set the fvGeometryScheme to basic (to avoid parallel communication).
static void writeProcAddressing(const fvMesh &procMesh, const mapDistributePolyMesh &map, const bool decompose, const fileName &writeInstance, refPtr< fileOperation > &writeHandler)
Write addressing if decomposing (1 to many) or reconstructing (many to 1).
static autoPtr< fvMesh > loadOrCreateMesh(const IOobject &io, const bool decompose, const bool verbose=false)
Read mesh if available, or create empty mesh with non-proc as per proc0 mesh.
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
static const word & zeroGradientType() noexcept
The type name for zeroGradient patch fields.
Calculates a non-overlapping list of offsets based on an input size (eg, number of cells) from differ...
Definition globalIndex.H:77
Various for reading/decomposing/reconstructing/distributing refinement data.
Definition hexRef8Data.H:57
An instant of time. Contains the time value and name. Uses Foam::Time when formatting the name.
Definition instant.H:56
Class containing mesh-to-mesh mapping information after a mesh distribution where we send parts of me...
void clear()
Reset to zero size, only retaining communicator(s).
static void removeFiles(const polyMesh &)
Helper: remove all relevant files from mesh instance.
Registry of regIOobjects.
const Time & time() const noexcept
Return time registry.
Simple container to manage read/write, redistribute finiteArea fields.
Finite volume reconstructor for volume and surface fields.
Distributor/redistributor for point fields, uses a two (or three) stage construction.
static int verbose_
Output verbosity when writing.
Mesh representing a set of points created from polyMesh.
Definition pointMesh.H:49
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
static const word & regionName(const word &region)
The mesh region name or word::null if polyMesh::defaultRegion.
Definition polyMesh.C:796
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
readUpdateState
Enumeration defining the state of the mesh after a read update.
Definition polyMesh.H:92
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh").
Definition polyMesh.H:411
Neighbour processor patch.
int neighbProcNo() const noexcept
Return neighbour processor number.
A class for managing references or pointers (no reference counting).
Definition refPtr.H:54
T * get() noexcept
Return pointer without nullptr checking.
Definition refPtr.H:249
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition refPtrI.H:230
static void unset(bool verbose=false)
Deactivate SIGFPE handler and NaN memory initialisation.
static void addOptions(const bool constant=true, const bool withZero=false)
Add timeSelector options to argList::validOptions.
static instantList selectIfPresent(Time &runTime, const argList &args)
If any time option provided return the set of times - as per select0() - otherwise return just the cu...
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 List of wordRe with additional matching capabilities.
Definition wordRes.H:56
A class for handling words, derived from Foam::string.
Definition word.H:66
static const word null
An empty word.
Definition word.H:84
A class representing the concept of 0 (zero) that can be used to avoid manipulating objects known to ...
Definition zero.H:58
fileName path(UMean.rootPath()/UMean.caseName()/"graphs"/UMean.instance())
dynamicFvMesh & mesh
engineTime & runTime
Foam::word regionName(args.getOrDefault< word >("region", Foam::polyMesh::defaultRegion))
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
const auto & io
label nPointFields
auto & name
wordList regionNames
Miscellaneous file handling for meshes.
#define WarningInFunction
Report a warning using Foam::Warning.
Namespace for handling debugging switches.
Definition debug.C:45
int debugSwitch(const char *name, const int deflt=0)
Lookup debug switch or add default value.
Definition debug.C:222
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > reconstruct(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
bool allowFaModels() noexcept
The enable/disable state for regionFaModel (default: true).
Namespace for OpenFOAM.
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
GeometricField< vector, fvPatchField, volMesh > volVectorField
const dimensionSet dimless
Dimensionless.
HashSet< word, Hash< word > > wordHashSet
A HashSet of words, uses string hasher.
Definition HashSet.H:80
List< labelList > labelListList
List of labelList.
Definition labelList.H:38
List< label > labelList
A List of labels.
Definition List.H:62
bool mkDir(const fileName &pathName, mode_t mode=0777)
Make a directory and return an error if it could not be created.
Definition POSIX.C:616
PtrList< unmappedPassivePositionParticleCloud > readLagrangian(const fvMesh &mesh, const UList< word > &cloudNames, const boolUList &haveClouds, const wordRes &selectedFields)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
IOList< label > labelIOList
IO for a List of label.
Definition labelIOList.H:32
List< fileName > fileNameList
List of fileName.
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
void reconstructLagrangian(autoPtr< parLagrangianDistributor > &distributorPtr, const fvMesh &baseMesh, const fvMesh &mesh, const mapDistributePolyMesh &distMap, const wordRes &selectedFields)
messageStream Info
Information stream (stdout output on master, null elsewhere).
List< instant > instantList
List of instants.
Definition instantList.H:41
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
GeometricField< tensor, fvPatchField, volMesh > volTensorField
void masterMeshInstance(const IOobject &io, fileName &facesInstance, fileName &pointsInstance)
Determine master faces instance.
GeometricField< sphericalTensor, fvPatchField, volMesh > volSphericalTensorField
bool returnReduceAnd(const bool value, const int communicator=UPstream::worldComm)
Perform logical (and) MPI Allreduce on a copy. Uses UPstream::reduceAnd.
void removeProcAddressing(const faMesh &mesh)
Remove procAddressing.
boolList haveMeshFile(const Time &runTime, const fileName &meshPath, const word &meshFile="faces", const bool verbose=true)
Check for availability of specified mesh file (default: "faces").
GeometricField< symmTensor, fvPatchField, volMesh > volSymmTensorField
List< bool > boolList
A List of bools.
Definition List.H:60
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
UList< bool > boolUList
A UList of bools.
Definition UList.H:73
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
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
labelList findIndices(const ListType &input, typename ListType::const_reference val, label start=0)
Linear search to find all occurrences of given element.
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
bool isDir(const fileName &name, const bool followLink=true)
Does the name exist as a DIRECTORY in the file system?
Definition POSIX.C:862
void redistributeLagrangian(autoPtr< parLagrangianDistributor > &distributorPtr, const fvMesh &mesh, const label nOldCells, const mapDistributePolyMesh &distMap, PtrList< unmappedPassivePositionParticleCloud > &clouds)
messageStream Warning
Warning stream (stdout output on master, null elsewhere), with additional 'FOAM Warning' header text.
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
Reading, reconstruct, redistribution of lagrangian fields.
Foam::argList args(argc, argv)
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
Foam::surfaceFields.