Loading...
Searching...
No Matches
decomposePar.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) 2016-2022,2024 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 decomposePar
29
30Group
31 grpParallelUtilities
32
33Description
34 Automatically decomposes a mesh and fields of a case for parallel
35 execution of OpenFOAM.
36
37Usage
38 \b decomposePar [OPTIONS]
39
40 Options:
41 - \par -allRegions
42 Decompose all regions in regionProperties. Does not check for
43 existence of processor*.
44
45 - \par -case <dir>
46 Specify case directory to use (instead of the cwd).
47
48 - \par -cellDist
49 Write the cell distribution as a labelList, for use with 'manual'
50 decomposition method and as a VTK or volScalarField for visualization.
51
52 - \par -constant
53 Include the 'constant/' dir in the times list.
54
55 - \par -copyUniform
56 Copy any \a uniform directories too.
57
58 - \par -copyZero
59 Copy \a 0 directory to processor* rather than decompose the fields.
60
61 - \par -debug-switch <name=val>
62 Specify the value of a registered debug switch. Default is 1
63 if the value is omitted. (Can be used multiple times)
64
65 - \par -decomposeParDict <file>
66 Use specified file for decomposePar dictionary.
67
68 - \par -dry-run
69 Test without writing the decomposition. Changes -cellDist to
70 only write VTK output.
71
72 - \par -fields
73 Use existing geometry decomposition and convert fields only.
74
75 - \par fileHandler <handler>
76 Override the file handler type.
77
78 - \par -force
79 Remove any existing \a processor subdirectories before decomposing the
80 geometry.
81
82 - \par -ifRequired
83 Only decompose the geometry if the number of domains has changed from a
84 previous decomposition. No \a processor subdirectories will be removed
85 unless the \a -force option is also specified. This option can be used
86 to avoid redundant geometry decomposition (eg, in scripts), but should
87 be used with caution when the underlying (serial) geometry or the
88 decomposition method etc. have been changed between decompositions.
89
90 - \par -info-switch <name=val>
91 Specify the value of a registered info switch. Default is 1
92 if the value is omitted. (Can be used multiple times)
93
94 - \par -latestTime
95 Select the latest time.
96
97 - \par -lib <name>
98 Additional library or library list to load (can be used multiple times).
99
100 - \par -noFunctionObjects
101 Do not execute function objects.
102
103 - \par -noSets
104 Skip decomposing cellSets, faceSets, pointSets.
105
106 - \par -noZero
107 Exclude the \a 0 dir from the times list.
108
109 - \par -opt-switch <name=val>
110 Specify the value of a registered optimisation switch (int/bool).
111 Default is 1 if the value is omitted. (Can be used multiple times)
112
113 - \par -region <regionName>
114 Decompose named region. Does not check for existence of processor*.
115
116 - \par -time <ranges>
117 Override controlDict settings and decompose selected times. Does not
118 re-decompose the mesh i.e. does not handle moving mesh or changing
119 mesh cases. Eg, ':10,20 40:70 1000:', 'none', etc.
120
121 - \par -verbose
122 Additional verbosity.
123
124 - \par -doc
125 Display documentation in browser.
126
127 - \par -doc-source
128 Display source code in browser.
129
130 - \par -help
131 Display short help and exit.
132
133 - \par -help-man
134 Display full help (manpage format) and exit.
135
136 - \par -help-notes
137 Display help notes (description) and exit.
138
139 - \par -help-full
140 Display full help and exit.
141
142\*---------------------------------------------------------------------------*/
143
144#include "argList.H"
145#include "timeSelector.H"
146#include "OSspecific.H"
147#include "IOobjectList.H"
148
149#include "decompositionModel.H"
150#include "domainDecomposition.H"
152
153#include "regionProperties.H"
154
155#include "fieldsDistributor.H"
156
157#include "fvFieldDecomposer.H"
158#include "pointFields.H"
159#include "pointFieldDecomposer.H"
160
162
163#include "emptyFaPatch.H"
164#include "faFieldDecomposer.H"
165#include "faMeshDecomposition.H"
166
167// * * * * * * * * * * * * * * * Local Functions * * * * * * * * * * * * * * //
168
169namespace Foam
170{
171
172// Read proc addressing at specific instance.
173// Uses polyMesh/fvMesh meshSubDir by default
174autoPtr<labelIOList> procAddressing
175(
176 const objectRegistry& procRegistry,
177 const word& name,
178 const word& instance,
180)
181{
183 (
185 (
186 name,
187 instance,
188 local,
189 procRegistry,
193 )
194 );
195}
196
197
198// Read proc addressing at specific instance.
199// Uses the finiteArea meshSubDir
200autoPtr<labelIOList> faProcAddressing
201(
202 const objectRegistry& procRegistry,
203 const word& name,
204 const word& instance,
206)
207{
208 return procAddressing(procRegistry, name, instance, local);
209}
210
211
212// Return cached or read proc addressing from facesInstance
214const labelIOList& procAddressing
215(
216 const PtrList<fvMesh>& procMeshList,
217 const label proci,
218 const word& name,
219 PtrList<labelIOList>& procAddressingList
220)
221{
222 const auto& procMesh = procMeshList[proci];
223
224 return procAddressingList.try_emplace
225 (
226 proci,
228 (
229 name,
230 procMesh.facesInstance(),
232 procMesh,
236 )
237 );
238}
239
240
241void decomposeUniform
242(
243 const bool copyUniform,
245 const Time& processorDb,
246 const word& regionDir = word::null
247)
248{
249 const Time& runTime = mesh.time();
250
251 // Any uniform data to copy/link?
252 const fileName uniformDir(regionDir/"uniform");
253
254 if (fileHandler().isDir(runTime.timePath()/uniformDir))
255 {
256 Info<< "Detected additional non-decomposed files in "
257 << runTime.timePath()/uniformDir
258 << endl;
259
260 // Bit of trickery to synthesise the correct directory base,
261 // e.g. processors4/0.01
262 const fileName timePath = fileHandler().objectPath
263 (
265 (
266 "dummy",
268 processorDb
269 ),
271 ).path();
272
273 // If no fields have been decomposed the destination
274 // directory will not have been created so make sure.
275 mkDir(timePath);
276
277 if (copyUniform || mesh.distributed())
278 {
279 if (!fileHandler().exists(timePath/uniformDir))
280 {
281 fileHandler().cp
282 (
283 runTime.timePath()/uniformDir,
284 timePath/uniformDir
285 );
286 }
287 }
288 else
289 {
290 // Link with relative paths
291 string parentPath = string("..")/"..";
292
293 if (!regionDir.empty())
294 {
295 parentPath = parentPath/"..";
296 }
297
298 fileName currentDir(cwd());
299 chDir(timePath);
300
301 if (!fileHandler().exists(uniformDir))
302 {
303 fileHandler().ln
304 (
305 parentPath/runTime.timeName()/uniformDir,
306 uniformDir
307 );
308 }
309 chDir(currentDir);
310 }
311 }
312}
313
314} // End namespace Foam
315
316
317using namespace Foam;
318
319// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
320
321int main(int argc, char *argv[])
322{
324 (
325 "Decompose a mesh and fields of a case for parallel execution"
326 );
327
330 (
331 "decomposeParDict",
332 "file",
333 "Alternative decomposePar dictionary file"
334 );
335
336 #include "addAllRegionOptions.H"
337 #include "addAllFaRegionOptions.H"
338
340 (
341 "Test without writing the decomposition. "
342 "Changes -cellDist to only write VTK output."
343 );
346 (
347 "domains",
348 "N",
349 "Override numberOfSubdomains (-dry-run only)",
350 true // Advanced option
351 );
353 (
354 "method",
355 "name",
356 "Override decomposition method (-dry-run only)",
357 true // Advanced option
358 );
359
361 (
362 "no-finite-area",
363 "Suppress finiteArea mesh/field decomposition",
364 true // Advanced option
365 );
366
368 (
369 "no-lagrangian",
370 "Suppress lagrangian (cloud) decomposition",
371 true // Advanced option
372 );
373
375 (
376 "cellDist",
377 "Write cell distribution as a labelList - for use with 'manual' "
378 "decomposition method and as a volScalarField for visualization."
379 );
381 (
382 "copyZero",
383 "Copy 0/ directory to processor*/ rather than decompose the fields"
384 );
386 (
387 "copyUniform",
388 "Copy any uniform/ directories too"
389 );
391 (
392 "fields",
393 "Use existing geometry decomposition and convert fields only"
394 );
396 (
397 "no-fields",
398 "Suppress conversion of fields (volume, finite-area, lagrangian)"
399 );
400
402 (
403 "no-sets",
404 "Skip decomposing cellSets, faceSets, pointSets"
405 );
406 argList::addOptionCompat("no-sets", {"noSets", 2106});
407
409 (
410 "force",
411 "Remove existing processor*/ subdirs before decomposing the geometry"
412 );
414 (
415 "ifRequired",
416 "Only decompose geometry if the number of domains has changed"
417 );
418
419 // Allow explicit -constant, have zero from time range
420 timeSelector::addOptions(true, false); // constant(true), zero(false)
421
422 // Prevent volume BCs from triggering finite-area
424
425 #include "setRootCase.H"
426
427 // ------------------------------------------------------------------------
428 // Configuration
429
430 const bool writeCellDist = args.found("cellDist");
431
432 // Most of these are ignored for dry-run (not triggered anywhere)
433 const bool copyZero = args.found("copyZero");
434 const bool copyUniform = args.found("copyUniform");
435 const bool decomposeSets = !args.found("no-sets");
436
437 const bool decomposeIfRequired = args.found("ifRequired");
438
439 const bool doDecompFields = !args.found("no-fields");
440 const bool doFiniteArea = !args.found("no-finite-area");
441 const bool doLagrangian = !args.found("no-lagrangian");
442
443 bool decomposeFieldsOnly = args.found("fields");
444 bool forceOverwrite = args.found("force");
445
446 // Set time from database
447 #include "createTime.H"
448
449 // Allow override of time (unless dry-run)
450 instantList times;
451 if (args.dryRun())
452 {
453 Info<< "\ndry-run: ignoring -copy*, -fields, -force, time selection"
454 << nl;
455 }
456 else
457 {
458 if (decomposeFieldsOnly && !doDecompFields)
459 {
460 FatalErrorIn(args.executable())
461 << "Options -fields and -no-fields are mutually exclusive"
462 << " ... giving up" << nl
463 << exit(FatalError);
464 }
465
466 if (!doDecompFields)
467 {
468 Info<< "Skip decompose of all fields" << nl;
469 }
470 if (!doFiniteArea)
471 {
472 Info<< "Skip decompose of finiteArea mesh/fields" << nl;
473 }
474 if (!doLagrangian)
475 {
476 Info<< "Skip decompose of lagrangian positions/fields" << nl;
477 }
478
480 }
481
482
483 // Allow override of decomposeParDict location
484 fileName decompDictFile;
485 if
486 (
487 args.readIfPresent("decomposeParDict", decompDictFile)
488 && !decompDictFile.empty() && !decompDictFile.isAbsolute()
489 )
490 {
491 decompDictFile = runTime.globalPath()/decompDictFile;
492 }
493
494 // Handle volume region selections
495 #include "getAllRegionOptions.H"
496
497 // Handle area region selections
498 #include "getAllFaRegionOptions.H"
499
500 if (!doFiniteArea)
501 {
502 areaRegionNames.clear(); // For consistency
503 }
504
505 const bool optRegions =
506 (regionNames.size() != 1 || regionNames[0] != polyMesh::defaultRegion);
507
508 if (regionNames.size() == 1 && regionNames[0] != polyMesh::defaultRegion)
509 {
510 Info<< "Using region: " << regionNames[0] << nl << endl;
511 }
512
513 forAll(regionNames, regioni)
514 {
515 const word& regionName = regionNames[regioni];
517
518 if (args.dryRun())
519 {
520 Info<< "dry-run: decomposing mesh " << regionName << nl << nl
521 << "Create mesh..." << flush;
522
524 (
526 (
528 runTime.timeName(),
529 runTime,
533 ),
534 decompDictFile,
535 args.getOrDefault<label>("domains", 0),
536 args.getOrDefault<word>("method", word::null)
537 );
538
539 decompTest.execute(writeCellDist, args.verbose());
540 continue;
541 }
542
543 Info<< "\n\nDecomposing mesh";
544 if (!regionDir.empty())
545 {
546 Info<< ' ' << regionName;
547 }
548 Info<< nl << endl;
549
550 // Determine the existing processor count directly
551 const label nProcsOld =
552 fileHandler().nProcs(runTime.path(), regionDir);
553
554 // Get requested numberOfSubdomains directly from the dictionary.
555 // Note: have no mesh yet so cannot use decompositionModel::New
556 const label nDomains = decompositionMethod::nDomains
557 (
559 (
561 (
563 (
565 runTime.time().system(),
566 regionDir, // region (if non-default)
567 runTime,
571 ),
572 decompDictFile
573 )
574 )
575 );
576
577 // Give file handler a chance to determine the output directory
578 const_cast<fileOperation&>(fileHandler()).nProcs(nDomains);
579
580 if (decomposeFieldsOnly)
581 {
582 // Sanity check on previously decomposed case
583 if (nProcsOld != nDomains)
584 {
585 FatalErrorIn(args.executable())
586 << "Specified -fields, but the case was decomposed with "
587 << nProcsOld << " domains"
588 << nl
589 << "instead of " << nDomains
590 << " domains as specified in decomposeParDict" << nl
591 << exit(FatalError);
592 }
593 }
594 else if (nProcsOld)
595 {
596 bool procDirsProblem = true;
597
598 if (decomposeIfRequired && nProcsOld == nDomains)
599 {
600 // We can reuse the decomposition
601 decomposeFieldsOnly = true;
602 procDirsProblem = false;
603 forceOverwrite = false;
604
605 Info<< "Using existing processor directories" << nl;
606 }
607
608 if (optRegions)
609 {
610 procDirsProblem = false;
611 forceOverwrite = false;
612 }
613
614 if (forceOverwrite)
615 {
616 Info<< "Removing " << nProcsOld
617 << " existing processor directories" << endl;
618
619 // Remove existing processors directory
620 fileNameList dirs
621 (
623 (
624 runTime.path(),
626 )
627 );
628 forAllReverse(dirs, diri)
629 {
630 const fileName& d = dirs[diri];
631
632 label proci = -1;
633
634 if
635 (
636 d.starts_with("processor")
637 &&
638 (
639 // Collated is "processors"
640 d[9] == 's'
641
642 // Uncollated has integer(s) after 'processor'
643 || Foam::read(d.substr(9), proci)
644 )
645 )
646 {
647 if (fileHandler().exists(d))
648 {
649 fileHandler().rmDir(d);
650 }
651 }
652 }
653
654 procDirsProblem = false;
655 }
656
657 if (procDirsProblem)
658 {
659 FatalErrorIn(args.executable())
660 << "Case is already decomposed with " << nProcsOld
661 << " domains, use the -force option or manually" << nl
662 << "remove processor directories before decomposing. e.g.,"
663 << nl
664 << " rm -rf " << runTime.path().c_str() << "/processor*"
665 << nl
666 << exit(FatalError);
667 }
668 }
669
670 Info<< "Create mesh" << endl;
672 (
674 (
676 runTime.timeName(),
677 runTime,
681 ),
682 decompDictFile
683 );
684 // Make sure pointMesh gets read as well
686
687
688 // Decompose the mesh
689 if (!decomposeFieldsOnly)
690 {
691 mesh.decomposeMesh();
692 mesh.writeDecomposition(decomposeSets);
693
694 if (writeCellDist)
695 {
696 const labelList& procIds = mesh.cellToProc();
697
698 // Write decomposition for visualization
699 mesh.writeVolField("cellDist");
700 //TBD: mesh.writeVTK("cellDist");
701
702 // Write decomposition as labelList for use with 'manual'
703 // decomposition method.
704 labelIOList cellDecomposition
705 (
707 (
708 "cellDecomposition",
709 mesh.facesInstance(),
710 mesh,
714 ),
715 procIds
716 );
717 cellDecomposition.write();
718
719 Info<< nl << "Wrote decomposition to "
720 << cellDecomposition.objectRelPath()
721 << " for use in manual decomposition." << endl;
722 }
723
724 fileHandler().flush();
725 }
726
727
728 if (copyZero)
729 {
730 // Copy the 0/ directory into each of the processor directories
731 // with fallback of 0.orig/ directory if necessary.
732
733 fileName inputDir(runTime.path()/"0");
734
735 bool canCopy = fileHandler().isDir(inputDir);
736 if (!canCopy)
737 {
738 // Try with "0.orig" instead
739 inputDir.ext("orig");
740 canCopy = fileHandler().isDir(inputDir);
741 }
742
743 if (canCopy)
744 {
745 // Avoid copying into the same directory multiple times
746 // (collated format). Don't need a hash here.
747 fileName prevOutputDir;
748 for (label proci = 0; proci < mesh.nProcs(); ++proci)
749 {
750 Time processorDb
751 (
753 args.rootPath(),
754 args.caseName()/("processor" + Foam::name(proci)),
755 false, // No function objects
756 false // No extra controlDict libs
757 );
758 // processorDb.setTime(runTime);
759
760 // Get corresponding directory name
761 // (to handle processors/)
762 const fileName outputDir
763 (
764 fileHandler().objectPath
765 (
767 (
768 word::null, // name
769 "0", // instance (time == 0)
770 processorDb
771 ),
773 )
774 );
775
776 if (outputDir != prevOutputDir)
777 {
778 Info<< "Processor " << proci
779 << ": copying \""
780 << inputDir.name() << "/\" to "
781 << runTime.relativePath(outputDir)
782 << endl;
783
784 fileHandler().cp(inputDir, outputDir);
785 prevOutputDir = outputDir;
786 }
787 }
788 }
789 else
790 {
791 Info<< "No 0/ or 0.orig/ directory to copy" << nl;
792 }
793 }
794 else
795 {
796 // Decompose field files, lagrangian, finite-area
797
798 // Cached processor meshes and maps. These are only preserved if
799 // running with multiple times.
800 PtrList<Time> processorDbList(mesh.nProcs());
801 PtrList<fvMesh> procMeshList(mesh.nProcs());
802 PtrList<labelIOList> faceProcAddressingList(mesh.nProcs());
803 PtrList<labelIOList> cellProcAddressingList(mesh.nProcs());
804 PtrList<labelIOList> boundaryProcAddressingList(mesh.nProcs());
805 PtrList<labelIOList> pointProcAddressingList(mesh.nProcs());
806 PtrList<labelIOList> pointBoundaryProcAddressingList(mesh.nProcs());
807
808 PtrList<fvFieldDecomposer> fieldDecomposerList(mesh.nProcs());
809 PtrList<pointFieldDecomposer> pointFieldDecomposerList
810 (
811 mesh.nProcs()
812 );
813
814
815 // Loop over all times
816 forAll(times, timei)
817 {
818 runTime.setTime(times[timei], timei);
819
820 Info<< "Time = " << runTime.timeName() << endl;
821
822 // Field objects at this time
823 IOobjectList objects;
824
825 // faMesh fields - can have multiple finite-area per volume
826 HashTable<IOobjectList> faObjects;
827
828 if (doDecompFields)
829 {
830 // List of volume mesh objects for this instance
831 objects = IOobjectList(mesh, runTime.timeName());
832
833 // Ignore generated fields: (cellDist)
834 objects.remove("cellDist");
835
836 // Lists of finite-area fields
837 faObjects.reserve(areaRegionNames.size());
838
839 for (const word& areaName : areaRegionNames)
840 {
841 // The finite-area objects for this area region
842 IOobjectList objs
843 (
845 runTime.timeName(),
846 polyMesh::regionName(areaName),
848 );
849
850 if (!objs.empty())
851 {
852 faObjects.emplace_set(areaName, std::move(objs));
853 }
854 }
855 }
856
857 // Finite area handling
858 // - all area regions use the same volume decomposition
859
860 HashPtrTable<faMeshDecomposition> faMeshDecompHashes;
861 if (doFiniteArea)
862 {
863 const word boundaryInst =
864 mesh.time().findInstance(mesh.meshDir(), "boundary");
865
866 for (const word& areaName : areaRegionNames)
867 {
869 (
870 "faBoundary",
871 boundaryInst,
872 faMesh::meshDir(mesh, areaName),
873 mesh.time(),
877 );
878
879 if (io.typeHeaderOk<faBoundaryMesh>(true))
880 {
881 // Always based on the volume decomposition!
882 faMeshDecompHashes.set
883 (
884 areaName,
886 (
887 areaName,
888 mesh,
889 mesh.nProcs(),
890 mesh.model()
891 )
892 );
893 }
894 }
895 }
896
897
898 // Volume/surface/internal fields
899 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
900
901 fvFieldDecomposer::fieldsCache volumeFieldCache;
902
903 if (doDecompFields)
904 {
905 volumeFieldCache.readAllFields(mesh, objects);
906 }
907
908
909 // Point fields
910 // ~~~~~~~~~~~~
911
912 // Read decomposed pointMesh
913 const pointMesh& pMesh =
915
916 pointFieldDecomposer::fieldsCache pointFieldCache;
917
918 if (doDecompFields)
919 {
920 pointFieldCache.readAllFields(pMesh, objects);
921 }
922
923
924 // Lagrangian fields
925 // ~~~~~~~~~~~~~~~~~
926
927 fileNameList cloudDirs;
928
929 if (doDecompFields && doLagrangian)
930 {
931 cloudDirs = fileHandler().readDir
932 (
933 runTime.timePath()/cloud::prefix,
935 );
936 }
937
938 // Particles
939 PtrList<Cloud<indexedParticle>> lagrangianPositions
940 (
941 cloudDirs.size()
942 );
943 // Particles per cell
945 (
946 cloudDirs.size()
947 );
948
949 lagrangianFieldDecomposer::fieldsCache lagrangianFieldCache
950 (
951 cloudDirs.size()
952 );
953
954 label cloudI = 0;
955
956 for (const fileName& cloudDir : cloudDirs)
957 {
958 IOobjectList cloudObjects
959 (
960 mesh,
961 runTime.timeName(),
962 cloud::prefix/cloudDir,
964 );
965
966 // Note: look up "positions" for backwards compatibility
967 if
968 (
969 cloudObjects.found("coordinates")
970 || cloudObjects.found("positions")
971 )
972 {
973 // Read lagrangian particles
974 // ~~~~~~~~~~~~~~~~~~~~~~~~~
975
976 Info<< "Identified lagrangian data set: "
977 << cloudDir << endl;
978
979 lagrangianPositions.set
980 (
981 cloudI,
983 (
984 mesh,
985 cloudDir,
986 false
987 )
988 );
989
990
991 // Sort particles per cell
992 // ~~~~~~~~~~~~~~~~~~~~~~~
993
994 cellParticles.set
995 (
996 cloudI,
997 new List<SLList<indexedParticle*>*>
998 (
999 mesh.nCells(),
1000 static_cast<SLList<indexedParticle*>*>(nullptr)
1001 )
1002 );
1003
1004 label i = 0;
1005
1006 for (indexedParticle& p : lagrangianPositions[cloudI])
1007 {
1008 p.index() = i++;
1009
1010 label celli = p.cell();
1011
1012 // Check
1013 if (celli < 0 || celli >= mesh.nCells())
1014 {
1015 FatalErrorIn(args.executable())
1016 << "Illegal cell number " << celli
1017 << " for particle with index "
1018 << p.index()
1019 << " at position "
1020 << p.position() << nl
1021 << "Cell number should be between 0 and "
1022 << mesh.nCells()-1 << nl
1023 << "On this mesh the particle should"
1024 << " be in cell "
1025 << mesh.findCell(p.position())
1026 << exit(FatalError);
1027 }
1028
1029 if (!cellParticles[cloudI][celli])
1030 {
1031 cellParticles[cloudI][celli] =
1033 }
1034
1035 cellParticles[cloudI][celli]->append(&p);
1036 }
1037
1038 // Read fields
1039 // ~~~~~~~~~~~
1040
1041 IOobjectList lagrangianObjects
1042 (
1043 mesh,
1044 runTime.timeName(),
1045 cloud::prefix/cloudDirs[cloudI],
1047 );
1048
1049 lagrangianFieldCache.readAllFields
1050 (
1051 cloudI,
1052 lagrangianObjects
1053 );
1054
1055 ++cloudI;
1056 }
1057 }
1058
1059 lagrangianPositions.resize(cloudI);
1060 cellParticles.resize(cloudI);
1061 lagrangianFieldCache.resize(cloudI);
1062
1063 Info<< endl;
1064
1065 // split the fields over processors
1066 for
1067 (
1068 label proci = 0;
1069 doDecompFields && proci < mesh.nProcs();
1070 ++proci
1071 )
1072 {
1073 Info<< "Processor " << proci << ": field transfer" << endl;
1074
1075 // open the database
1076 if (!processorDbList.set(proci))
1077 {
1078 processorDbList.set
1079 (
1080 proci,
1081 new Time
1082 (
1084 args.rootPath(),
1085 args.caseName()
1086 / ("processor" + Foam::name(proci)),
1087 args.allowFunctionObjects(),
1088 args.allowLibs()
1089 )
1090 );
1091 }
1092 Time& processorDb = processorDbList[proci];
1093
1094
1095 processorDb.setTime(runTime);
1096
1097 // Read the mesh
1098 if (!procMeshList.set(proci))
1099 {
1100 procMeshList.set
1101 (
1102 proci,
1103 new fvMesh
1104 (
1105 IOobject
1106 (
1107 regionName,
1108 processorDb.timeName(),
1109 processorDb
1110 )
1111 )
1112 );
1113 }
1114 const fvMesh& procMesh = procMeshList[proci];
1115
1116 const labelIOList& faceProcAddressing = procAddressing
1117 (
1118 procMeshList,
1119 proci,
1120 "faceProcAddressing",
1121 faceProcAddressingList
1122 );
1123
1124 const labelIOList& cellProcAddressing = procAddressing
1125 (
1126 procMeshList,
1127 proci,
1128 "cellProcAddressing",
1129 cellProcAddressingList
1130 );
1131
1132 const labelIOList& boundaryProcAddressing = procAddressing
1133 (
1134 procMeshList,
1135 proci,
1136 "boundaryProcAddressing",
1137 boundaryProcAddressingList
1138 );
1139
1140
1141 // FV fields: volume, surface, internal
1142 {
1143 if (!fieldDecomposerList.set(proci))
1144 {
1145 fieldDecomposerList.set
1146 (
1147 proci,
1149 (
1150 mesh,
1151 procMesh,
1152 faceProcAddressing,
1153 cellProcAddressing,
1154 boundaryProcAddressing
1155 )
1156 );
1157 }
1158
1159 volumeFieldCache.decomposeAllFields
1160 (
1161 fieldDecomposerList[proci]
1162 );
1163
1164 if (times.size() == 1)
1165 {
1166 // Clear cached decomposer
1167 fieldDecomposerList.set(proci, nullptr);
1168 }
1169 }
1170
1171
1172 // Point fields
1173 if (!pointFieldCache.empty())
1174 {
1175 const labelIOList& pointProcAddressing = procAddressing
1176 (
1177 procMeshList,
1178 proci,
1179 "pointProcAddressing",
1180 pointProcAddressingList
1181 );
1182
1183 const pointMesh& procPMesh =
1185
1186 if (!pointBoundaryProcAddressingList.set(proci))
1187 {
1188 pointBoundaryProcAddressingList.set
1189 (
1190 proci,
1192 (
1193 IOobject
1194 (
1195 "boundaryProcAddressing",
1196 procMesh.facesInstance(),
1199 procPMesh.thisDb(),
1203 ),
1204 boundaryProcAddressing
1205 )
1206 );
1207 }
1208 const auto& pointBoundaryProcAddressing =
1209 pointBoundaryProcAddressingList[proci];
1210
1211
1212 if (!pointFieldDecomposerList.set(proci))
1213 {
1214 pointFieldDecomposerList.set
1215 (
1216 proci,
1218 (
1219 pMesh,
1220 procPMesh,
1221 pointProcAddressing,
1222 pointBoundaryProcAddressing
1223 )
1224 );
1225 }
1226
1227 pointFieldCache.decomposeAllFields
1228 (
1229 pointFieldDecomposerList[proci]
1230 );
1231
1232 if (times.size() == 1)
1233 {
1234 // Early deletion
1235 pointBoundaryProcAddressingList.set
1236 (
1237 proci,
1238 nullptr
1239 );
1240 pointProcAddressingList.set(proci, nullptr);
1241 pointFieldDecomposerList.set(proci, nullptr);
1242 }
1243 }
1244
1245
1246 // If there is lagrangian data write it out
1247 forAll(lagrangianPositions, cloudi)
1248 {
1249 if (lagrangianPositions[cloudi].size())
1250 {
1251 lagrangianFieldDecomposer fieldDecomposer
1252 (
1253 mesh,
1254 procMesh,
1255 faceProcAddressing,
1256 cellProcAddressing,
1257 cloudDirs[cloudi],
1258 lagrangianPositions[cloudi],
1259 cellParticles[cloudi]
1260 );
1261
1262 // Lagrangian fields
1263 lagrangianFieldCache.decomposeAllFields
1264 (
1265 cloudi,
1266 cloudDirs[cloudi],
1267 fieldDecomposer
1268 );
1269 }
1270 }
1271
1272 if (doDecompFields)
1273 {
1274 // Decompose "uniform" directory in the time region
1275 // directory
1276 decomposeUniform
1277 (
1278 copyUniform, mesh, processorDb, regionDir
1279 );
1280
1281 // For a multi-region case, also decompose "uniform"
1282 // directory in the time directory
1283 if (regionNames.size() > 1 && regioni == 0)
1284 {
1285 decomposeUniform(copyUniform, mesh, processorDb);
1286 }
1287 }
1288
1289
1290 // We have cached all the constant mesh data for the current
1291 // processor. This is only important if running with
1292 // multiple times, otherwise it is just extra storage.
1293 if (times.size() == 1)
1294 {
1295 boundaryProcAddressingList.set(proci, nullptr);
1296 cellProcAddressingList.set(proci, nullptr);
1297 faceProcAddressingList.set(proci, nullptr);
1298 procMeshList.set(proci, nullptr);
1299 processorDbList.set(proci, nullptr);
1300 }
1301 }
1302
1303
1304 // Finite-area mesh and field decomposition
1305 for (auto& iter : faMeshDecompHashes.sorted())
1306 {
1307 const word& areaName = iter.key();
1308
1309 faMeshDecomposition& aMesh = *(iter.val());
1310
1311 Info<< "\nFinite area mesh decomposition: "
1312 << areaName << endl;
1313
1314 aMesh.decomposeMesh();
1315 aMesh.writeDecomposition();
1316
1317
1318 // Area/edge fields
1319 // ~~~~~~~~~~~~~~~~
1320
1321 faFieldDecomposer::fieldsCache areaFieldCache;
1322
1323 if
1324 (
1325 const auto objs = faObjects.cfind(areaName);
1326 doDecompFields && objs.good()
1327 )
1328 {
1329 areaFieldCache.readAllFields(aMesh, objs.val());
1330 }
1331
1332 const label nAreaFields = areaFieldCache.size();
1333
1334 Info<< endl;
1335 Info<< "Finite area field transfer: "
1336 << nAreaFields << " fields" << endl;
1337
1338 // Split the fields over processors
1339 for
1340 (
1341 label proci = 0;
1342 nAreaFields && proci < mesh.nProcs();
1343 ++proci
1344 )
1345 {
1346 Info<< " Processor " << proci << endl;
1347
1348 // open the database
1349 Time processorDb
1350 (
1352 args.rootPath(),
1353 args.caseName()/("processor" + Foam::name(proci)),
1354 false, // No function objects
1355 false // No extra controlDict libs
1356 );
1357
1358 processorDb.setTime(runTime);
1359
1360 // Read the volume mesh
1361 fvMesh procFvMesh
1362 (
1363 IOobject
1364 (
1365 regionName,
1366 processorDb.timeName(),
1367 processorDb
1368 )
1369 );
1370
1371 faMesh procMesh(areaName, procFvMesh);
1372
1373 // // Does not work. HJ, 15/Aug/2017
1374 // const labelIOList& faceProcAddressing =
1375 // procAddressing
1376 // (
1377 // procMeshList,
1378 // proci,
1379 // "faceProcAddressing",
1380 // faceProcAddressingList
1381 // );
1382
1383 // const labelIOList& boundaryProcAddressing =
1384 // procAddressing
1385 // (
1386 // procMeshList,
1387 // proci,
1388 // "boundaryProcAddressing",
1389 // boundaryProcAddressingList
1390 // );
1391
1392 // Addressing from faMesh (not polyMesh) meshSubDir
1393
1394 autoPtr<labelIOList> tfaceProcAddr =
1395 faProcAddressing
1396 (
1397 procMesh,
1398 "faceProcAddressing",
1399 runTime.constant()
1400 );
1401 auto& faceProcAddressing = *tfaceProcAddr;
1402
1403 autoPtr<labelIOList> tboundaryProcAddr =
1404 faProcAddressing
1405 (
1406 procMesh,
1407 "boundaryProcAddressing",
1408 runTime.constant()
1409 );
1410 auto& boundaryProcAddressing = *tboundaryProcAddr;
1411
1412 autoPtr<labelIOList> tedgeProcAddr =
1413 faProcAddressing
1414 (
1415 procMesh,
1416 "edgeProcAddressing",
1417 runTime.constant()
1418 );
1419 const auto& edgeProcAddressing = *tedgeProcAddr;
1420
1421 faFieldDecomposer fieldDecomposer
1422 (
1423 aMesh,
1424 procMesh,
1425 edgeProcAddressing,
1426 faceProcAddressing,
1427 boundaryProcAddressing
1428 );
1429
1430 areaFieldCache.decomposeAllFields
1431 (
1432 fieldDecomposer,
1433 args.verbose() // report
1434 );
1435 }
1436 }
1437 }
1438 }
1439 }
1440
1441 Info<< "\nEnd\n" << endl;
1442
1443 return 0;
1444}
1445
1446
1447// ************************************************************************* //
Functions used by OpenFOAM that are specific to POSIX compliant operating systems and need to be repl...
Required Classes.
Required Classes.
Base cloud calls templated on particle type.
Definition Cloud.H:64
A HashTable of pointers to objects of type <T>, with deallocation management of the pointers.
bool set(const Key &key, T *ptr)
Assign a new entry, overwrites existing.
bool good() const noexcept
True if iterator points to an entry.
Definition HashTable.H:990
A HashTable similar to std::unordered_map.
Definition HashTable.H:124
const_iterator cfind(const Key &key) const
Find and return an const_iterator set at the hashed entry.
Definition HashTableI.H:113
void reserve(label numEntries)
Reserve space for at least the specified number of elements (not the number of buckets) and regenerat...
Definition HashTable.C:729
bool emplace_set(const Key &key, Args &&... args)
Emplace set an entry, overwriting any existing entries.
Definition HashTableI.H:141
UPtrList< node_type > sorted()
Non-const access to the hash-table contents in sorted order (sorted by keys).
Definition HashTable.C:200
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
List of IOobjects with searching and retrieving facilities. Implemented as a HashTable,...
autoPtr< IOobject > remove(const IOobject &io)
Remove object from the list by its IOobject::name().
@ NO_REGISTER
Do not request registration (bool: false).
@ NO_READ
Nothing to be read.
@ READ_IF_PRESENT
Reading is optional [identical to LAZY_READ].
@ MUST_READ
Reading required.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
static IOobject selectIO(const IOobject &io, const fileName &altFile, const word &ioName="")
Return the IOobject, but also consider an alternative file name.
Definition IOobject.C:256
static FOAM_NO_DANGLING_REFERENCE const pointMesh & New(const polyMesh &mesh, Args &&... args)
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition PtrList.H:67
const T * set(const label i) const
Return const pointer to element (can be nullptr), or nullptr for out-of-range access (ie,...
Definition PtrList.H:171
T & try_emplace(const label i, Args &&... args)
Like emplace_set() but will not overwrite an occupied (non-null) location.
Definition PtrListI.H:210
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition Time.H:75
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
fileName timePath() const
Return current time path = path/timeName.
Definition Time.H:512
virtual void setTime(const Time &t)
Reset the time and time-index to those of the given time.
Definition Time.C:905
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
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 noParallel()
Remove the parallel options.
Definition argList.C:599
static void addOption(const word &optName, const string &param="", const string &usage="", bool advanced=false)
Add an option to validOptions with usage information.
Definition argList.C:400
static void addNote(const string &note)
Add extra notes for the usage information.
Definition argList.C:477
static void addOptionCompat(const word &optName, std::pair< const char *, int > compat)
Specify an alias for the option name.
Definition argList.C:433
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
static autoPtr< T > New(Args &&... args)
Construct autoPtr with forwarding arguments.
Definition autoPtr.H:178
static const word prefix
The prefix to local: lagrangian.
Definition cloud.H:79
static label nDomains(const dictionary &decompDict, const word &regionName="")
Return region-specific or top-level numberOfSubdomains entry.
static const word canonicalName
The canonical name ("decomposeParDict") under which the MeshObject is registered.
Testing of domain decomposition for finite-volume meshes.
Automatic domain decomposition class for finite-volume meshes.
Finite area boundary mesh, which is a faPatch list with registered IO, a reference to the associated ...
label size() const
Number of fields.
void readAllFields(const faMesh &mesh, const IOobjectList &objects)
Read all fields given mesh and objects.
void decomposeAllFields(const faFieldDecomposer &decomposer, bool report=false) const
Decompose and write all fields.
Finite Area area and edge field decomposer.
Automatic faMesh decomposition class.
bool writeDecomposition() const
Write decomposition.
void decomposeMesh()
Decompose mesh.
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 const objectRegistry & Registry(const polyMesh &pMesh)
Return the singleton parent registry (on the polyMesh) that contains all objects related to finite-ar...
Definition faMesh.C:157
static word meshSubDir
The mesh sub-directory name (usually "faMesh").
Definition faMesh.H:750
A class for handling file names.
Definition fileName.H:75
@ DIRECTORY
A directory.
Definition fileName.H:85
static std::string path(const std::string &str)
Return directory path name (part before last /).
Definition fileNameI.H:169
static bool isAbsolute(const std::string &str)
Return true if filename starts with a '/' or '\' or (windows-only) with a filesystem-root.
Definition fileNameI.H:129
An encapsulation of filesystem-related operations.
void readAllFields(const fvMesh &mesh, const IOobjectList &objects)
Read all fields given mesh and objects.
void decomposeAllFields(const fvFieldDecomposer &decomposer, bool report=false) const
Decompose and write all fields.
Finite Volume volume and surface field decomposer.
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
const Time & time() const
Return the top-level database.
Definition fvMesh.H:360
Adds label index to base particle.
Registry of regIOobjects.
void readAllFields(const pointMesh &mesh, const IOobjectList &objects)
Read all fields given mesh and objects.
bool empty() const noexcept
No fields.
void decomposeAllFields(const pointFieldDecomposer &decomposer, bool report=false) const
Decompose and write all fields.
Point field decomposer.
Mesh representing a set of points created from polyMesh.
Definition pointMesh.H:49
const objectRegistry & thisDb() const
Return database. For now is its polyMesh.
Definition pointMesh.H:201
static word meshSubDir
Return the mesh sub-directory name (usually "pointMesh").
Definition pointMesh.H:107
static const word & regionName(const word &region)
The mesh region name or word::null if polyMesh::defaultRegion.
Definition polyMesh.C:796
const fileName & facesInstance() const
Return the current instance directory for faces.
Definition polyMesh.C:844
static word defaultRegion
Return the default region name.
Definition polyMesh.H:406
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh").
Definition polyMesh.H:411
A class for handling character strings derived from std::string.
Definition string.H:76
bool starts_with(char c) const
True if string starts with given character (cf. C++20).
Definition string.H:436
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...
A class for handling words, derived from Foam::string.
Definition word.H:66
static const word null
An empty word.
Definition word.H:84
volScalarField & p
bool local
Definition EEqn.H:20
dynamicFvMesh & mesh
engineTime & runTime
Foam::word regionName(args.getOrDefault< word >("region", Foam::polyMesh::defaultRegion))
#define FatalErrorIn(functionName)
Report an error message using Foam::FatalError.
Definition error.H:592
const auto & io
const word & regionDir
wordList areaRegionNames
wordList regionNames
bool allowFaModels() noexcept
The enable/disable state for regionFaModel (default: true).
Namespace for OpenFOAM.
fileName cwd()
The physical or logical current working directory path name.
Definition POSIX.C:592
LList< SLListBase, T > SLList
Definition SLListFwd.H:41
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition int32.H:127
bool exists(const fileName &name, const bool checkGzip=true, const bool followLink=true)
Does the name exist (as DIRECTORY or FILE) in the file system?
Definition POSIX.C:837
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
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().
messageStream Info
Information stream (stdout output on master, null elsewhere).
List< instant > instantList
List of instants.
Definition instantList.H:41
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
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
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
bool isDir(const fileName &name, const bool followLink=true)
Does the name exist as a DIRECTORY in the file system?
Definition POSIX.C:862
Ostream & flush(Ostream &os)
Flush stream.
Definition Ostream.H:509
fileNameList readDir(const fileName &directory, const fileName::Type type=fileName::Type::FILE, const bool filtergz=true, const bool followLink=true)
Read a directory and return the entries as a fileName List.
Definition POSIX.C:965
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
bool chDir(const fileName &dir)
Change current directory to the one specified and return true on success.
Definition POSIX.C:609
Foam::argList args(argc, argv)
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
#define FOAM_NO_DANGLING_REFERENCE
Definition stdFoam.H:80
#define forAllReverse(list, i)
Reverse loop across all elements in list.
Definition stdFoam.H:315