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 -------------------------------------------------------------------------------
11 License
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 
27 Application
28  decomposePar
29 
30 Group
31  grpParallelUtilities
32 
33 Description
34  Automatically decomposes a mesh and fields of a case for parallel
35  execution of OpenFOAM.
36 
37 Usage
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 
169 namespace Foam
170 {
171 
172 // Read proc addressing at specific instance.
173 // Uses polyMesh/fvMesh meshSubDir by default
174 autoPtr<labelIOList> procAddressing
175 (
176  const objectRegistry& procRegistry,
177  const word& name,
178  const word& instance,
180 )
181 {
183  (
184  IOobject
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
200 autoPtr<labelIOList> faProcAddressing
201 (
202  const objectRegistry& procRegistry,
203  const word& name,
204  const word& instance,
205  const word& local = faMesh::meshSubDir
206 )
207 {
208  return procAddressing(procRegistry, name, instance, local);
209 }
210 
211 
212 // Return cached or read proc addressing from facesInstance
214 const 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,
227  IOobject
228  (
229  name,
230  procMesh.facesInstance(),
232  procMesh,
236  )
237  );
238 }
239 
240 
241 void decomposeUniform
242 (
243  const bool copyUniform,
244  const domainDecomposition& mesh,
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  (
264  IOobject
265  (
266  "dummy",
267  runTime.timeName(),
268  processorDb
269  ),
270  word::null
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 
317 using namespace Foam;
318 
319 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
320 
321 int 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 
339  (
340  "Test without writing the decomposition. "
341  "Changes -cellDist to only write VTK output."
342  );
345  (
346  "domains",
347  "N",
348  "Override numberOfSubdomains (-dry-run only)",
349  true // Advanced option
350  );
352  (
353  "method",
354  "name",
355  "Override decomposition method (-dry-run only)",
356  true // Advanced option
357  );
358 
360  (
361  "no-finite-area",
362  "Suppress finiteArea mesh/field decomposition",
363  true // Advanced option
364  );
365 
367  (
368  "no-lagrangian",
369  "Suppress lagrangian (cloud) decomposition",
370  true // Advanced option
371  );
372 
374  (
375  "cellDist",
376  "Write cell distribution as a labelList - for use with 'manual' "
377  "decomposition method and as a volScalarField for visualization."
378  );
380  (
381  "copyZero",
382  "Copy 0/ directory to processor*/ rather than decompose the fields"
383  );
385  (
386  "copyUniform",
387  "Copy any uniform/ directories too"
388  );
390  (
391  "fields",
392  "Use existing geometry decomposition and convert fields only"
393  );
395  (
396  "no-fields",
397  "Suppress conversion of fields (volume, finite-area, lagrangian)"
398  );
399 
401  (
402  "no-sets",
403  "Skip decomposing cellSets, faceSets, pointSets"
404  );
405  argList::addOptionCompat("no-sets", {"noSets", 2106});
406 
408  (
409  "force",
410  "Remove existing processor*/ subdirs before decomposing the geometry"
411  );
413  (
414  "ifRequired",
415  "Only decompose geometry if the number of domains has changed"
416  );
417 
418  // Allow explicit -constant, have zero from time range
419  timeSelector::addOptions(true, false); // constant(true), zero(false)
420 
421  #include "setRootCase.H"
422 
423  const bool writeCellDist = args.found("cellDist");
424 
425  // Most of these are ignored for dry-run (not triggered anywhere)
426  const bool copyZero = args.found("copyZero");
427  const bool copyUniform = args.found("copyUniform");
428  const bool decomposeSets = !args.found("no-sets");
429 
430  const bool decomposeIfRequired = args.found("ifRequired");
431 
432  const bool doDecompFields = !args.found("no-fields");
433  const bool doFiniteArea = !args.found("no-finite-area");
434  const bool doLagrangian = !args.found("no-lagrangian");
435 
436  bool decomposeFieldsOnly = args.found("fields");
437  bool forceOverwrite = args.found("force");
438 
439 
440  // Set time from database
441  #include "createTime.H"
442 
443  // Allow override of time (unless dry-run)
444  instantList times;
445  if (args.dryRun())
446  {
447  Info<< "\ndry-run: ignoring -copy*, -fields, -force, time selection"
448  << nl;
449  }
450  else
451  {
452  if (decomposeFieldsOnly && !doDecompFields)
453  {
455  << "Options -fields and -no-fields are mutually exclusive"
456  << " ... giving up" << nl
457  << exit(FatalError);
458  }
459 
460  if (!doDecompFields)
461  {
462  Info<< "Skip decompose of all fields" << nl;
463  }
464  if (!doFiniteArea)
465  {
466  Info<< "Skip decompose of finiteArea mesh/fields" << nl;
467  }
468  if (!doLagrangian)
469  {
470  Info<< "Skip decompose of lagrangian positions/fields" << nl;
471  }
472 
474  }
475 
476 
477  // Allow override of decomposeParDict location
478  fileName decompDictFile;
479  if
480  (
481  args.readIfPresent("decomposeParDict", decompDictFile)
482  && !decompDictFile.empty() && !decompDictFile.isAbsolute()
483  )
484  {
485  decompDictFile = runTime.globalPath()/decompDictFile;
486  }
487 
488  // Get region names
489  #include "getAllRegionOptions.H"
490 
491  const bool optRegions =
493 
495  {
496  Info<< "Using region: " << regionNames[0] << nl << endl;
497  }
498 
499  forAll(regionNames, regioni)
500  {
501  const word& regionName = regionNames[regioni];
503 
504  if (args.dryRun())
505  {
506  Info<< "dry-run: decomposing mesh " << regionName << nl << nl
507  << "Create mesh..." << flush;
508 
509  domainDecompositionDryRun decompTest
510  (
511  IOobject
512  (
513  regionName,
514  runTime.timeName(),
515  runTime,
519  ),
520  decompDictFile,
521  args.getOrDefault<label>("domains", 0),
522  args.getOrDefault<word>("method", word::null)
523  );
524 
525  decompTest.execute(writeCellDist, args.verbose());
526  continue;
527  }
528 
529  Info<< "\n\nDecomposing mesh";
530  if (!regionDir.empty())
531  {
532  Info<< ' ' << regionName;
533  }
534  Info<< nl << endl;
535 
536  // Determine the existing processor count directly
537  const label nProcsOld =
538  fileHandler().nProcs(runTime.path(), regionDir);
539 
540  // Get requested numberOfSubdomains directly from the dictionary.
541  // Note: have no mesh yet so cannot use decompositionModel::New
542  const label nDomains = decompositionMethod::nDomains
543  (
545  (
547  (
548  IOobject
549  (
551  runTime.time().system(),
552  regionDir, // region (if non-default)
553  runTime,
557  ),
558  decompDictFile
559  )
560  )
561  );
562 
563  // Give file handler a chance to determine the output directory
564  const_cast<fileOperation&>(fileHandler()).nProcs(nDomains);
565 
566  if (decomposeFieldsOnly)
567  {
568  // Sanity check on previously decomposed case
569  if (nProcsOld != nDomains)
570  {
572  << "Specified -fields, but the case was decomposed with "
573  << nProcsOld << " domains"
574  << nl
575  << "instead of " << nDomains
576  << " domains as specified in decomposeParDict" << nl
577  << exit(FatalError);
578  }
579  }
580  else if (nProcsOld)
581  {
582  bool procDirsProblem = true;
583 
584  if (decomposeIfRequired && nProcsOld == nDomains)
585  {
586  // We can reuse the decomposition
587  decomposeFieldsOnly = true;
588  procDirsProblem = false;
589  forceOverwrite = false;
590 
591  Info<< "Using existing processor directories" << nl;
592  }
593 
594  if (optRegions)
595  {
596  procDirsProblem = false;
597  forceOverwrite = false;
598  }
599 
600  if (forceOverwrite)
601  {
602  Info<< "Removing " << nProcsOld
603  << " existing processor directories" << endl;
604 
605  // Remove existing processors directory
606  fileNameList dirs
607  (
609  (
610  runTime.path(),
611  fileName::Type::DIRECTORY
612  )
613  );
614  forAllReverse(dirs, diri)
615  {
616  const fileName& d = dirs[diri];
617 
618  label proci = -1;
619 
620  if
621  (
622  d.starts_with("processor")
623  &&
624  (
625  // Collated is "processors"
626  d[9] == 's'
627 
628  // Uncollated has integer(s) after 'processor'
629  || Foam::read(d.substr(9), proci)
630  )
631  )
632  {
633  if (fileHandler().exists(d))
634  {
635  fileHandler().rmDir(d);
636  }
637  }
638  }
639 
640  procDirsProblem = false;
641  }
642 
643  if (procDirsProblem)
644  {
646  << "Case is already decomposed with " << nProcsOld
647  << " domains, use the -force option or manually" << nl
648  << "remove processor directories before decomposing. e.g.,"
649  << nl
650  << " rm -rf " << runTime.path().c_str() << "/processor*"
651  << nl
652  << exit(FatalError);
653  }
654  }
655 
656  Info<< "Create mesh" << endl;
658  (
659  IOobject
660  (
661  regionName,
662  runTime.timeName(),
663  runTime,
667  ),
668  decompDictFile
669  );
670  // Make sure pointMesh gets read as well
672 
673 
674  // Decompose the mesh
675  if (!decomposeFieldsOnly)
676  {
677  mesh.decomposeMesh();
678  mesh.writeDecomposition(decomposeSets);
679 
680  if (writeCellDist)
681  {
682  const labelList& procIds = mesh.cellToProc();
683 
684  // Write decomposition for visualization
685  mesh.writeVolField("cellDist");
686  //TBD: mesh.writeVTK("cellDist");
687 
688  // Write decomposition as labelList for use with 'manual'
689  // decomposition method.
690  labelIOList cellDecomposition
691  (
692  IOobject
693  (
694  "cellDecomposition",
696  mesh,
700  ),
701  procIds
702  );
703  cellDecomposition.write();
704 
705  Info<< nl << "Wrote decomposition to "
706  << cellDecomposition.objectRelPath()
707  << " for use in manual decomposition." << endl;
708  }
709 
710  fileHandler().flush();
711  }
712 
713 
714  if (copyZero)
715  {
716  // Copy the 0/ directory into each of the processor directories
717  // with fallback of 0.orig/ directory if necessary.
718 
719  fileName inputDir(runTime.path()/"0");
720 
721  bool canCopy = fileHandler().isDir(inputDir);
722  if (!canCopy)
723  {
724  // Try with "0.orig" instead
725  inputDir.ext("orig");
726  canCopy = fileHandler().isDir(inputDir);
727  }
728 
729  if (canCopy)
730  {
731  // Avoid copying into the same directory multiple times
732  // (collated format). Don't need a hash here.
733  fileName prevOutputDir;
734  for (label proci = 0; proci < mesh.nProcs(); ++proci)
735  {
736  Time processorDb
737  (
739  args.rootPath(),
740  args.caseName()/("processor" + Foam::name(proci)),
741  false, // No function objects
742  false // No extra controlDict libs
743  );
744  // processorDb.setTime(runTime);
745 
746  // Get corresponding directory name
747  // (to handle processors/)
748  const fileName outputDir
749  (
750  fileHandler().objectPath
751  (
752  IOobject
753  (
754  word::null, // name
755  "0", // instance (time == 0)
756  processorDb
757  ),
758  word::null
759  )
760  );
761 
762  if (outputDir != prevOutputDir)
763  {
764  Info<< "Processor " << proci
765  << ": copying \""
766  << inputDir.name() << "/\" to "
767  << runTime.relativePath(outputDir)
768  << endl;
769 
770  fileHandler().cp(inputDir, outputDir);
771  prevOutputDir = outputDir;
772  }
773  }
774  }
775  else
776  {
777  Info<< "No 0/ or 0.orig/ directory to copy" << nl;
778  }
779  }
780  else
781  {
782  // Decompose field files, lagrangian, finite-area
783 
784  // Cached processor meshes and maps. These are only preserved if
785  // running with multiple times.
786  PtrList<Time> processorDbList(mesh.nProcs());
787  PtrList<fvMesh> procMeshList(mesh.nProcs());
788  PtrList<labelIOList> faceProcAddressingList(mesh.nProcs());
789  PtrList<labelIOList> cellProcAddressingList(mesh.nProcs());
790  PtrList<labelIOList> boundaryProcAddressingList(mesh.nProcs());
791  PtrList<labelIOList> pointProcAddressingList(mesh.nProcs());
792  PtrList<labelIOList> pointBoundaryProcAddressingList(mesh.nProcs());
793 
794  PtrList<fvFieldDecomposer> fieldDecomposerList(mesh.nProcs());
795  PtrList<pointFieldDecomposer> pointFieldDecomposerList
796  (
797  mesh.nProcs()
798  );
799 
800 
801  // Loop over all times
802  forAll(times, timei)
803  {
804  runTime.setTime(times[timei], timei);
805 
806  Info<< "Time = " << runTime.timeName() << endl;
807 
808  // Field objects at this time
809  IOobjectList objects;
810  IOobjectList faObjects;
811 
812  if (doDecompFields)
813  {
814  // List of volume mesh objects for this instance
815  objects = IOobjectList(mesh, runTime.timeName());
816 
817  // List of area mesh objects (assuming single region)
818  faObjects = IOobjectList
819  (
820  mesh.time(),
821  runTime.timeName(),
824  );
825 
826  // Ignore generated fields: (cellDist)
827  objects.remove("cellDist");
828  }
829 
830  // Finite area handling
831  autoPtr<faMeshDecomposition> faMeshDecompPtr;
832  if (doFiniteArea)
833  {
834  const word boundaryInst =
835  mesh.time().findInstance(mesh.meshDir(), "boundary");
836 
837  IOobject io
838  (
839  "faBoundary",
840  boundaryInst,
842  mesh.time(),
846  );
847 
848  if (io.typeHeaderOk<faBoundaryMesh>(true))
849  {
850  // Always based on the volume decomposition!
851  faMeshDecompPtr.reset
852  (
854  (
855  mesh,
856  mesh.nProcs(),
857  mesh.model()
858  )
859  );
860  }
861  }
862 
863 
864  // Volume/surface/internal fields
865  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
866 
867  fvFieldDecomposer::fieldsCache volumeFieldCache;
868 
869  if (doDecompFields)
870  {
871  volumeFieldCache.readAllFields(mesh, objects);
872  }
873 
874 
875  // Point fields
876  // ~~~~~~~~~~~~
877 
878  // Read decomposed pointMesh
879  const pointMesh& pMesh =
881 
882  pointFieldDecomposer::fieldsCache pointFieldCache;
883 
884  if (doDecompFields)
885  {
886  pointFieldCache.readAllFields(pMesh, objects);
887  }
888 
889 
890  // Lagrangian fields
891  // ~~~~~~~~~~~~~~~~~
892 
893  fileNameList cloudDirs;
894 
895  if (doDecompFields && doLagrangian)
896  {
897  cloudDirs = fileHandler().readDir
898  (
901  );
902  }
903 
904  // Particles
905  PtrList<Cloud<indexedParticle>> lagrangianPositions
906  (
907  cloudDirs.size()
908  );
909  // Particles per cell
911  (
912  cloudDirs.size()
913  );
914 
915  lagrangianFieldDecomposer::fieldsCache lagrangianFieldCache
916  (
917  cloudDirs.size()
918  );
919 
920  label cloudI = 0;
921 
922  for (const fileName& cloudDir : cloudDirs)
923  {
924  IOobjectList cloudObjects
925  (
926  mesh,
927  runTime.timeName(),
928  cloud::prefix/cloudDir,
930  );
931 
932  // Note: look up "positions" for backwards compatibility
933  if
934  (
935  cloudObjects.found("coordinates")
936  || cloudObjects.found("positions")
937  )
938  {
939  // Read lagrangian particles
940  // ~~~~~~~~~~~~~~~~~~~~~~~~~
941 
942  Info<< "Identified lagrangian data set: "
943  << cloudDir << endl;
944 
945  lagrangianPositions.set
946  (
947  cloudI,
949  (
950  mesh,
951  cloudDir,
952  false
953  )
954  );
955 
956 
957  // Sort particles per cell
958  // ~~~~~~~~~~~~~~~~~~~~~~~
959 
960  cellParticles.set
961  (
962  cloudI,
964  (
965  mesh.nCells(),
966  static_cast<SLList<indexedParticle*>*>(nullptr)
967  )
968  );
969 
970  label i = 0;
971 
972  for (indexedParticle& p : lagrangianPositions[cloudI])
973  {
974  p.index() = i++;
975 
976  label celli = p.cell();
977 
978  // Check
979  if (celli < 0 || celli >= mesh.nCells())
980  {
982  << "Illegal cell number " << celli
983  << " for particle with index "
984  << p.index()
985  << " at position "
986  << p.position() << nl
987  << "Cell number should be between 0 and "
988  << mesh.nCells()-1 << nl
989  << "On this mesh the particle should"
990  << " be in cell "
991  << mesh.findCell(p.position())
992  << exit(FatalError);
993  }
994 
995  if (!cellParticles[cloudI][celli])
996  {
997  cellParticles[cloudI][celli] =
999  }
1000 
1001  cellParticles[cloudI][celli]->append(&p);
1002  }
1003 
1004  // Read fields
1005  // ~~~~~~~~~~~
1006 
1007  IOobjectList lagrangianObjects
1008  (
1009  mesh,
1010  runTime.timeName(),
1011  cloud::prefix/cloudDirs[cloudI],
1013  );
1014 
1015  lagrangianFieldCache.readAllFields
1016  (
1017  cloudI,
1018  lagrangianObjects
1019  );
1020 
1021  ++cloudI;
1022  }
1023  }
1024 
1025  lagrangianPositions.resize(cloudI);
1026  cellParticles.resize(cloudI);
1027  lagrangianFieldCache.resize(cloudI);
1028 
1029  Info<< endl;
1030 
1031  // split the fields over processors
1032  for
1033  (
1034  label proci = 0;
1035  doDecompFields && proci < mesh.nProcs();
1036  ++proci
1037  )
1038  {
1039  Info<< "Processor " << proci << ": field transfer" << endl;
1040 
1041  // open the database
1042  if (!processorDbList.set(proci))
1043  {
1044  processorDbList.set
1045  (
1046  proci,
1047  new Time
1048  (
1050  args.rootPath(),
1051  args.caseName()
1052  / ("processor" + Foam::name(proci)),
1054  args.allowLibs()
1055  )
1056  );
1057  }
1058  Time& processorDb = processorDbList[proci];
1059 
1060 
1061  processorDb.setTime(runTime);
1062 
1063  // read the mesh
1064  if (!procMeshList.set(proci))
1065  {
1066  procMeshList.set
1067  (
1068  proci,
1069  new fvMesh
1070  (
1071  IOobject
1072  (
1073  regionName,
1074  processorDb.timeName(),
1075  processorDb
1076  )
1077  )
1078  );
1079  }
1080  const fvMesh& procMesh = procMeshList[proci];
1081 
1082  const labelIOList& faceProcAddressing = procAddressing
1083  (
1084  procMeshList,
1085  proci,
1086  "faceProcAddressing",
1087  faceProcAddressingList
1088  );
1089 
1090  const labelIOList& cellProcAddressing = procAddressing
1091  (
1092  procMeshList,
1093  proci,
1094  "cellProcAddressing",
1095  cellProcAddressingList
1096  );
1097 
1098  const labelIOList& boundaryProcAddressing = procAddressing
1099  (
1100  procMeshList,
1101  proci,
1102  "boundaryProcAddressing",
1103  boundaryProcAddressingList
1104  );
1105 
1106 
1107  // FV fields: volume, surface, internal
1108  {
1109  if (!fieldDecomposerList.set(proci))
1110  {
1111  fieldDecomposerList.set
1112  (
1113  proci,
1114  new fvFieldDecomposer
1115  (
1116  mesh,
1117  procMesh,
1118  faceProcAddressing,
1119  cellProcAddressing,
1120  boundaryProcAddressing
1121  )
1122  );
1123  }
1124 
1125  volumeFieldCache.decomposeAllFields
1126  (
1127  fieldDecomposerList[proci]
1128  );
1129 
1130  if (times.size() == 1)
1131  {
1132  // Clear cached decomposer
1133  fieldDecomposerList.set(proci, nullptr);
1134  }
1135  }
1136 
1137 
1138  // Point fields
1139  if (!pointFieldCache.empty())
1140  {
1141  const labelIOList& pointProcAddressing = procAddressing
1142  (
1143  procMeshList,
1144  proci,
1145  "pointProcAddressing",
1146  pointProcAddressingList
1147  );
1148 
1149  const pointMesh& procPMesh =
1151 
1152  if (!pointBoundaryProcAddressingList.set(proci))
1153  {
1154  pointBoundaryProcAddressingList.set
1155  (
1156  proci,
1158  (
1159  IOobject
1160  (
1161  "boundaryProcAddressing",
1162  procMesh.facesInstance(),
1165  procPMesh.thisDb(),
1169  ),
1170  boundaryProcAddressing
1171  )
1172  );
1173  }
1174  const auto& pointBoundaryProcAddressing =
1175  pointBoundaryProcAddressingList[proci];
1176 
1177 
1178  if (!pointFieldDecomposerList.set(proci))
1179  {
1180  pointFieldDecomposerList.set
1181  (
1182  proci,
1184  (
1185  pMesh,
1186  procPMesh,
1187  pointProcAddressing,
1188  pointBoundaryProcAddressing
1189  )
1190  );
1191  }
1192 
1193  pointFieldCache.decomposeAllFields
1194  (
1195  pointFieldDecomposerList[proci]
1196  );
1197 
1198  if (times.size() == 1)
1199  {
1200  // Early deletion
1201  pointBoundaryProcAddressingList.set
1202  (
1203  proci,
1204  nullptr
1205  );
1206  pointProcAddressingList.set(proci, nullptr);
1207  pointFieldDecomposerList.set(proci, nullptr);
1208  }
1209  }
1210 
1211 
1212  // If there is lagrangian data write it out
1213  forAll(lagrangianPositions, cloudi)
1214  {
1215  if (lagrangianPositions[cloudi].size())
1216  {
1217  lagrangianFieldDecomposer fieldDecomposer
1218  (
1219  mesh,
1220  procMesh,
1221  faceProcAddressing,
1222  cellProcAddressing,
1223  cloudDirs[cloudi],
1224  lagrangianPositions[cloudi],
1225  cellParticles[cloudi]
1226  );
1227 
1228  // Lagrangian fields
1229  lagrangianFieldCache.decomposeAllFields
1230  (
1231  cloudi,
1232  cloudDirs[cloudi],
1233  fieldDecomposer
1234  );
1235  }
1236  }
1237 
1238  if (doDecompFields)
1239  {
1240  // Decompose "uniform" directory in the time region
1241  // directory
1242  decomposeUniform
1243  (
1244  copyUniform, mesh, processorDb, regionDir
1245  );
1246 
1247  // For a multi-region case, also decompose "uniform"
1248  // directory in the time directory
1249  if (regionNames.size() > 1 && regioni == 0)
1250  {
1251  decomposeUniform(copyUniform, mesh, processorDb);
1252  }
1253  }
1254 
1255 
1256  // We have cached all the constant mesh data for the current
1257  // processor. This is only important if running with
1258  // multiple times, otherwise it is just extra storage.
1259  if (times.size() == 1)
1260  {
1261  boundaryProcAddressingList.set(proci, nullptr);
1262  cellProcAddressingList.set(proci, nullptr);
1263  faceProcAddressingList.set(proci, nullptr);
1264  procMeshList.set(proci, nullptr);
1265  processorDbList.set(proci, nullptr);
1266  }
1267  }
1268 
1269 
1270  // Finite area mesh and field decomposition
1271  if (faMeshDecompPtr)
1272  {
1273  Info<< "\nFinite area mesh decomposition" << endl;
1274 
1275  faMeshDecomposition& aMesh = faMeshDecompPtr();
1276 
1277  aMesh.decomposeMesh();
1278  aMesh.writeDecomposition();
1279 
1280 
1281  // Area/edge fields
1282  // ~~~~~~~~~~~~~~~~
1283 
1284  faFieldDecomposer::fieldsCache areaFieldCache;
1285 
1286  if (doDecompFields)
1287  {
1288  areaFieldCache.readAllFields(aMesh, faObjects);
1289  }
1290 
1291  const label nAreaFields = areaFieldCache.size();
1292 
1293  Info<< endl;
1294  Info<< "Finite area field transfer: "
1295  << nAreaFields << " fields" << endl;
1296 
1297  // Split the fields over processors
1298  for
1299  (
1300  label proci = 0;
1301  nAreaFields && proci < mesh.nProcs();
1302  ++proci
1303  )
1304  {
1305  Info<< " Processor " << proci << endl;
1306 
1307  // open the database
1308  Time processorDb
1309  (
1311  args.rootPath(),
1312  args.caseName()/("processor" + Foam::name(proci)),
1313  false, // No function objects
1314  false // No extra controlDict libs
1315  );
1316 
1317  processorDb.setTime(runTime);
1318 
1319  // Read the mesh
1320  fvMesh procFvMesh
1321  (
1322  IOobject
1323  (
1324  regionName,
1325  processorDb.timeName(),
1326  processorDb
1327  )
1328  );
1329 
1330  faMesh procMesh(procFvMesh);
1331 
1332  // // Does not work. HJ, 15/Aug/2017
1333  // const labelIOList& faceProcAddressing =
1334  // procAddressing
1335  // (
1336  // procMeshList,
1337  // proci,
1338  // "faceProcAddressing",
1339  // faceProcAddressingList
1340  // );
1341 
1342  // const labelIOList& boundaryProcAddressing =
1343  // procAddressing
1344  // (
1345  // procMeshList,
1346  // proci,
1347  // "boundaryProcAddressing",
1348  // boundaryProcAddressingList
1349  // );
1350 
1351  // Addressing from faMesh (not polyMesh) meshSubDir
1352 
1353  autoPtr<labelIOList> tfaceProcAddr =
1354  faProcAddressing
1355  (
1356  procMesh,
1357  "faceProcAddressing",
1358  runTime.constant()
1359  );
1360  auto& faceProcAddressing = *tfaceProcAddr;
1361 
1362  autoPtr<labelIOList> tboundaryProcAddr =
1363  faProcAddressing
1364  (
1365  procMesh,
1366  "boundaryProcAddressing",
1367  runTime.constant()
1368  );
1369  auto& boundaryProcAddressing = *tboundaryProcAddr;
1370 
1371  autoPtr<labelIOList> tedgeProcAddr =
1372  faProcAddressing
1373  (
1374  procMesh,
1375  "edgeProcAddressing",
1376  runTime.constant()
1377  );
1378  const auto& edgeProcAddressing = *tedgeProcAddr;
1379 
1380  faFieldDecomposer fieldDecomposer
1381  (
1382  aMesh,
1383  procMesh,
1384  edgeProcAddressing,
1385  faceProcAddressing,
1386  boundaryProcAddressing
1387  );
1388 
1389  areaFieldCache.decomposeAllFields(fieldDecomposer);
1390  }
1391  }
1392  }
1393  }
1394  }
1395 
1396  Info<< "\nEnd\n" << endl;
1397 
1398  return 0;
1399 }
1400 
1401 
1402 // ************************************************************************* //
Finite area mesh (used for 2-D non-Euclidian finite area method) defined using a patch of faces on a ...
Definition: faMesh.H:133
static void addNote(const string &note)
Add extra notes for the usage information.
Definition: argList.C:476
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
fileName path() const
Return path = rootPath/caseName. Same as TimePaths::path()
Definition: Time.H:503
A class for handling file names.
Definition: fileName.H:72
Finite Area area and edge field decomposer.
List of IOobjects with searching and retrieving facilities. Implemented as a HashTable, so the various sorted methods should be used if traversing in parallel.
Definition: IOobjectList.H:55
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
const fileName & facesInstance() const
Return the current instance directory for faces.
Definition: polyMesh.C:844
void append(const T &elem)
Add copy at back of list.
Definition: LList.H:744
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
fileName timePath() const
Return current time path = path/timeName.
Definition: Time.H:520
const word & regionName() const
The mesh region name or word::null if polyMesh::defaultRegion.
Definition: polyMesh.C:832
void readAllFields(const faMesh &mesh, const IOobjectList &objects)
Read all fields given mesh and objects.
Template class for non-intrusive linked lists.
Definition: LList.H:46
static FOAM_NO_DANGLING_REFERENCE const pointMesh & New(const polyMesh &mesh, Args &&... args)
Get existing or create MeshObject registered with typeName.
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
Definition: polyMesh.H:411
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:529
bool chDir(const fileName &dir)
Change current directory to the one specified and return true on success.
Definition: POSIX.C:609
fileName relativePath(const fileName &input, const bool caseTag=false) const
Return the input relative to the globalPath by stripping off a leading value of the globalPath...
Definition: TimePathsI.H:103
static void addBoolOption(const word &optName, const string &usage="", bool advanced=false)
Add a bool option to validOptions with usage information.
Definition: argList.C:388
static void noParallel()
Remove the parallel options.
Definition: argList.C:598
static void addOptionCompat(const word &optName, std::pair< const char *, int > compat)
Specify an alias for the option name.
Definition: argList.C:432
refPtr< fileOperation > fileHandler(std::nullptr_t)
Delete current file handler - forwards to fileOperation::handler()
wordList regionNames
static bool isAbsolute(const std::string &str)
Return true if filename starts with a &#39;/&#39; or &#39;\&#39; or (windows-only) with a filesystem-root.
Definition: fileNameI.H:129
Ignore writing from objectRegistry::writeObject()
void decomposeAllFields(const fvFieldDecomposer &decomposer, bool report=false) const
Decompose and write all fields.
T getOrDefault(const word &optName, const T &deflt) const
Get a value from the named option if present, or return default.
Definition: argListI.H:300
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:360
Automatic domain decomposition class for finite-volume meshes.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:69
An encapsulation of filesystem-related operations.
bool allowFunctionObjects() const
The controlDict &#39;functions&#39; entry is allowed to be used.
Definition: argList.C:2252
#define FOAM_NO_DANGLING_REFERENCE
Definition: stdFoam.H:80
void readAllFields(const pointMesh &mesh, const IOobjectList &objects)
Read all fields given mesh and objects.
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition: int32.H:127
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:50
Mesh representing a set of points created from polyMesh.
Definition: pointMesh.H:45
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:286
void reset(T *p=nullptr) noexcept
Delete managed object and set to new given pointer.
Definition: autoPtrI.H:37
bool isDir(const fileName &name, const bool followLink=true)
Does the name exist as a DIRECTORY in the file system?
Definition: POSIX.C:862
virtual const fileName & dbDir() const
Local directory path of the objectRegistry relative to Time with override for the single-region case...
Definition: faMesh.H:1212
Functions used by OpenFOAM that are specific to POSIX compliant operating systems and need to be repl...
bool allowLibs() const
The controlDict &#39;libs&#39; entry is allowed to be used. (eg, has not been disabled by the -no-libs option...
Definition: argList.C:2270
void readAllFields(const fvMesh &mesh, const IOobjectList &objects)
Read all fields given mesh and objects.
dynamicFvMesh & mesh
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
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
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
fileName globalPath() const
Return global path for the case = rootPath/globalCaseName. Same as TimePaths::globalPath() ...
Definition: Time.H:512
A class for handling words, derived from Foam::string.
Definition: word.H:63
static void addDryRunOption(const string &usage, bool advanced=false)
Enable a &#39;dry-run&#39; bool option, with usage information.
Definition: argList.C:518
fileName meshDir() const
Return the local mesh directory (dbDir()/meshSubDir)
Definition: faMesh.C:1015
const Time & time() const noexcept
Return time registry.
const word & executable() const noexcept
Name of executable without the path.
Definition: argListI.H:44
int dryRun() const noexcept
Return the dry-run flag.
Definition: argListI.H:109
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:406
bool writeDecomposition()
Write decomposition.
const word & system() const noexcept
Return system name.
Definition: TimePathsI.H:118
const objectRegistry & thisDb() const
Return database. For now is its polyMesh.
Definition: pointMesh.H:201
Reading is optional [identical to LAZY_READ].
static const word null
An empty word.
Definition: word.H:84
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
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:399
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...
Definition: timeSelector.C:291
virtual void setTime(const Time &t)
Reset the time and time-index to those of the given time.
Definition: Time.C:902
bool local
Definition: EEqn.H:20
static word controlDictName
The default control dictionary name (normally "controlDict")
Definition: Time.H:268
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:159
const fileName & caseName() const noexcept
Return case name (parallel run) or global case (serial run)
Definition: argListI.H:62
label nDomains() const noexcept
Number of domains.
bool typeHeaderOk(const bool checkType=true, const bool search=true, const bool verbose=true)
Read header (respects is_globalIOobject trait) and check its info.
void resize(label newCapacity)
Rehash the hash table with new number of buckets. Currently identical to setCapacity() ...
Definition: HashTable.C:722
static void addVerboseOption(const string &usage="", bool advanced=false)
Enable a &#39;verbose&#39; bool option, with usage information.
Definition: argList.C:534
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
const word & constant() const noexcept
Return constant name.
Definition: TimePathsI.H:112
Point field decomposer.
const fileName & rootPath() const noexcept
Return root path.
Definition: argListI.H:56
T & try_emplace(const label i, Args &&... args)
Like emplace_set() but will not overwrite an occupied (non-null) location.
Definition: PtrListI.H:193
Finite Volume volume and surface field decomposer.
static const word canonicalName
The canonical name ("decomposeParDict") under which the MeshObject is registered. ...
void decomposeAllFields(const pointFieldDecomposer &decomposer, bool report=false) const
Decompose and write all fields.
static word meshSubDir
Return the mesh sub-directory name (usually "pointMesh")
Definition: pointMesh.H:107
static word meshSubDir
The mesh sub-directory name (usually "faMesh")
Definition: faMesh.H:750
bool starts_with(char c) const
True if string starts with given character (cf. C++20)
Definition: string.H:435
void decomposeMesh()
Decompose mesh.
Ostream & flush(Ostream &os)
Flush stream.
Definition: Ostream.H:519
Foam::word regionName(args.getOrDefault< word >("region", Foam::polyMesh::defaultRegion))
label nCells() const noexcept
Number of mesh cells.
autoPtr< IOobject > remove(const IOobject &io)
Remove object from the list by its IOobject::name().
Definition: IOobjectList.H:302
const word & regionDir
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers...
Definition: List.H:54
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
fileName meshDir() const
Return the local mesh directory (dbDir()/meshSubDir)
Definition: polyMesh.C:826
Nothing to be read.
Required Classes.
fileName cwd()
The physical or logical current working directory path name.
Definition: POSIX.C:592
Finite area boundary mesh, which is a faPatch list with registered IO, a reference to the associated ...
label findCell(const point &p, const cellDecomposition=CELL_TETS) const
Find cell enclosing this location and return index.
Definition: polyMesh.C:1496
messageStream Info
Information stream (stdout output on master, null elsewhere)
Automatic faMesh decomposition class.
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
#define forAllReverse(list, i)
Reverse loop across all elements in list.
Definition: stdFoam.H:302
bool readIfPresent(const word &optName, T &val) const
Read a value from the named option if present.
Definition: argListI.H:316
volScalarField & p
IOobject io("surfaceFilmProperties", mesh.time().constant(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE, IOobject::NO_REGISTER)
Registry of regIOobjects.
static autoPtr< T > New(Args &&... args)
Construct autoPtr with forwarding arguments.
Definition: autoPtr.H:178
Testing of domain decomposition for finite-volume meshes.
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
Foam::argList args(argc, argv)
label size() const
Number of fields.
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:180
A class for handling character strings derived from std::string.
Definition: string.H:72
Adds label index to base particle.
#define FatalErrorIn(functionName)
Report an error message using Foam::FatalError.
Definition: error.H:592
int verbose() const noexcept
Return the verbose flag.
Definition: argListI.H:121
Do not request registration (bool: false)
bool found(const word &optName) const
Return true if the named option is found.
Definition: argListI.H:171
static void addOptions(const bool constant=true, const bool withZero=false)
Add timeSelector options to argList::validOptions.
Definition: timeSelector.C:101
bool set(const Key &key, const T &obj)
Copy assign a new entry, overwriting existing entries.
Definition: HashTableI.H:174
Namespace for OpenFOAM.
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
void decomposeAllFields(const faFieldDecomposer &decomposer, bool report=false) const
Decompose and write all fields.
Lagrangian field decomposer.
static const word prefix
The prefix to local: lagrangian.
Definition: cloud.H:79