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 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 fvMesh& procMesh,
177 const word& name,
178 const word& instance,
179 const word& local = fvMesh::meshSubDir
180)
181{
183 (
185 (
186 name,
187 instance,
188 local,
189 procMesh,
192 false // do not register
193 )
194 );
195}
196
197
198// Read proc addressing at specific instance.
199// Uses the finiteArea meshSubDir
200autoPtr<labelIOList> faProcAddressing
201(
202 const fvMesh& procMesh,
203 const word& name,
204 const word& instance,
205 const word& local = faMesh::meshSubDir
206)
207{
208 return procAddressing(procMesh, name, instance, local);
209}
210
211
212// Return cached or read proc addressing from facesInstance
213const labelIOList& procAddressing
214(
215 const PtrList<fvMesh>& procMeshList,
216 const label proci,
217 const word& name,
218 PtrList<labelIOList>& procAddressingList
219)
220{
221 const fvMesh& procMesh = procMeshList[proci];
222
223 if (!procAddressingList.set(proci))
224 {
225 procAddressingList.set
226 (
227 proci,
228 procAddressing(procMesh, name, procMesh.facesInstance())
229 );
230 }
231 return procAddressingList[proci];
232}
233
234
235void decomposeUniform
236(
237 const bool copyUniform,
239 const Time& processorDb,
240 const word& regionDir = word::null
241)
242{
243 const Time& runTime = mesh.time();
244
245 // Any uniform data to copy/link?
246 const fileName uniformDir(regionDir/"uniform");
247
248 if (fileHandler().isDir(runTime.timePath()/uniformDir))
249 {
250 Info<< "Detected additional non-decomposed files in "
251 << runTime.timePath()/uniformDir
252 << endl;
253
254 const fileName timePath =
255 fileHandler().filePath(processorDb.timePath());
256
257 // If no fields have been decomposed the destination
258 // directory will not have been created so make sure.
259 mkDir(timePath);
260
261 if (copyUniform || mesh.distributed())
262 {
263 if (!fileHandler().exists(timePath/uniformDir))
264 {
266 (
267 runTime.timePath()/uniformDir,
268 timePath/uniformDir
269 );
270 }
271 }
272 else
273 {
274 // Link with relative paths
275 string parentPath = string("..")/"..";
276
277 if (!regionDir.empty())
278 {
279 parentPath = parentPath/"..";
280 }
281
282 fileName currentDir(cwd());
283 chDir(timePath);
284
285 if (!fileHandler().exists(uniformDir))
286 {
288 (
289 parentPath/runTime.timeName()/uniformDir,
290 uniformDir
291 );
292 }
293 chDir(currentDir);
294 }
295 }
296}
297
298} // End namespace Foam
299
300
301using namespace Foam;
302
303// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
304
305int main(int argc, char *argv[])
306{
307 argList::addNote
308 (
309 "Decompose a mesh and fields of a case for parallel execution"
310 );
311
312 argList::noParallel();
313 argList::addOption
314 (
315 "decomposeParDict",
316 "file",
317 "Use specified file for decomposePar dictionary"
318 );
319
320 #include "addAllRegionOptions.H"
321
322 argList::addDryRunOption
323 (
324 "Test without writing the decomposition. "
325 "Changes -cellDist to only write VTK output."
326 );
327 argList::addVerboseOption("Additional verbosity");
328 argList::addOption
329 (
330 "domains",
331 "N",
332 "Override numberOfSubdomains (-dry-run only)",
333 true // Advanced option
334 );
335 argList::addOption
336 (
337 "method",
338 "name",
339 "Override decomposition method (-dry-run only)",
340 true // Advanced option
341 );
342
343 argList::addBoolOption
344 (
345 "no-finite-area",
346 "Suppress finiteArea mesh/field decomposition",
347 true // Advanced option
348 );
349
350 argList::addBoolOption
351 (
352 "no-lagrangian",
353 "Suppress lagrangian (cloud) decomposition",
354 true // Advanced option
355 );
356
357 argList::addBoolOption
358 (
359 "cellDist",
360 "Write cell distribution as a labelList - for use with 'manual' "
361 "decomposition method and as a volScalarField for visualization."
362 );
363 argList::addBoolOption
364 (
365 "copyZero",
366 "Copy 0/ directory to processor*/ rather than decompose the fields"
367 );
368 argList::addBoolOption
369 (
370 "copyUniform",
371 "Copy any uniform/ directories too"
372 );
373 argList::addBoolOption
374 (
375 "fields",
376 "Use existing geometry decomposition and convert fields only"
377 );
378 argList::addBoolOption
379 (
380 "no-fields",
381 "Suppress conversion of fields (volume, finite-area, lagrangian)"
382 );
383
384 argList::addBoolOption
385 (
386 "no-sets",
387 "Skip decomposing cellSets, faceSets, pointSets"
388 );
389 argList::addOptionCompat("no-sets", {"noSets", 2106});
390
391 argList::addBoolOption
392 (
393 "force",
394 "Remove existing processor*/ subdirs before decomposing the geometry"
395 );
396 argList::addBoolOption
397 (
398 "ifRequired",
399 "Only decompose geometry if the number of domains has changed"
400 );
401
402 // Allow explicit -constant, have zero from time range
403 timeSelector::addOptions(true, false); // constant(true), zero(false)
404
405 #include "setRootCase.H"
406
407 const bool writeCellDist = args.found("cellDist");
408
409 // Most of these are ignored for dry-run (not triggered anywhere)
410 const bool copyZero = args.found("copyZero");
411 const bool copyUniform = args.found("copyUniform");
412 const bool decomposeSets = !args.found("no-sets");
413
414 const bool decomposeIfRequired = args.found("ifRequired");
415
416 const bool doDecompFields = !args.found("no-fields");
417 const bool doFiniteArea = !args.found("no-finite-area");
418 const bool doLagrangian = !args.found("no-lagrangian");
419
420 bool decomposeFieldsOnly = args.found("fields");
421 bool forceOverwrite = args.found("force");
422
423
424 // Set time from database
425 #include "createTime.H"
426
427 // Allow override of time (unless dry-run)
428 instantList times;
429 if (args.dryRun())
430 {
431 Info<< "\ndry-run: ignoring -copy*, -fields, -force, time selection"
432 << nl;
433 }
434 else
435 {
436 if (decomposeFieldsOnly && !doDecompFields)
437 {
439 << "Options -fields and -no-fields are mutually exclusive"
440 << " ... giving up" << nl
441 << exit(FatalError);
442 }
443
444 if (!doDecompFields)
445 {
446 Info<< "Skip decompose of all fields" << nl;
447 }
448 if (!doFiniteArea)
449 {
450 Info<< "Skip decompose of finiteArea mesh/fields" << nl;
451 }
452 if (!doLagrangian)
453 {
454 Info<< "Skip decompose of lagrangian positions/fields" << nl;
455 }
456
457 times = timeSelector::selectIfPresent(runTime, args);
458 }
459
460
461 // Allow override of decomposeParDict location
462 fileName decompDictFile(args.get<fileName>("decomposeParDict", ""));
463 if (!decompDictFile.empty() && !decompDictFile.isAbsolute())
464 {
465 decompDictFile = runTime.globalPath()/decompDictFile;
466 }
467
468 // Get region names
469 #include "getAllRegionOptions.H"
470
471 const bool optRegions =
472 (regionNames.size() != 1 || regionNames[0] != polyMesh::defaultRegion);
473
474 if (regionNames.size() == 1 && regionNames[0] != polyMesh::defaultRegion)
475 {
476 Info<< "Using region: " << regionNames[0] << nl << endl;
477 }
478
479 forAll(regionNames, regioni)
480 {
481 const word& regionName = regionNames[regioni];
482 const word& regionDir = polyMesh::regionName(regionName);
483
484 if (args.dryRun())
485 {
486 Info<< "dry-run: decomposing mesh " << regionName << nl << nl
487 << "Create mesh..." << flush;
488
490 (
492 (
494 runTime.timeName(),
495 runTime,
496 IOobject::MUST_READ,
497 IOobject::NO_WRITE,
498 false
499 ),
500 decompDictFile,
501 args.getOrDefault<label>("domains", 0),
502 args.getOrDefault<word>("method", word::null)
503 );
504
505 decompTest.execute(writeCellDist, args.verbose());
506 continue;
507 }
508
509 Info<< "\n\nDecomposing mesh";
510 if (!regionDir.empty())
511 {
512 Info<< ' ' << regionName;
513 }
514 Info<< nl << endl;
515
516 // Determine the existing processor count directly
517 const label nProcsOld =
519
520 // Get requested numberOfSubdomains directly from the dictionary.
521 // Note: have no mesh yet so cannot use decompositionModel::New
522 const label nDomains = decompositionMethod::nDomains
523 (
525 (
526 IOobject::selectIO
527 (
529 (
530 decompositionModel::canonicalName,
531 runTime.time().system(),
532 regionDir, // region (if non-default)
533 runTime,
534 IOobject::MUST_READ,
535 IOobject::NO_WRITE,
536 false // do not register
537 ),
538 decompDictFile
539 )
540 )
541 );
542
543 // Give file handler a chance to determine the output directory
544 const_cast<fileOperation&>(fileHandler()).setNProcs(nDomains);
545
546 if (decomposeFieldsOnly)
547 {
548 // Sanity check on previously decomposed case
549 if (nProcsOld != nDomains)
550 {
552 << "Specified -fields, but the case was decomposed with "
553 << nProcsOld << " domains"
554 << nl
555 << "instead of " << nDomains
556 << " domains as specified in decomposeParDict" << nl
557 << exit(FatalError);
558 }
559 }
560 else if (nProcsOld)
561 {
562 bool procDirsProblem = true;
563
564 if (decomposeIfRequired && nProcsOld == nDomains)
565 {
566 // We can reuse the decomposition
567 decomposeFieldsOnly = true;
568 procDirsProblem = false;
569 forceOverwrite = false;
570
571 Info<< "Using existing processor directories" << nl;
572 }
573
574 if (optRegions)
575 {
576 procDirsProblem = false;
577 forceOverwrite = false;
578 }
579
580 if (forceOverwrite)
581 {
582 Info<< "Removing " << nProcsOld
583 << " existing processor directories" << endl;
584
585 // Remove existing processors directory
586 fileNameList dirs
587 (
588 fileHandler().readDir
589 (
590 runTime.path(),
591 fileName::Type::DIRECTORY
592 )
593 );
594 forAllReverse(dirs, diri)
595 {
596 const fileName& d = dirs[diri];
597
598 label proci = -1;
599
600 if
601 (
602 d.starts_with("processor")
603 &&
604 (
605 // Collated is "processors"
606 d[9] == 's'
607
608 // Uncollated has integer(s) after 'processor'
609 || Foam::read(d.substr(9), proci)
610 )
611 )
612 {
613 if (fileHandler().exists(d))
614 {
615 fileHandler().rmDir(d);
616 }
617 }
618 }
619
620 procDirsProblem = false;
621 }
622
623 if (procDirsProblem)
624 {
626 << "Case is already decomposed with " << nProcsOld
627 << " domains, use the -force option or manually" << nl
628 << "remove processor directories before decomposing. e.g.,"
629 << nl
630 << " rm -rf " << runTime.path().c_str() << "/processor*"
631 << nl
632 << exit(FatalError);
633 }
634 }
635
636 Info<< "Create mesh" << endl;
638 (
640 (
642 runTime.timeName(),
643 runTime,
644 IOobject::NO_READ,
645 IOobject::NO_WRITE,
646 false
647 ),
648 decompDictFile
649 );
650
651 // Decompose the mesh
652 if (!decomposeFieldsOnly)
653 {
654 mesh.decomposeMesh();
655 mesh.writeDecomposition(decomposeSets);
656
657 if (writeCellDist)
658 {
659 const labelList& procIds = mesh.cellToProc();
660
661 // Write decomposition for visualization
662 mesh.writeVolField("cellDist");
663 //TBD: mesh.writeVTK("cellDist");
664
665 // Write decomposition as labelList for use with 'manual'
666 // decomposition method.
667 labelIOList cellDecomposition
668 (
670 (
671 "cellDecomposition",
672 mesh.facesInstance(),
673 mesh,
674 IOobject::NO_READ,
675 IOobject::NO_WRITE,
676 false
677 ),
678 procIds
679 );
680 cellDecomposition.write();
681
682 Info<< nl << "Wrote decomposition to "
683 << cellDecomposition.objectRelPath()
684 << " for use in manual decomposition." << endl;
685 }
686
687 fileHandler().flush();
688 }
689
690
691 if (copyZero)
692 {
693 // Copy the 0/ directory into each of the processor directories
694 // with fallback of 0.orig/ directory if necessary.
695
696 fileName inputDir(runTime.path()/"0");
697
698 bool canCopy = fileHandler().isDir(inputDir);
699 if (!canCopy)
700 {
701 // Try with "0.orig" instead
702 inputDir.ext("orig");
703 canCopy = fileHandler().isDir(inputDir);
704 }
705
706 if (canCopy)
707 {
708 // Avoid copying into the same directory multiple times
709 // (collated format). Don't need a hash here.
710 fileName prevOutputDir;
711 for (label proci = 0; proci < mesh.nProcs(); ++proci)
712 {
713 Time processorDb
714 (
715 Time::controlDictName,
716 args.rootPath(),
717 args.caseName()/("processor" + Foam::name(proci))
718 );
719 // processorDb.setTime(runTime);
720
721 // Get corresponding directory name
722 // (to handle processors/)
723 const fileName outputDir
724 (
725 fileHandler().objectPath
726 (
728 (
729 word::null, // name
730 "0", // instance (time == 0)
731 processorDb
732 ),
734 )
735 );
736
737 if (outputDir != prevOutputDir)
738 {
739 Info<< "Processor " << proci
740 << ": copying \""
741 << inputDir.name() << "/\" to "
742 << runTime.relativePath(outputDir)
743 << endl;
744
745 fileHandler().cp(inputDir, outputDir);
746 prevOutputDir = outputDir;
747 }
748 }
749 }
750 else
751 {
752 Info<< "No 0/ or 0.orig/ directory to copy" << nl;
753 }
754 }
755 else
756 {
757 // Decompose field files, lagrangian, finite-area
758
759 // Cached processor meshes and maps. These are only preserved if
760 // running with multiple times.
761 PtrList<Time> processorDbList(mesh.nProcs());
762 PtrList<fvMesh> procMeshList(mesh.nProcs());
763 PtrList<labelIOList> faceProcAddressingList(mesh.nProcs());
764 PtrList<labelIOList> cellProcAddressingList(mesh.nProcs());
765 PtrList<labelIOList> boundaryProcAddressingList(mesh.nProcs());
766 PtrList<labelIOList> pointProcAddressingList(mesh.nProcs());
767
768 PtrList<fvFieldDecomposer> fieldDecomposerList(mesh.nProcs());
769 PtrList<pointFieldDecomposer> pointFieldDecomposerList
770 (
771 mesh.nProcs()
772 );
773
774
775 // Loop over all times
776 forAll(times, timei)
777 {
778 runTime.setTime(times[timei], timei);
779
780 Info<< "Time = " << runTime.timeName() << endl;
781
782 // Field objects at this time
783 IOobjectList objects;
784
785 if (doDecompFields)
786 {
787 objects = IOobjectList(mesh, runTime.timeName());
788
789 // Ignore generated fields: (cellDist)
790 objects.remove("cellDist");
791 }
792
793 // Finite area handling
794 autoPtr<faMeshDecomposition> faMeshDecompPtr;
795 if (doFiniteArea)
796 {
798 (
799 "faBoundary",
800 mesh.time().findInstance
801 (
802 mesh.dbDir()/polyMesh::meshSubDir,
803 "boundary"
804 ),
805 faMesh::meshSubDir,
806 mesh,
807 IOobject::READ_IF_PRESENT,
808 IOobject::NO_WRITE,
809 false // not registered
810 );
811
812 if (io.typeHeaderOk<faBoundaryMesh>(true))
813 {
814 // Always based on the volume decomposition!
815 faMeshDecompPtr.reset
816 (
818 (
819 mesh,
820 mesh.nProcs(),
821 mesh.model()
822 )
823 );
824 }
825 }
826
827
828 // Volume/surface/internal fields
829 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
830
831 fvFieldDecomposer::fieldsCache volumeFieldCache;
832
833 if (doDecompFields)
834 {
835 volumeFieldCache.readAllFields(mesh, objects);
836 }
837
838
839 // Point fields
840 // ~~~~~~~~~~~~
841 const pointMesh& pMesh = pointMesh::New(mesh);
842
843 pointFieldDecomposer::fieldsCache pointFieldCache;
844
845 if (doDecompFields)
846 {
847 pointFieldCache.readAllFields(pMesh, objects);
848 }
849
850
851 // Lagrangian fields
852 // ~~~~~~~~~~~~~~~~~
853
854 fileNameList cloudDirs;
855
856 if (doDecompFields && doLagrangian)
857 {
858 cloudDirs = fileHandler().readDir
859 (
860 runTime.timePath()/cloud::prefix,
861 fileName::DIRECTORY
862 );
863 }
864
865 // Particles
866 PtrList<Cloud<indexedParticle>> lagrangianPositions
867 (
868 cloudDirs.size()
869 );
870 // Particles per cell
872 (
873 cloudDirs.size()
874 );
875
876 lagrangianFieldDecomposer::fieldsCache lagrangianFieldCache
877 (
878 cloudDirs.size()
879 );
880
881 label cloudI = 0;
882
883 for (const fileName& cloudDir : cloudDirs)
884 {
885 IOobjectList cloudObjects
886 (
887 mesh,
888 runTime.timeName(),
889 cloud::prefix/cloudDir,
890 IOobject::MUST_READ,
891 IOobject::NO_WRITE,
892 false
893 );
894
895 // Note: look up "positions" for backwards compatibility
896 if
897 (
898 cloudObjects.found("coordinates")
899 || cloudObjects.found("positions")
900 )
901 {
902 // Read lagrangian particles
903 // ~~~~~~~~~~~~~~~~~~~~~~~~~
904
905 Info<< "Identified lagrangian data set: "
906 << cloudDir << endl;
907
908 lagrangianPositions.set
909 (
910 cloudI,
912 (
913 mesh,
914 cloudDir,
915 false
916 )
917 );
918
919
920 // Sort particles per cell
921 // ~~~~~~~~~~~~~~~~~~~~~~~
922
923 cellParticles.set
924 (
925 cloudI,
927 (
928 mesh.nCells(),
929 static_cast<SLList<indexedParticle*>*>(nullptr)
930 )
931 );
932
933 label i = 0;
934
935 for (indexedParticle& p : lagrangianPositions[cloudI])
936 {
937 p.index() = i++;
938
939 label celli = p.cell();
940
941 // Check
942 if (celli < 0 || celli >= mesh.nCells())
943 {
945 << "Illegal cell number " << celli
946 << " for particle with index "
947 << p.index()
948 << " at position "
949 << p.position() << nl
950 << "Cell number should be between 0 and "
951 << mesh.nCells()-1 << nl
952 << "On this mesh the particle should"
953 << " be in cell "
954 << mesh.findCell(p.position())
955 << exit(FatalError);
956 }
957
958 if (!cellParticles[cloudI][celli])
959 {
960 cellParticles[cloudI][celli] =
962 }
963
964 cellParticles[cloudI][celli]->append(&p);
965 }
966
967 // Read fields
968 // ~~~~~~~~~~~
969
970 IOobjectList lagrangianObjects
971 (
972 mesh,
973 runTime.timeName(),
974 cloud::prefix/cloudDirs[cloudI],
975 IOobject::MUST_READ,
976 IOobject::NO_WRITE,
977 false
978 );
979
980 lagrangianFieldCache.readAllFields
981 (
982 cloudI,
983 lagrangianObjects
984 );
985
986 ++cloudI;
987 }
988 }
989
990 lagrangianPositions.resize(cloudI);
991 cellParticles.resize(cloudI);
992 lagrangianFieldCache.resize(cloudI);
993
994 Info<< endl;
995
996 // split the fields over processors
997 for
998 (
999 label proci = 0;
1000 doDecompFields && proci < mesh.nProcs();
1001 ++proci
1002 )
1003 {
1004 Info<< "Processor " << proci << ": field transfer" << endl;
1005
1006 // open the database
1007 if (!processorDbList.set(proci))
1008 {
1009 processorDbList.set
1010 (
1011 proci,
1012 new Time
1013 (
1014 Time::controlDictName,
1015 args.rootPath(),
1016 args.caseName()
1017 / ("processor" + Foam::name(proci))
1018 )
1019 );
1020 }
1021 Time& processorDb = processorDbList[proci];
1022
1023
1024 processorDb.setTime(runTime);
1025
1026 // read the mesh
1027 if (!procMeshList.set(proci))
1028 {
1029 procMeshList.set
1030 (
1031 proci,
1032 new fvMesh
1033 (
1034 IOobject
1035 (
1036 regionName,
1037 processorDb.timeName(),
1038 processorDb
1039 )
1040 )
1041 );
1042 }
1043 const fvMesh& procMesh = procMeshList[proci];
1044
1045 const labelIOList& faceProcAddressing = procAddressing
1046 (
1047 procMeshList,
1048 proci,
1049 "faceProcAddressing",
1050 faceProcAddressingList
1051 );
1052
1053 const labelIOList& cellProcAddressing = procAddressing
1054 (
1055 procMeshList,
1056 proci,
1057 "cellProcAddressing",
1058 cellProcAddressingList
1059 );
1060
1061 const labelIOList& boundaryProcAddressing = procAddressing
1062 (
1063 procMeshList,
1064 proci,
1065 "boundaryProcAddressing",
1066 boundaryProcAddressingList
1067 );
1068
1069
1070 // FV fields: volume, surface, internal
1071 {
1072 if (!fieldDecomposerList.set(proci))
1073 {
1074 fieldDecomposerList.set
1075 (
1076 proci,
1078 (
1079 mesh,
1080 procMesh,
1081 faceProcAddressing,
1082 cellProcAddressing,
1083 boundaryProcAddressing
1084 )
1085 );
1086 }
1087
1088 volumeFieldCache.decomposeAllFields
1089 (
1090 fieldDecomposerList[proci]
1091 );
1092
1093 if (times.size() == 1)
1094 {
1095 // Clear cached decomposer
1096 fieldDecomposerList.set(proci, nullptr);
1097 }
1098 }
1099
1100
1101 // Point fields
1102 if (!pointFieldCache.empty())
1103 {
1104 const labelIOList& pointProcAddressing = procAddressing
1105 (
1106 procMeshList,
1107 proci,
1108 "pointProcAddressing",
1109 pointProcAddressingList
1110 );
1111
1112 const pointMesh& procPMesh = pointMesh::New(procMesh);
1113
1114 if (!pointFieldDecomposerList.set(proci))
1115 {
1116 pointFieldDecomposerList.set
1117 (
1118 proci,
1120 (
1121 pMesh,
1122 procPMesh,
1123 pointProcAddressing,
1124 boundaryProcAddressing
1125 )
1126 );
1127 }
1128
1129 pointFieldCache.decomposeAllFields
1130 (
1131 pointFieldDecomposerList[proci]
1132 );
1133
1134 if (times.size() == 1)
1135 {
1136 pointProcAddressingList.set(proci, nullptr);
1137 pointFieldDecomposerList.set(proci, nullptr);
1138 }
1139 }
1140
1141
1142 // If there is lagrangian data write it out
1143 forAll(lagrangianPositions, cloudi)
1144 {
1145 if (lagrangianPositions[cloudi].size())
1146 {
1147 lagrangianFieldDecomposer fieldDecomposer
1148 (
1149 mesh,
1150 procMesh,
1151 faceProcAddressing,
1152 cellProcAddressing,
1153 cloudDirs[cloudi],
1154 lagrangianPositions[cloudi],
1155 cellParticles[cloudi]
1156 );
1157
1158 // Lagrangian fields
1159 lagrangianFieldCache.decomposeAllFields
1160 (
1161 cloudi,
1162 cloudDirs[cloudi],
1163 fieldDecomposer
1164 );
1165 }
1166 }
1167
1168 if (doDecompFields)
1169 {
1170 // Decompose "uniform" directory in the time region
1171 // directory
1172 decomposeUniform
1173 (
1174 copyUniform, mesh, processorDb, regionDir
1175 );
1176
1177 // For a multi-region case, also decompose "uniform"
1178 // directory in the time directory
1179 if (regionNames.size() > 1 && regioni == 0)
1180 {
1181 decomposeUniform(copyUniform, mesh, processorDb);
1182 }
1183 }
1184
1185
1186 // We have cached all the constant mesh data for the current
1187 // processor. This is only important if running with
1188 // multiple times, otherwise it is just extra storage.
1189 if (times.size() == 1)
1190 {
1191 boundaryProcAddressingList.set(proci, nullptr);
1192 cellProcAddressingList.set(proci, nullptr);
1193 faceProcAddressingList.set(proci, nullptr);
1194 procMeshList.set(proci, nullptr);
1195 processorDbList.set(proci, nullptr);
1196 }
1197 }
1198
1199
1200 // Finite area mesh and field decomposition
1201 if (faMeshDecompPtr)
1202 {
1203 Info<< "\nFinite area mesh decomposition" << endl;
1204
1205 faMeshDecomposition& aMesh = faMeshDecompPtr();
1206
1207 aMesh.decomposeMesh();
1208 aMesh.writeDecomposition();
1209
1210
1211 // Area/edge fields
1212 // ~~~~~~~~~~~~~~~~
1213
1214 faFieldDecomposer::fieldsCache areaFieldCache;
1215
1216 if (doDecompFields)
1217 {
1218 areaFieldCache.readAllFields(aMesh, objects);
1219 }
1220
1221 const label nAreaFields = areaFieldCache.size();
1222
1223 Info<< endl;
1224 Info<< "Finite area field transfer: "
1225 << nAreaFields << " fields" << endl;
1226
1227 // Split the fields over processors
1228 for
1229 (
1230 label proci = 0;
1231 nAreaFields && proci < mesh.nProcs();
1232 ++proci
1233 )
1234 {
1235 Info<< " Processor " << proci << endl;
1236
1237 // open the database
1238 Time processorDb
1239 (
1240 Time::controlDictName,
1241 args.rootPath(),
1242 args.caseName()
1243 / ("processor" + Foam::name(proci))
1244 );
1245
1246 processorDb.setTime(runTime);
1247
1248 // Read the mesh
1249 fvMesh procFvMesh
1250 (
1251 IOobject
1252 (
1253 regionName,
1254 processorDb.timeName(),
1255 processorDb
1256 )
1257 );
1258
1259 faMesh procMesh(procFvMesh);
1260
1261 // // Does not work. HJ, 15/Aug/2017
1262 // const labelIOList& faceProcAddressing =
1263 // procAddressing
1264 // (
1265 // procMeshList,
1266 // proci,
1267 // "faceProcAddressing",
1268 // faceProcAddressingList
1269 // );
1270
1271 // const labelIOList& boundaryProcAddressing =
1272 // procAddressing
1273 // (
1274 // procMeshList,
1275 // proci,
1276 // "boundaryProcAddressing",
1277 // boundaryProcAddressingList
1278 // );
1279
1280 // Addressing from faMesh (not polyMesh) meshSubDir
1281
1282 autoPtr<labelIOList> tfaceProcAddr =
1283 faProcAddressing
1284 (
1285 procFvMesh,
1286 "faceProcAddressing",
1287 runTime.constant()
1288 );
1289 auto& faceProcAddressing = *tfaceProcAddr;
1290
1291 autoPtr<labelIOList> tboundaryProcAddr =
1292 faProcAddressing
1293 (
1294 procFvMesh,
1295 "boundaryProcAddressing",
1296 runTime.constant()
1297 );
1298 auto& boundaryProcAddressing = *tboundaryProcAddr;
1299
1300 autoPtr<labelIOList> tedgeProcAddr =
1301 faProcAddressing
1302 (
1303 procFvMesh,
1304 "edgeProcAddressing",
1305 runTime.constant()
1306 );
1307 const auto& edgeProcAddressing = *tedgeProcAddr;
1308
1309 faFieldDecomposer fieldDecomposer
1310 (
1311 aMesh,
1312 procMesh,
1313 edgeProcAddressing,
1314 faceProcAddressing,
1315 boundaryProcAddressing
1316 );
1317
1318 areaFieldCache.decomposeAllFields(fieldDecomposer);
1319 }
1320 }
1321 }
1322 }
1323 }
1324
1325 Info<< "\nEnd\n" << endl;
1326
1327 return 0;
1328}
1329
1330
1331// ************************************************************************* //
Functions used by OpenFOAM that are specific to POSIX compliant operating systems and need to be repl...
Base cloud calls templated on particle type.
Definition: Cloud.H:68
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:57
List of IOobjects with searching and retrieving facilities.
Definition: IOobjectList.H:59
bool remove(const IOobject &io)
Remove IOobject from the list.
Definition: IOobjectList.C:259
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
Template class for non-intrusive linked lists.
Definition: LList.H:79
std::enable_if< std::is_same< bool, TypeT >::value, bool >::type set(const label i, bool val=true)
A bitSet::set() method for a list of bool.
Definition: List.H:330
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: PtrList.H:73
const T * set(const label i) const
Definition: PtrList.H:138
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:80
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:780
fileName timePath() const
Return current time path.
Definition: Time.H:375
virtual void setTime(const Time &t)
Reset the time and time-index to those of the given time.
Definition: Time.C:997
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
T get(const label index) const
Get a value from the argument at index.
Definition: argListI.H:278
int dryRun() const noexcept
Return the dry-run flag.
Definition: argListI.H:116
int verbose() const noexcept
Return the verbose flag.
Definition: argListI.H:128
const fileName & rootPath() const noexcept
Return root path.
Definition: argListI.H:63
bool found(const word &optName) const
Return true if the named option is found.
Definition: argListI.H:178
T getOrDefault(const word &optName, const T &deflt) const
Get a value from the named option if present, or return default.
Definition: argListI.H:307
const word & executable() const noexcept
Name of executable without the path.
Definition: argListI.H:51
const fileName & caseName() const noexcept
Return case name (parallel run) or global case (serial run)
Definition: argListI.H:69
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
void reset(autoPtr< T > &&other) noexcept
Delete managed object and set to new given pointer.
Definition: autoPtrI.H:117
Testing of domain decomposition for finite-volume meshes.
Automatic domain decomposition class for finite-volume meshes.
Finite area boundary mesh.
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.
Finite area mesh (used for 2-D non-Euclidian finite area method) defined using a patch of faces on a ...
Definition: faMesh.H:100
static word meshSubDir
The mesh sub-directory name (usually "faMesh")
Definition: faMesh.H:504
A class for handling file names.
Definition: fileName.H:76
An encapsulation of filesystem-related operations.
Definition: fileOperation.H:69
virtual void flush() const
Forcibly wait until all output done. Flush any cached data.
virtual bool cp(const fileName &src, const fileName &dst, const bool followLink=true) const =0
Copy, recursively if necessary, the source to the destination.
virtual bool rmDir(const fileName &dir, const bool silent=false) const =0
Remove a directory and its contents.
virtual bool ln(const fileName &src, const fileName &dst) const =0
Create a softlink. dst should not exist. Returns true if.
virtual label nProcs(const fileName &dir, const fileName &local="") const
Get number of processor directories/results. Used for e.g.
virtual bool isDir(const fileName &, const bool followLink=true) const =0
Does the name exist as a DIRECTORY in the file system?
virtual fileName filePath(const bool checkGlobal, const IOobject &, const word &typeName, const bool search=true) const =0
Search for an object. checkGlobal : also check undecomposed case.
virtual fileNameList readDir(const fileName &, const fileName::Type=fileName::FILE, const bool filtergz=true, const bool followLink=true) const =0
Read a directory and return the entries as a string list.
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:91
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:290
Adds label index to base particle.
void readAllFields(const pointMesh &mesh, const IOobjectList &objects)
Read all fields given mesh and objects.
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:55
const fileName & facesInstance() const
Return the current instance directory for faces.
Definition: polyMesh.C:866
A class for handling character strings derived from std::string.
Definition: string.H:79
bool starts_with(const std::string &s) const
True if string starts with the given prefix (cf. C++20)
Definition: string.H:297
A class for handling words, derived from Foam::string.
Definition: word.H:68
static const word null
An empty word.
Definition: word.H:80
volScalarField & p
dynamicFvMesh & mesh
engineTime & runTime
Foam::word regionName(Foam::polyMesh::defaultRegion)
#define FatalErrorIn(functionName)
Report an error message using Foam::FatalError.
Definition: error.H:448
const word & regionDir
wordList regionNames
IOobject io("surfaceFilmProperties", mesh.time().constant(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE, false)
Namespace for OpenFOAM.
const fileOperation & fileHandler()
Get current file handler.
fileName cwd()
The physical or logical current working directory path name.
Definition: MSwindows.C:476
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition: int32.H:108
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: MSwindows.C:633
bool mkDir(const fileName &pathName, mode_t mode=0777)
Make a directory and return an error if it could not be created.
Definition: MSwindows.C:515
messageStream Info
Information stream (stdout output on master, null elsewhere)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
IOList< label > labelIOList
Label container classes.
Definition: labelIOList.H:44
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
bool isDir(const fileName &name, const bool followLink=true)
Does the name exist as a DIRECTORY in the file system?
Definition: MSwindows.C:651
Ostream & flush(Ostream &os)
Flush stream.
Definition: Ostream.H:364
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
bool chDir(const fileName &dir)
Change current directory to the one specified and return true on success.
Definition: MSwindows.C:508
Foam::argList args(argc, argv)
faMesh aMesh(mesh)
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
#define forAllReverse(list, i)
Reverse loop across all elements in list.
Definition: stdFoam.H:346