reconstructPar.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2017 OpenFOAM Foundation
9  Copyright (C) 2015-2023 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  reconstructPar
29 
30 Group
31  grpParallelUtilities
32 
33 Description
34  Reconstructs fields of a case that is decomposed for parallel
35  execution of OpenFOAM.
36 
37 \*---------------------------------------------------------------------------*/
38 
39 #include "argList.H"
40 #include "timeSelector.H"
41 
42 #include "fvCFD.H"
43 #include "IOobjectList.H"
44 #include "processorMeshes.H"
45 #include "regionProperties.H"
46 #include "fvFieldReconstructor.H"
49 
50 #include "faCFD.H"
51 #include "faMesh.H"
52 #include "processorFaMeshes.H"
53 #include "faFieldReconstructor.H"
54 
55 #include "cellSet.H"
56 #include "faceSet.H"
57 #include "pointSet.H"
58 
59 #include "hexRef8Data.H"
60 
61 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
62 
63 bool haveAllTimes
64 (
65  const wordHashSet& masterTimeDirSet,
66  const instantList& timeDirs
67 )
68 {
69  // Loop over all times
70  for (const instant& t : timeDirs)
71  {
72  if (!masterTimeDirSet.found(t.name()))
73  {
74  return false;
75  }
76  }
77  return true;
78 }
79 
80 
81 int main(int argc, char *argv[])
82 {
83  argList::addNote
84  (
85  "Reconstruct fields of a parallel case"
86  );
87 
88  // Enable -constant ... if someone really wants it
89  // Enable -withZero to prevent accidentally trashing the initial fields
90  timeSelector::addOptions(true, true); // constant(true), zero(true)
91  argList::noParallel();
92 
93  #include "addAllRegionOptions.H"
94 
95  argList::addVerboseOption();
96  argList::addOption
97  (
98  "fields",
99  "wordRes",
100  "Specify single or multiple fields to reconstruct (all by default)."
101  " Eg, 'T' or '(p T U \"alpha.*\")'"
102  );
103  argList::addBoolOption
104  (
105  "no-fields", // noFields
106  "Skip reconstructing fields"
107  );
108  argList::addOptionCompat("no-fields", {"noFields", 2106});
109  argList::addOption
110  (
111  "lagrangianFields",
112  "wordRes",
113  "Specify single or multiple lagrangian fields to reconstruct"
114  " (all by default)."
115  " Eg, '(U d)'"
116  " - Positions are always included."
117  );
118  argList::addBoolOption
119  (
120  "no-lagrangian", // noLagrangian
121  "Skip reconstructing lagrangian positions and fields"
122  );
123  argList::addOptionCompat("no-lagrangian", {"noLagrangian", 2106});
124 
125  argList::addBoolOption
126  (
127  "no-sets",
128  "Skip reconstructing cellSets, faceSets, pointSets"
129  );
130  argList::addOptionCompat("no-sets", {"noSets", 2106});
131 
132  argList::addBoolOption
133  (
134  "newTimes",
135  "Only reconstruct new times (i.e. that do not exist already)"
136  );
137 
138  #include "setRootCase.H"
139  #include "createTime.H"
140 
141 
142  const bool doFields = !args.found("no-fields");
143  wordRes selectedFields;
144 
145  if (doFields)
146  {
147  args.readListIfPresent<wordRe>("fields", selectedFields);
148  }
149  else
150  {
151  Info<< "Skipping reconstructing fields";
152  if (args.found("fields"))
153  {
154  Info<< ". Ignore -fields option";
155  }
156  Info<< nl << endl;
157  }
158 
159 
160  const bool doFiniteArea = !args.found("no-finite-area");
161  if (!doFiniteArea)
162  {
163  Info<< "Skipping reconstructing finiteArea mesh/fields"
164  << nl << endl;
165  }
166 
167 
168  const bool doLagrangian = !args.found("no-lagrangian");
169  wordRes selectedLagrangianFields;
170 
171  if (doLagrangian)
172  {
173  args.readListIfPresent<wordRe>
174  (
175  "lagrangianFields", selectedLagrangianFields
176  );
177  }
178  else
179  {
180  Info<< "Skipping reconstructing lagrangian positions/fields";
181  if (args.found("lagrangianFields"))
182  {
183  Info<< ". Ignore -lagrangianFields option";
184  }
185  Info<< nl << endl;
186  }
187 
188 
189  const bool doReconstructSets = !args.found("no-sets");
190 
191  if (!doReconstructSets)
192  {
193  Info<< "Skipping reconstructing cellSets, faceSets and pointSets"
194  << nl << endl;
195  }
196 
197  const bool newTimes = args.found("newTimes");
198 
199  // Get region names
200  #include "getAllRegionOptions.H"
201 
202  // Determine the processor count
203  label nProcs{0};
204 
205  if (regionNames.empty())
206  {
208  << "No regions specified or detected."
209  << exit(FatalError);
210  }
211  else if (regionNames[0] == polyMesh::defaultRegion)
212  {
213  nProcs = fileHandler().nProcs(args.path());
214  }
215  else
216  {
217  nProcs = fileHandler().nProcs(args.path(), regionNames[0]);
218 
219  if (regionNames.size() == 1)
220  {
221  Info<< "Using region: " << regionNames[0] << nl << endl;
222  }
223  }
224 
225  if (!nProcs)
226  {
228  << "No processor* directories found"
229  << exit(FatalError);
230  }
231 
232  // Warn fileHandler of number of processors
233  const_cast<fileOperation&>(fileHandler()).nProcs(nProcs);
234 
235  // Create the processor databases
236  PtrList<Time> databases(nProcs);
237 
238  forAll(databases, proci)
239  {
240  databases.set
241  (
242  proci,
243  new Time
244  (
245  Time::controlDictName,
246  args.rootPath(),
247  args.caseName()/("processor" + Foam::name(proci)),
249  args.allowLibs()
250  )
251  );
252  }
253 
254  // Use the times list from the master processor
255  // and select a subset based on the command-line options
257  (
258  databases[0].times(),
259  args
260  );
261 
262  // Note that we do not set the runTime time so it is still the
263  // one set through the controlDict. The -time option
264  // only affects the selected set of times from processor0.
265  // - can be illogical
266  // + any point motion handled through mesh.readUpdate
267 
268 
269  if (timeDirs.empty())
270  {
271  WarningInFunction << "No times selected";
272  exit(1);
273  }
274 
275 
276  // Get current times if -newTimes
277  instantList masterTimeDirs;
278  if (newTimes)
279  {
280  masterTimeDirs = runTime.times();
281  }
282  wordHashSet masterTimeDirSet(2*masterTimeDirs.size());
283  for (const instant& t : masterTimeDirs)
284  {
285  masterTimeDirSet.insert(t.name());
286  }
287 
288 
289  // Set all times on processor meshes equal to reconstructed mesh
290  forAll(databases, proci)
291  {
292  databases[proci].setTime(runTime);
293  }
294 
295 
296  forAll(regionNames, regioni)
297  {
298  const word& regionName = regionNames[regioni];
300 
301  Info<< "\n\nReconstructing fields" << nl
302  << "region=" << regionName << nl << endl;
303 
304  if
305  (
306  newTimes
307  && regionNames.size() == 1
308  && regionDir.empty()
309  && haveAllTimes(masterTimeDirSet, timeDirs)
310  )
311  {
312  Info<< "Skipping region " << regionName
313  << " since already have all times"
314  << endl << endl;
315  continue;
316  }
317 
318 
319  fvMesh mesh
320  (
321  IOobject
322  (
323  regionName,
324  runTime.timeName(),
325  runTime,
327  )
328  );
329 
330 
331  // Read all meshes and addressing to reconstructed mesh
332  processorMeshes procMeshes(databases, regionName);
333 
334  // Loop over all times
335  forAll(timeDirs, timei)
336  {
337  if (newTimes && masterTimeDirSet.found(timeDirs[timei].name()))
338  {
339  Info<< "Skipping time " << timeDirs[timei].name()
340  << endl << endl;
341  continue;
342  }
343 
344 
345  // Set time for global database
346  runTime.setTime(timeDirs[timei], timei);
347 
348  Info<< "Time = " << runTime.timeName() << endl << endl;
349 
350  // Set time for all databases
351  forAll(databases, proci)
352  {
353  databases[proci].setTime(timeDirs[timei], timei);
354  }
355 
356  // Check if any new meshes need to be read.
357  polyMesh::readUpdateState meshStat = mesh.readUpdate();
358 
359  polyMesh::readUpdateState procStat = procMeshes.readUpdate();
360 
361  if (procStat == polyMesh::POINTS_MOVED)
362  {
363  // Reconstruct the points for moving mesh cases and write
364  // them out
365  procMeshes.reconstructPoints(mesh);
366  }
367  else if (meshStat != procStat)
368  {
370  << "readUpdate for the reconstructed mesh:"
371  << meshStat << nl
372  << "readUpdate for the processor meshes :"
373  << procStat << nl
374  << "These should be equal or your addressing"
375  << " might be incorrect."
376  << " Please check your time directories for any "
377  << "mesh directories." << endl;
378  }
379 
380 
381  // Get list of objects from processor0 database
382  IOobjectList objects
383  (
384  procMeshes.meshes()[0],
385  databases[0].timeName(),
386  IOobjectOption::NO_REGISTER
387  );
388 
389  IOobjectList faObjects;
390 
391  if (doFiniteArea && doFields)
392  {
393  // List of area mesh objects (assuming single region)
394  // - scan on processor0
395  faObjects = IOobjectList
396  (
397  procMeshes.meshes()[0],
398  databases[0].timeName(),
399  faMesh::dbDir(word::null), // local relative to mesh
400  IOobjectOption::NO_REGISTER
401  );
402  }
403 
404  if (doFields)
405  {
406  // If there are any FV fields, reconstruct them
407  Info<< "Reconstructing FV fields" << nl << endl;
408 
409  fvFieldReconstructor reconstructor
410  (
411  mesh,
412  procMeshes.meshes(),
413  procMeshes.faceProcAddressing(),
414  procMeshes.cellProcAddressing(),
415  procMeshes.boundaryProcAddressing()
416  );
417 
418  reconstructor.reconstructAllFields(objects, selectedFields);
419 
420  if (reconstructor.nReconstructed() == 0)
421  {
422  Info<< "No FV fields" << nl << endl;
423  }
424  }
425 
426  if (doFields)
427  {
428  Info<< "Reconstructing point fields" << nl << endl;
429 
430  const pointMesh& pMesh = pointMesh::New
431  (
432  mesh,
433  IOobject::READ_IF_PRESENT
434  );
435 
436  pointFieldReconstructor reconstructor
437  (
438  pMesh,
439  procMeshes.pointMeshes(),
440  procMeshes.pointProcAddressing(),
441  procMeshes.pointMeshBoundaryProcAddressing()
442  );
443 
444  reconstructor.reconstructAllFields(objects, selectedFields);
445 
446  if (reconstructor.nReconstructed() == 0)
447  {
448  Info<< "No point fields" << nl << endl;
449  }
450  }
451 
452 
453  // If there are any clouds, reconstruct them.
454  // The problem is that a cloud of size zero will not get written so
455  // in pass 1 we determine the cloud names and per cloud name the
456  // fields. Note that the fields are stored as IOobjectList from
457  // the first processor that has them. They are in pass2 only used
458  // for name and type (scalar, vector etc).
459 
460  if (doLagrangian)
461  {
462  HashTable<IOobjectList> allCloudObjects;
463 
464  forAll(databases, proci)
465  {
466  fileName lagrangianDir
467  (
468  fileHandler().filePath
469  (
470  databases[proci].timePath()
471  / regionDir
472  / cloud::prefix
473  )
474  );
475 
476  fileNameList cloudDirs;
477  if (!lagrangianDir.empty())
478  {
479  cloudDirs = fileHandler().readDir
480  (
481  lagrangianDir,
482  fileName::DIRECTORY
483  );
484  }
485 
486  for (const fileName& cloudDir : cloudDirs)
487  {
488  // Check if we already have cloud objects for this
489  // cloudname
490  if (!allCloudObjects.found(cloudDir))
491  {
492  // Do local scan for valid cloud objects
493  IOobjectList localObjs
494  (
495  procMeshes.meshes()[proci],
496  databases[proci].timeName(),
497  cloud::prefix/cloudDir
498  );
499 
500  if
501  (
502  localObjs.found("coordinates")
503  || localObjs.found("positions")
504  )
505  {
506  allCloudObjects.insert(cloudDir, localObjs);
507  }
508  }
509  }
510  }
511 
512 
513  if (allCloudObjects.size())
514  {
515  lagrangianReconstructor reconstructor
516  (
517  mesh,
518  procMeshes.meshes(),
519  procMeshes.faceProcAddressing(),
520  procMeshes.cellProcAddressing()
521  );
522 
523  // Pass2: reconstruct the cloud
524  forAllConstIters(allCloudObjects, iter)
525  {
526  const word cloudName = word::validate(iter.key());
527 
528  // Objects (on arbitrary processor)
529  const IOobjectList& cloudObjs = iter.val();
530 
531  Info<< "Reconstructing lagrangian fields for cloud "
532  << cloudName << nl << endl;
533 
534  reconstructor.reconstructPositions(cloudName);
535 
536  reconstructor.reconstructAllFields
537  (
538  cloudName,
539  cloudObjs,
540  selectedLagrangianFields
541  );
542  }
543  }
544  else
545  {
546  Info<< "No lagrangian fields" << nl << endl;
547  }
548  }
549 
550 
551  // If there are any FA fields, reconstruct them
552 
553  if (!doFiniteArea)
554  {
555  }
556  else if
557  (
558  faObjects.count<areaScalarField>()
559  || faObjects.count<areaVectorField>()
560  || faObjects.count<areaSphericalTensorField>()
561  || faObjects.count<areaSymmTensorField>()
562  || faObjects.count<areaTensorField>()
563  || faObjects.count<edgeScalarField>()
564  )
565  {
566  Info << "Reconstructing FA fields" << nl << endl;
567 
568  faMesh aMesh(mesh);
569 
570  processorFaMeshes procFaMeshes(procMeshes.meshes());
571 
572  faFieldReconstructor reconstructor
573  (
574  aMesh,
575  procFaMeshes.meshes(),
576  procFaMeshes.edgeProcAddressing(),
577  procFaMeshes.faceProcAddressing(),
578  procFaMeshes.boundaryProcAddressing()
579  );
580 
581  reconstructor.reconstructAllFields(faObjects);
582  }
583  else
584  {
585  Info << "No FA fields" << nl << endl;
586  }
587 
588  if (doReconstructSets)
589  {
590  // Scan to find all sets
591  HashTable<label> cSetNames;
592  HashTable<label> fSetNames;
593  HashTable<label> pSetNames;
594 
595  forAll(procMeshes.meshes(), proci)
596  {
597  const fvMesh& procMesh = procMeshes.meshes()[proci];
598 
599  // Note: look at sets in current time only or between
600  // mesh and current time?. For now current time. This will
601  // miss out on sets in intermediate times that have not
602  // been reconstructed.
603  IOobjectList objects
604  (
605  procMesh,
606  databases[0].timeName(), //procMesh.facesInstance()
607  polyMesh::meshSubDir/"sets"
608  );
609 
610  for (const IOobject& io : objects.csorted<cellSet>())
611  {
612  cSetNames.insert(io.name(), cSetNames.size());
613  }
614 
615  for (const IOobject& io : objects.csorted<faceSet>())
616  {
617  fSetNames.insert(io.name(), fSetNames.size());
618  }
619 
620  for (const IOobject& io : objects.csorted<pointSet>())
621  {
622  pSetNames.insert(io.name(), pSetNames.size());
623  }
624  }
625 
626  if (cSetNames.size() || fSetNames.size() || pSetNames.size())
627  {
628  // Construct all sets
629  PtrList<cellSet> cellSets(cSetNames.size());
630  PtrList<faceSet> faceSets(fSetNames.size());
631  PtrList<pointSet> pointSets(pSetNames.size());
632 
633  Info<< "Reconstructing sets:" << endl;
634  if (cSetNames.size())
635  {
636  Info<< " cellSets "
637  << cSetNames.sortedToc() << endl;
638  }
639  if (fSetNames.size())
640  {
641  Info<< " faceSets "
642  << fSetNames.sortedToc() << endl;
643  }
644  if (pSetNames.size())
645  {
646  Info<< " pointSets "
647  << pSetNames.sortedToc() << endl;
648  }
649 
650  // Load sets
651  forAll(procMeshes.meshes(), proci)
652  {
653  const fvMesh& procMesh = procMeshes.meshes()[proci];
654 
655  IOobjectList objects
656  (
657  procMesh,
658  databases[0].timeName(),
659  polyMesh::meshSubDir/"sets"
660  );
661 
662  // cellSets
663  const labelList& cellMap =
664  procMeshes.cellProcAddressing()[proci];
665 
666  for (const IOobject& io : objects.csorted<cellSet>())
667  {
668  // Load cellSet
669  const cellSet procSet(io);
670  const label seti = cSetNames[io.name()];
671  if (!cellSets.set(seti))
672  {
673  cellSets.set
674  (
675  seti,
676  new cellSet
677  (
678  mesh,
679  io.name(),
680  procSet.size()
681  )
682  );
683  }
684  cellSet& cSet = cellSets[seti];
685  cSet.instance() = runTime.timeName();
686 
687  for (const label celli : procSet)
688  {
689  cSet.insert(cellMap[celli]);
690  }
691  }
692 
693  // faceSets
694  const labelList& faceMap =
695  procMeshes.faceProcAddressing()[proci];
696 
697  for (const IOobject& io : objects.csorted<faceSet>())
698  {
699  // Load faceSet
700  const faceSet procSet(io);
701  const label seti = fSetNames[io.name()];
702  if (!faceSets.set(seti))
703  {
704  faceSets.set
705  (
706  seti,
707  new faceSet
708  (
709  mesh,
710  io.name(),
711  procSet.size()
712  )
713  );
714  }
715  faceSet& fSet = faceSets[seti];
716  fSet.instance() = runTime.timeName();
717 
718  for (const label facei : procSet)
719  {
720  fSet.insert(mag(faceMap[facei])-1);
721  }
722  }
723  // pointSets
724  const labelList& pointMap =
725  procMeshes.pointProcAddressing()[proci];
726 
727  for (const IOobject& io : objects.csorted<pointSet>())
728  {
729  // Load pointSet
730  const pointSet procSet(io);
731  const label seti = pSetNames[io.name()];
732  if (!pointSets.set(seti))
733  {
734  pointSets.set
735  (
736  seti,
737  new pointSet
738  (
739  mesh,
740  io.name(),
741  procSet.size()
742  )
743  );
744  }
745  pointSet& pSet = pointSets[seti];
746  pSet.instance() = runTime.timeName();
747 
748  for (const label pointi : procSet)
749  {
750  pSet.insert(pointMap[pointi]);
751  }
752  }
753  }
754 
755  // Write sets
756 
757  for (const auto& set : cellSets)
758  {
759  set.write();
760  }
761  for (const auto& set : faceSets)
762  {
763  set.write();
764  }
765  for (const auto& set : pointSets)
766  {
767  set.write();
768  }
769  }
770 
771 
772  // Reconstruct refinement data
773  {
774  PtrList<hexRef8Data> procData(procMeshes.meshes().size());
775 
776  forAll(procMeshes.meshes(), procI)
777  {
778  const fvMesh& procMesh = procMeshes.meshes()[procI];
779 
780  procData.set
781  (
782  procI,
783  new hexRef8Data
784  (
785  IOobject
786  (
787  "dummy",
788  procMesh.time().timeName(),
789  polyMesh::meshSubDir,
790  procMesh,
791  IOobject::READ_IF_PRESENT,
792  IOobject::NO_WRITE,
793  IOobject::NO_REGISTER
794  )
795  )
796  );
797  }
798 
799  // Combine individual parts
800 
801  const PtrList<labelIOList>& cellAddr =
802  procMeshes.cellProcAddressing();
803 
804  UPtrList<const labelList> cellMaps(cellAddr.size());
805  forAll(cellAddr, i)
806  {
807  cellMaps.set(i, &cellAddr[i]);
808  }
809 
810  const PtrList<labelIOList>& pointAddr =
811  procMeshes.pointProcAddressing();
812 
813  UPtrList<const labelList> pointMaps(pointAddr.size());
814  forAll(pointAddr, i)
815  {
816  pointMaps.set(i, &pointAddr[i]);
817  }
818 
819  UPtrList<const hexRef8Data> procRefs(procData.size());
820  forAll(procData, i)
821  {
822  procRefs.set(i, &procData[i]);
823  }
824 
825  hexRef8Data
826  (
827  IOobject
828  (
829  "dummy",
830  mesh.time().timeName(),
831  polyMesh::meshSubDir,
832  mesh,
833  IOobject::NO_READ,
834  IOobject::NO_WRITE,
835  IOobject::NO_REGISTER
836  ),
837  cellMaps,
838  pointMaps,
839  procRefs
840  ).write();
841  }
842  }
843 
844 
845  // Reconstruct refinement data
846  {
847  PtrList<hexRef8Data> procData(procMeshes.meshes().size());
848 
849  forAll(procMeshes.meshes(), procI)
850  {
851  const fvMesh& procMesh = procMeshes.meshes()[procI];
852 
853  procData.set
854  (
855  procI,
856  new hexRef8Data
857  (
858  IOobject
859  (
860  "dummy",
861  procMesh.time().timeName(),
862  polyMesh::meshSubDir,
863  procMesh,
864  IOobject::READ_IF_PRESENT,
865  IOobject::NO_WRITE,
866  IOobject::NO_REGISTER
867  )
868  )
869  );
870  }
871 
872  // Combine individual parts
873 
874  const PtrList<labelIOList>& cellAddr =
875  procMeshes.cellProcAddressing();
876 
877  UPtrList<const labelList> cellMaps(cellAddr.size());
878  forAll(cellAddr, i)
879  {
880  cellMaps.set(i, &cellAddr[i]);
881  }
882 
883  const PtrList<labelIOList>& pointAddr =
884  procMeshes.pointProcAddressing();
885 
886  UPtrList<const labelList> pointMaps(pointAddr.size());
887  forAll(pointAddr, i)
888  {
889  pointMaps.set(i, &pointAddr[i]);
890  }
891 
892  UPtrList<const hexRef8Data> procRefs(procData.size());
893  forAll(procData, i)
894  {
895  procRefs.set(i, &procData[i]);
896  }
897 
898  hexRef8Data
899  (
900  IOobject
901  (
902  "dummy",
903  mesh.time().timeName(),
904  polyMesh::meshSubDir,
905  mesh,
906  IOobject::NO_READ,
907  IOobject::NO_WRITE,
908  IOobject::NO_REGISTER
909  ),
910  cellMaps,
911  pointMaps,
912  procRefs
913  ).write();
914  }
915 
916  // If there is a "uniform" directory in the time region
917  // directory copy from the master processor
918  {
919  fileName uniformDir0
920  (
921  fileHandler().filePath
922  (
923  databases[0].timePath()/regionDir/"uniform"
924  )
925  );
926 
927  if (!uniformDir0.empty() && fileHandler().isDir(uniformDir0))
928  {
929  fileHandler().cp(uniformDir0, runTime.timePath()/regionDir);
930  }
931  }
932 
933  // For the first region of a multi-region case additionally
934  // copy the "uniform" directory in the time directory
935  if (regioni == 0 && !regionDir.empty())
936  {
937  fileName uniformDir0
938  (
939  fileHandler().filePath
940  (
941  databases[0].timePath()/"uniform"
942  )
943  );
944 
945  if (!uniformDir0.empty() && fileHandler().isDir(uniformDir0))
946  {
947  fileHandler().cp(uniformDir0, runTime.timePath());
948  }
949  }
950  }
951  }
952 
953  Info<< "\nEnd\n" << endl;
954 
955  return 0;
956 }
957 
958 
959 // ************************************************************************* //
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:600
List< bool > select(const label n, const labelUList &locations)
Construct a selection list of bools (all false) with the given pre-size, subsequently add specified l...
Definition: BitOps.C:134
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:529
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
refPtr< fileOperation > fileHandler(std::nullptr_t)
Delete current file handler - forwards to fileOperation::handler()
GeometricField< symmTensor, faPatchField, areaMesh > areaSymmTensorField
Definition: areaFieldsFwd.H:84
wordList regionNames
bool insert(const Key &key)
Insert a new entry, not overwriting existing entries.
Definition: HashSet.H:229
bool allowFunctionObjects() const
The controlDict &#39;functions&#39; entry is allowed to be used.
Definition: argList.C:2252
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
StringType validate(const std::string &str, const UnaryPredicate &accept, const bool invert=false)
Return a copy of the input string with validated characters.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:286
word timeName
Definition: getTimeIndex.H:3
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
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
const word cloudName(propsDict.get< word >("cloud"))
HashSet< word, Hash< word > > wordHashSet
A HashSet of words, uses string hasher.
Definition: HashSet.H:73
GeometricField< scalar, faePatchField, edgeMesh > edgeScalarField
Definition: edgeFieldsFwd.H:43
GeometricField< tensor, faPatchField, areaMesh > areaTensorField
Definition: areaFieldsFwd.H:88
const fileName & caseName() const noexcept
Return case name (parallel run) or global case (serial run)
Definition: argListI.H:62
GeometricField< sphericalTensor, faPatchField, areaMesh > areaSphericalTensorField
Definition: areaFieldsFwd.H:80
const fileName & rootPath() const noexcept
Return root path.
Definition: argListI.H:56
bool readListIfPresent(const word &optName, List< T > &list) const
If named option is present, get a List of values treating a single entry like a list of size 1...
Definition: argListI.H:387
fileName path() const
Return the full path to the (processor local) case.
Definition: argListI.H:74
#define WarningInFunction
Report a warning using Foam::Warning.
Foam::word regionName(args.getOrDefault< word >("region", Foam::polyMesh::defaultRegion))
const word & regionDir
Required Classes.
messageStream Info
Information stream (stdout output on master, null elsewhere)
List< label > labelList
A List of labels.
Definition: List.H:61
IOobject io("surfaceFilmProperties", mesh.time().constant(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE, IOobject::NO_REGISTER)
List< instant > instantList
List of instants.
Definition: instantList.H:41
List< fileName > fileNameList
List of fileName.
Definition: fileNameList.H:32
Foam::argList args(argc, argv)
GeometricField< vector, faPatchField, areaMesh > areaVectorField
Definition: areaFieldsFwd.H:76
bool found(const word &optName) const
Return true if the named option is found.
Definition: argListI.H:171
forAllConstIters(mixture.phases(), phase)
Definition: pEqn.H:28
GeometricField< scalar, faPatchField, areaMesh > areaScalarField
Definition: areaFieldsFwd.H:72