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 -------------------------------------------------------------------------------
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 fvMesh& procMesh,
177  const word& name,
178  const word& instance,
179  const word& local = fvMesh::meshSubDir
180 )
181 {
183  (
184  IOobject
185  (
186  name,
187  instance,
188  local,
189  procMesh,
193  )
194  );
195 }
196 
197 
198 // Read proc addressing at specific instance.
199 // Uses the finiteArea meshSubDir
200 autoPtr<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
213 const 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 
235 void decomposeUniform
236 (
237  const bool copyUniform,
238  const domainDecomposition& mesh,
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  // Bit of trickery to synthesise the correct directory base,
255  // e.g. processors4/0.01
256  const fileName timePath = fileHandler().objectPath
257  (
258  IOobject
259  (
260  "dummy",
261  runTime.timeName(),
262  processorDb
263  ),
264  word::null
265  ).path();
266 
267  // If no fields have been decomposed the destination
268  // directory will not have been created so make sure.
269  mkDir(timePath);
270 
271  if (copyUniform || mesh.distributed())
272  {
273  if (!fileHandler().exists(timePath/uniformDir))
274  {
275  fileHandler().cp
276  (
277  runTime.timePath()/uniformDir,
278  timePath/uniformDir
279  );
280  }
281  }
282  else
283  {
284  // Link with relative paths
285  string parentPath = string("..")/"..";
286 
287  if (!regionDir.empty())
288  {
289  parentPath = parentPath/"..";
290  }
291 
292  fileName currentDir(cwd());
293  chDir(timePath);
294 
295  if (!fileHandler().exists(uniformDir))
296  {
297  fileHandler().ln
298  (
299  parentPath/runTime.timeName()/uniformDir,
300  uniformDir
301  );
302  }
303  chDir(currentDir);
304  }
305  }
306 }
307 
308 } // End namespace Foam
309 
310 
311 using namespace Foam;
312 
313 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
314 
315 int main(int argc, char *argv[])
316 {
318  (
319  "Decompose a mesh and fields of a case for parallel execution"
320  );
321 
324  (
325  "decomposeParDict",
326  "file",
327  "Alternative decomposePar dictionary file"
328  );
329 
330  #include "addAllRegionOptions.H"
331 
333  (
334  "Test without writing the decomposition. "
335  "Changes -cellDist to only write VTK output."
336  );
339  (
340  "domains",
341  "N",
342  "Override numberOfSubdomains (-dry-run only)",
343  true // Advanced option
344  );
346  (
347  "method",
348  "name",
349  "Override decomposition method (-dry-run only)",
350  true // Advanced option
351  );
352 
354  (
355  "no-finite-area",
356  "Suppress finiteArea mesh/field decomposition",
357  true // Advanced option
358  );
359 
361  (
362  "no-lagrangian",
363  "Suppress lagrangian (cloud) decomposition",
364  true // Advanced option
365  );
366 
368  (
369  "cellDist",
370  "Write cell distribution as a labelList - for use with 'manual' "
371  "decomposition method and as a volScalarField for visualization."
372  );
374  (
375  "copyZero",
376  "Copy 0/ directory to processor*/ rather than decompose the fields"
377  );
379  (
380  "copyUniform",
381  "Copy any uniform/ directories too"
382  );
384  (
385  "fields",
386  "Use existing geometry decomposition and convert fields only"
387  );
389  (
390  "no-fields",
391  "Suppress conversion of fields (volume, finite-area, lagrangian)"
392  );
393 
395  (
396  "no-sets",
397  "Skip decomposing cellSets, faceSets, pointSets"
398  );
399  argList::addOptionCompat("no-sets", {"noSets", 2106});
400 
402  (
403  "force",
404  "Remove existing processor*/ subdirs before decomposing the geometry"
405  );
407  (
408  "ifRequired",
409  "Only decompose geometry if the number of domains has changed"
410  );
411 
412  // Allow explicit -constant, have zero from time range
413  timeSelector::addOptions(true, false); // constant(true), zero(false)
414 
415  #include "setRootCase.H"
416 
417  const bool writeCellDist = args.found("cellDist");
418 
419  // Most of these are ignored for dry-run (not triggered anywhere)
420  const bool copyZero = args.found("copyZero");
421  const bool copyUniform = args.found("copyUniform");
422  const bool decomposeSets = !args.found("no-sets");
423 
424  const bool decomposeIfRequired = args.found("ifRequired");
425 
426  const bool doDecompFields = !args.found("no-fields");
427  const bool doFiniteArea = !args.found("no-finite-area");
428  const bool doLagrangian = !args.found("no-lagrangian");
429 
430  bool decomposeFieldsOnly = args.found("fields");
431  bool forceOverwrite = args.found("force");
432 
433 
434  // Set time from database
435  #include "createTime.H"
436 
437  // Allow override of time (unless dry-run)
438  instantList times;
439  if (args.dryRun())
440  {
441  Info<< "\ndry-run: ignoring -copy*, -fields, -force, time selection"
442  << nl;
443  }
444  else
445  {
446  if (decomposeFieldsOnly && !doDecompFields)
447  {
449  << "Options -fields and -no-fields are mutually exclusive"
450  << " ... giving up" << nl
451  << exit(FatalError);
452  }
453 
454  if (!doDecompFields)
455  {
456  Info<< "Skip decompose of all fields" << nl;
457  }
458  if (!doFiniteArea)
459  {
460  Info<< "Skip decompose of finiteArea mesh/fields" << nl;
461  }
462  if (!doLagrangian)
463  {
464  Info<< "Skip decompose of lagrangian positions/fields" << nl;
465  }
466 
468  }
469 
470 
471  // Allow override of decomposeParDict location
472  fileName decompDictFile;
473  if
474  (
475  args.readIfPresent("decomposeParDict", decompDictFile)
476  && !decompDictFile.empty() && !decompDictFile.isAbsolute()
477  )
478  {
479  decompDictFile = runTime.globalPath()/decompDictFile;
480  }
481 
482  // Get region names
483  #include "getAllRegionOptions.H"
484 
485  const bool optRegions =
487 
489  {
490  Info<< "Using region: " << regionNames[0] << nl << endl;
491  }
492 
493  forAll(regionNames, regioni)
494  {
495  const word& regionName = regionNames[regioni];
497 
498  if (args.dryRun())
499  {
500  Info<< "dry-run: decomposing mesh " << regionName << nl << nl
501  << "Create mesh..." << flush;
502 
503  domainDecompositionDryRun decompTest
504  (
505  IOobject
506  (
507  regionName,
508  runTime.timeName(),
509  runTime,
513  ),
514  decompDictFile,
515  args.getOrDefault<label>("domains", 0),
516  args.getOrDefault<word>("method", word::null)
517  );
518 
519  decompTest.execute(writeCellDist, args.verbose());
520  continue;
521  }
522 
523  Info<< "\n\nDecomposing mesh";
524  if (!regionDir.empty())
525  {
526  Info<< ' ' << regionName;
527  }
528  Info<< nl << endl;
529 
530  // Determine the existing processor count directly
531  const label nProcsOld =
532  fileHandler().nProcs(runTime.path(), regionDir);
533 
534  // Get requested numberOfSubdomains directly from the dictionary.
535  // Note: have no mesh yet so cannot use decompositionModel::New
536  const label nDomains = decompositionMethod::nDomains
537  (
539  (
541  (
542  IOobject
543  (
545  runTime.time().system(),
546  regionDir, // region (if non-default)
547  runTime,
551  ),
552  decompDictFile
553  )
554  )
555  );
556 
557  // Give file handler a chance to determine the output directory
558  const_cast<fileOperation&>(fileHandler()).nProcs(nDomains);
559 
560  if (decomposeFieldsOnly)
561  {
562  // Sanity check on previously decomposed case
563  if (nProcsOld != nDomains)
564  {
566  << "Specified -fields, but the case was decomposed with "
567  << nProcsOld << " domains"
568  << nl
569  << "instead of " << nDomains
570  << " domains as specified in decomposeParDict" << nl
571  << exit(FatalError);
572  }
573  }
574  else if (nProcsOld)
575  {
576  bool procDirsProblem = true;
577 
578  if (decomposeIfRequired && nProcsOld == nDomains)
579  {
580  // We can reuse the decomposition
581  decomposeFieldsOnly = true;
582  procDirsProblem = false;
583  forceOverwrite = false;
584 
585  Info<< "Using existing processor directories" << nl;
586  }
587 
588  if (optRegions)
589  {
590  procDirsProblem = false;
591  forceOverwrite = false;
592  }
593 
594  if (forceOverwrite)
595  {
596  Info<< "Removing " << nProcsOld
597  << " existing processor directories" << endl;
598 
599  // Remove existing processors directory
600  fileNameList dirs
601  (
603  (
604  runTime.path(),
605  fileName::Type::DIRECTORY
606  )
607  );
608  forAllReverse(dirs, diri)
609  {
610  const fileName& d = dirs[diri];
611 
612  label proci = -1;
613 
614  if
615  (
616  d.starts_with("processor")
617  &&
618  (
619  // Collated is "processors"
620  d[9] == 's'
621 
622  // Uncollated has integer(s) after 'processor'
623  || Foam::read(d.substr(9), proci)
624  )
625  )
626  {
627  if (fileHandler().exists(d))
628  {
629  fileHandler().rmDir(d);
630  }
631  }
632  }
633 
634  procDirsProblem = false;
635  }
636 
637  if (procDirsProblem)
638  {
640  << "Case is already decomposed with " << nProcsOld
641  << " domains, use the -force option or manually" << nl
642  << "remove processor directories before decomposing. e.g.,"
643  << nl
644  << " rm -rf " << runTime.path().c_str() << "/processor*"
645  << nl
646  << exit(FatalError);
647  }
648  }
649 
650  Info<< "Create mesh" << endl;
652  (
653  IOobject
654  (
655  regionName,
656  runTime.timeName(),
657  runTime,
661  ),
662  decompDictFile
663  );
664 
665  // Decompose the mesh
666  if (!decomposeFieldsOnly)
667  {
668  mesh.decomposeMesh();
669  mesh.writeDecomposition(decomposeSets);
670 
671  if (writeCellDist)
672  {
673  const labelList& procIds = mesh.cellToProc();
674 
675  // Write decomposition for visualization
676  mesh.writeVolField("cellDist");
677  //TBD: mesh.writeVTK("cellDist");
678 
679  // Write decomposition as labelList for use with 'manual'
680  // decomposition method.
681  labelIOList cellDecomposition
682  (
683  IOobject
684  (
685  "cellDecomposition",
687  mesh,
691  ),
692  procIds
693  );
694  cellDecomposition.write();
695 
696  Info<< nl << "Wrote decomposition to "
697  << cellDecomposition.objectRelPath()
698  << " for use in manual decomposition." << endl;
699  }
700 
701  fileHandler().flush();
702  }
703 
704 
705  if (copyZero)
706  {
707  // Copy the 0/ directory into each of the processor directories
708  // with fallback of 0.orig/ directory if necessary.
709 
710  fileName inputDir(runTime.path()/"0");
711 
712  bool canCopy = fileHandler().isDir(inputDir);
713  if (!canCopy)
714  {
715  // Try with "0.orig" instead
716  inputDir.ext("orig");
717  canCopy = fileHandler().isDir(inputDir);
718  }
719 
720  if (canCopy)
721  {
722  // Avoid copying into the same directory multiple times
723  // (collated format). Don't need a hash here.
724  fileName prevOutputDir;
725  for (label proci = 0; proci < mesh.nProcs(); ++proci)
726  {
727  Time processorDb
728  (
730  args.rootPath(),
731  args.caseName()/("processor" + Foam::name(proci)),
732  false, // No function objects
733  false // No extra controlDict libs
734  );
735  // processorDb.setTime(runTime);
736 
737  // Get corresponding directory name
738  // (to handle processors/)
739  const fileName outputDir
740  (
741  fileHandler().objectPath
742  (
743  IOobject
744  (
745  word::null, // name
746  "0", // instance (time == 0)
747  processorDb
748  ),
749  word::null
750  )
751  );
752 
753  if (outputDir != prevOutputDir)
754  {
755  Info<< "Processor " << proci
756  << ": copying \""
757  << inputDir.name() << "/\" to "
758  << runTime.relativePath(outputDir)
759  << endl;
760 
761  fileHandler().cp(inputDir, outputDir);
762  prevOutputDir = outputDir;
763  }
764  }
765  }
766  else
767  {
768  Info<< "No 0/ or 0.orig/ directory to copy" << nl;
769  }
770  }
771  else
772  {
773  // Decompose field files, lagrangian, finite-area
774 
775  // Cached processor meshes and maps. These are only preserved if
776  // running with multiple times.
777  PtrList<Time> processorDbList(mesh.nProcs());
778  PtrList<fvMesh> procMeshList(mesh.nProcs());
779  PtrList<labelIOList> faceProcAddressingList(mesh.nProcs());
780  PtrList<labelIOList> cellProcAddressingList(mesh.nProcs());
781  PtrList<labelIOList> boundaryProcAddressingList(mesh.nProcs());
782  PtrList<labelIOList> pointProcAddressingList(mesh.nProcs());
783 
784  PtrList<fvFieldDecomposer> fieldDecomposerList(mesh.nProcs());
785  PtrList<pointFieldDecomposer> pointFieldDecomposerList
786  (
787  mesh.nProcs()
788  );
789 
790 
791  // Loop over all times
792  forAll(times, timei)
793  {
794  runTime.setTime(times[timei], timei);
795 
796  Info<< "Time = " << runTime.timeName() << endl;
797 
798  // Field objects at this time
799  IOobjectList objects;
800 
801  if (doDecompFields)
802  {
803  objects = IOobjectList(mesh, runTime.timeName());
804 
805  // Ignore generated fields: (cellDist)
806  objects.remove("cellDist");
807  }
808 
809  // Finite area handling
810  autoPtr<faMeshDecomposition> faMeshDecompPtr;
811  if (doFiniteArea)
812  {
813  IOobject io
814  (
815  "faBoundary",
816  mesh.time().findInstance(mesh.meshDir(), "boundary"),
818  mesh,
822  );
823 
824  if (io.typeHeaderOk<faBoundaryMesh>(true))
825  {
826  // Always based on the volume decomposition!
827  faMeshDecompPtr.reset
828  (
830  (
831  mesh,
832  mesh.nProcs(),
833  mesh.model()
834  )
835  );
836  }
837  }
838 
839 
840  // Volume/surface/internal fields
841  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
842 
843  fvFieldDecomposer::fieldsCache volumeFieldCache;
844 
845  if (doDecompFields)
846  {
847  volumeFieldCache.readAllFields(mesh, objects);
848  }
849 
850 
851  // Point fields
852  // ~~~~~~~~~~~~
853  const pointMesh& pMesh = pointMesh::New(mesh);
854 
855  pointFieldDecomposer::fieldsCache pointFieldCache;
856 
857  if (doDecompFields)
858  {
859  pointFieldCache.readAllFields(pMesh, objects);
860  }
861 
862 
863  // Lagrangian fields
864  // ~~~~~~~~~~~~~~~~~
865 
866  fileNameList cloudDirs;
867 
868  if (doDecompFields && doLagrangian)
869  {
870  cloudDirs = fileHandler().readDir
871  (
874  );
875  }
876 
877  // Particles
878  PtrList<Cloud<indexedParticle>> lagrangianPositions
879  (
880  cloudDirs.size()
881  );
882  // Particles per cell
884  (
885  cloudDirs.size()
886  );
887 
888  lagrangianFieldDecomposer::fieldsCache lagrangianFieldCache
889  (
890  cloudDirs.size()
891  );
892 
893  label cloudI = 0;
894 
895  for (const fileName& cloudDir : cloudDirs)
896  {
897  IOobjectList cloudObjects
898  (
899  mesh,
900  runTime.timeName(),
901  cloud::prefix/cloudDir,
903  );
904 
905  // Note: look up "positions" for backwards compatibility
906  if
907  (
908  cloudObjects.found("coordinates")
909  || cloudObjects.found("positions")
910  )
911  {
912  // Read lagrangian particles
913  // ~~~~~~~~~~~~~~~~~~~~~~~~~
914 
915  Info<< "Identified lagrangian data set: "
916  << cloudDir << endl;
917 
918  lagrangianPositions.set
919  (
920  cloudI,
922  (
923  mesh,
924  cloudDir,
925  false
926  )
927  );
928 
929 
930  // Sort particles per cell
931  // ~~~~~~~~~~~~~~~~~~~~~~~
932 
933  cellParticles.set
934  (
935  cloudI,
937  (
938  mesh.nCells(),
939  static_cast<SLList<indexedParticle*>*>(nullptr)
940  )
941  );
942 
943  label i = 0;
944 
945  for (indexedParticle& p : lagrangianPositions[cloudI])
946  {
947  p.index() = i++;
948 
949  label celli = p.cell();
950 
951  // Check
952  if (celli < 0 || celli >= mesh.nCells())
953  {
955  << "Illegal cell number " << celli
956  << " for particle with index "
957  << p.index()
958  << " at position "
959  << p.position() << nl
960  << "Cell number should be between 0 and "
961  << mesh.nCells()-1 << nl
962  << "On this mesh the particle should"
963  << " be in cell "
964  << mesh.findCell(p.position())
965  << exit(FatalError);
966  }
967 
968  if (!cellParticles[cloudI][celli])
969  {
970  cellParticles[cloudI][celli] =
972  }
973 
974  cellParticles[cloudI][celli]->append(&p);
975  }
976 
977  // Read fields
978  // ~~~~~~~~~~~
979 
980  IOobjectList lagrangianObjects
981  (
982  mesh,
983  runTime.timeName(),
984  cloud::prefix/cloudDirs[cloudI],
986  );
987 
988  lagrangianFieldCache.readAllFields
989  (
990  cloudI,
991  lagrangianObjects
992  );
993 
994  ++cloudI;
995  }
996  }
997 
998  lagrangianPositions.resize(cloudI);
999  cellParticles.resize(cloudI);
1000  lagrangianFieldCache.resize(cloudI);
1001 
1002  Info<< endl;
1003 
1004  // split the fields over processors
1005  for
1006  (
1007  label proci = 0;
1008  doDecompFields && proci < mesh.nProcs();
1009  ++proci
1010  )
1011  {
1012  Info<< "Processor " << proci << ": field transfer" << endl;
1013 
1014  // open the database
1015  if (!processorDbList.set(proci))
1016  {
1017  processorDbList.set
1018  (
1019  proci,
1020  new Time
1021  (
1023  args.rootPath(),
1024  args.caseName()
1025  / ("processor" + Foam::name(proci)),
1027  args.allowLibs()
1028  )
1029  );
1030  }
1031  Time& processorDb = processorDbList[proci];
1032 
1033 
1034  processorDb.setTime(runTime);
1035 
1036  // read the mesh
1037  if (!procMeshList.set(proci))
1038  {
1039  procMeshList.set
1040  (
1041  proci,
1042  new fvMesh
1043  (
1044  IOobject
1045  (
1046  regionName,
1047  processorDb.timeName(),
1048  processorDb
1049  )
1050  )
1051  );
1052  }
1053  const fvMesh& procMesh = procMeshList[proci];
1054 
1055  const labelIOList& faceProcAddressing = procAddressing
1056  (
1057  procMeshList,
1058  proci,
1059  "faceProcAddressing",
1060  faceProcAddressingList
1061  );
1062 
1063  const labelIOList& cellProcAddressing = procAddressing
1064  (
1065  procMeshList,
1066  proci,
1067  "cellProcAddressing",
1068  cellProcAddressingList
1069  );
1070 
1071  const labelIOList& boundaryProcAddressing = procAddressing
1072  (
1073  procMeshList,
1074  proci,
1075  "boundaryProcAddressing",
1076  boundaryProcAddressingList
1077  );
1078 
1079 
1080  // FV fields: volume, surface, internal
1081  {
1082  if (!fieldDecomposerList.set(proci))
1083  {
1084  fieldDecomposerList.set
1085  (
1086  proci,
1087  new fvFieldDecomposer
1088  (
1089  mesh,
1090  procMesh,
1091  faceProcAddressing,
1092  cellProcAddressing,
1093  boundaryProcAddressing
1094  )
1095  );
1096  }
1097 
1098  volumeFieldCache.decomposeAllFields
1099  (
1100  fieldDecomposerList[proci]
1101  );
1102 
1103  if (times.size() == 1)
1104  {
1105  // Clear cached decomposer
1106  fieldDecomposerList.set(proci, nullptr);
1107  }
1108  }
1109 
1110 
1111  // Point fields
1112  if (!pointFieldCache.empty())
1113  {
1114  const labelIOList& pointProcAddressing = procAddressing
1115  (
1116  procMeshList,
1117  proci,
1118  "pointProcAddressing",
1119  pointProcAddressingList
1120  );
1121 
1122  const pointMesh& procPMesh = pointMesh::New(procMesh);
1123 
1124  if (!pointFieldDecomposerList.set(proci))
1125  {
1126  pointFieldDecomposerList.set
1127  (
1128  proci,
1130  (
1131  pMesh,
1132  procPMesh,
1133  pointProcAddressing,
1134  boundaryProcAddressing
1135  )
1136  );
1137  }
1138 
1139  pointFieldCache.decomposeAllFields
1140  (
1141  pointFieldDecomposerList[proci]
1142  );
1143 
1144  if (times.size() == 1)
1145  {
1146  pointProcAddressingList.set(proci, nullptr);
1147  pointFieldDecomposerList.set(proci, nullptr);
1148  }
1149  }
1150 
1151 
1152  // If there is lagrangian data write it out
1153  forAll(lagrangianPositions, cloudi)
1154  {
1155  if (lagrangianPositions[cloudi].size())
1156  {
1157  lagrangianFieldDecomposer fieldDecomposer
1158  (
1159  mesh,
1160  procMesh,
1161  faceProcAddressing,
1162  cellProcAddressing,
1163  cloudDirs[cloudi],
1164  lagrangianPositions[cloudi],
1165  cellParticles[cloudi]
1166  );
1167 
1168  // Lagrangian fields
1169  lagrangianFieldCache.decomposeAllFields
1170  (
1171  cloudi,
1172  cloudDirs[cloudi],
1173  fieldDecomposer
1174  );
1175  }
1176  }
1177 
1178  if (doDecompFields)
1179  {
1180  // Decompose "uniform" directory in the time region
1181  // directory
1182  decomposeUniform
1183  (
1184  copyUniform, mesh, processorDb, regionDir
1185  );
1186 
1187  // For a multi-region case, also decompose "uniform"
1188  // directory in the time directory
1189  if (regionNames.size() > 1 && regioni == 0)
1190  {
1191  decomposeUniform(copyUniform, mesh, processorDb);
1192  }
1193  }
1194 
1195 
1196  // We have cached all the constant mesh data for the current
1197  // processor. This is only important if running with
1198  // multiple times, otherwise it is just extra storage.
1199  if (times.size() == 1)
1200  {
1201  boundaryProcAddressingList.set(proci, nullptr);
1202  cellProcAddressingList.set(proci, nullptr);
1203  faceProcAddressingList.set(proci, nullptr);
1204  procMeshList.set(proci, nullptr);
1205  processorDbList.set(proci, nullptr);
1206  }
1207  }
1208 
1209 
1210  // Finite area mesh and field decomposition
1211  if (faMeshDecompPtr)
1212  {
1213  Info<< "\nFinite area mesh decomposition" << endl;
1214 
1215  faMeshDecomposition& aMesh = faMeshDecompPtr();
1216 
1217  aMesh.decomposeMesh();
1218  aMesh.writeDecomposition();
1219 
1220 
1221  // Area/edge fields
1222  // ~~~~~~~~~~~~~~~~
1223 
1224  faFieldDecomposer::fieldsCache areaFieldCache;
1225 
1226  if (doDecompFields)
1227  {
1228  areaFieldCache.readAllFields(aMesh, objects);
1229  }
1230 
1231  const label nAreaFields = areaFieldCache.size();
1232 
1233  Info<< endl;
1234  Info<< "Finite area field transfer: "
1235  << nAreaFields << " fields" << endl;
1236 
1237  // Split the fields over processors
1238  for
1239  (
1240  label proci = 0;
1241  nAreaFields && proci < mesh.nProcs();
1242  ++proci
1243  )
1244  {
1245  Info<< " Processor " << proci << endl;
1246 
1247  // open the database
1248  Time processorDb
1249  (
1251  args.rootPath(),
1252  args.caseName()/("processor" + Foam::name(proci)),
1253  false, // No function objects
1254  false // No extra controlDict libs
1255  );
1256 
1257  processorDb.setTime(runTime);
1258 
1259  // Read the mesh
1260  fvMesh procFvMesh
1261  (
1262  IOobject
1263  (
1264  regionName,
1265  processorDb.timeName(),
1266  processorDb
1267  )
1268  );
1269 
1270  faMesh procMesh(procFvMesh);
1271 
1272  // // Does not work. HJ, 15/Aug/2017
1273  // const labelIOList& faceProcAddressing =
1274  // procAddressing
1275  // (
1276  // procMeshList,
1277  // proci,
1278  // "faceProcAddressing",
1279  // faceProcAddressingList
1280  // );
1281 
1282  // const labelIOList& boundaryProcAddressing =
1283  // procAddressing
1284  // (
1285  // procMeshList,
1286  // proci,
1287  // "boundaryProcAddressing",
1288  // boundaryProcAddressingList
1289  // );
1290 
1291  // Addressing from faMesh (not polyMesh) meshSubDir
1292 
1293  autoPtr<labelIOList> tfaceProcAddr =
1294  faProcAddressing
1295  (
1296  procFvMesh,
1297  "faceProcAddressing",
1298  runTime.constant()
1299  );
1300  auto& faceProcAddressing = *tfaceProcAddr;
1301 
1302  autoPtr<labelIOList> tboundaryProcAddr =
1303  faProcAddressing
1304  (
1305  procFvMesh,
1306  "boundaryProcAddressing",
1307  runTime.constant()
1308  );
1309  auto& boundaryProcAddressing = *tboundaryProcAddr;
1310 
1311  autoPtr<labelIOList> tedgeProcAddr =
1312  faProcAddressing
1313  (
1314  procFvMesh,
1315  "edgeProcAddressing",
1316  runTime.constant()
1317  );
1318  const auto& edgeProcAddressing = *tedgeProcAddr;
1319 
1320  faFieldDecomposer fieldDecomposer
1321  (
1322  aMesh,
1323  procMesh,
1324  edgeProcAddressing,
1325  faceProcAddressing,
1326  boundaryProcAddressing
1327  );
1328 
1329  areaFieldCache.decomposeAllFields(fieldDecomposer);
1330  }
1331  }
1332  }
1333  }
1334  }
1335 
1336  Info<< "\nEnd\n" << endl;
1337 
1338  return 0;
1339 }
1340 
1341 
1342 // ************************************************************************* //
Finite area mesh (used for 2-D non-Euclidian finite area method) defined using a patch of faces on a ...
Definition: faMesh.H:87
static void addNote(const string &note)
Add extra notes for the usage information.
Definition: argList.C:462
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
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:859
void append(const T &elem)
Add copy at back of list.
Definition: LList.H:744
word findInstance(const fileName &dir, const word &name=word::null, IOobjectOption::readOption rOpt=IOobjectOption::MUST_READ, const word &stopInstance=word::null) const
Return time instance (location) of dir that contains the file name (eg, used in reading mesh data)...
Definition: Time.C:725
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:847
static const pointMesh & New(const polyMesh &mesh, Args &&... args)
Get existing or create a new MeshObject. Registered with typeName.
Definition: MeshObject.C:53
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
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:410
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
bool chDir(const fileName &dir)
Change current directory to the one specified and return true on success.
Definition: POSIX.C:607
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:374
static void noParallel()
Remove the parallel options.
Definition: argList.C:584
static void addOptionCompat(const word &optName, std::pair< const char *, int > compat)
Specify an alias for the option name.
Definition: argList.C:418
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:2088
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:421
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:860
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:2106
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:614
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
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:504
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:405
bool writeDecomposition()
Write decomposition.
const word & system() const noexcept
Return system name.
Definition: TimePathsI.H:118
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:835
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:385
static instantList selectIfPresent(Time &runTime, const argList &args)
If any time option provided return the set of times (as select0) otherwise return just the current ti...
Definition: timeSelector.C:265
virtual void setTime(const Time &t)
Reset the time and time-index to those of the given time.
Definition: Time.C:935
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:705
static void addVerboseOption(const string &usage="", bool advanced=false)
Enable a &#39;verbose&#39; bool option, with usage information.
Definition: argList.C:520
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
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
The mesh sub-directory name (usually "faMesh")
Definition: faMesh.H:701
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:521
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:297
const word & regionDir
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers...
Definition: List.H:55
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:841
Nothing to be read.
Required Classes.
fileName cwd()
The physical or logical current working directory path name.
Definition: POSIX.C:590
Finite area boundary mesh.
label findCell(const point &p, const cellDecomposition=CELL_TETS) const
Find cell enclosing this location and return index.
Definition: polyMesh.C:1511
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:437
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)
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:963
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:172
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:590
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