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  );
387 
388  if (doFields)
389  {
390  // If there are any FV fields, reconstruct them
391  Info<< "Reconstructing FV fields" << nl << endl;
392 
393  fvFieldReconstructor reconstructor
394  (
395  mesh,
396  procMeshes.meshes(),
397  procMeshes.faceProcAddressing(),
398  procMeshes.cellProcAddressing(),
399  procMeshes.boundaryProcAddressing()
400  );
401 
402  reconstructor.reconstructAllFields(objects, selectedFields);
403 
404  if (reconstructor.nReconstructed() == 0)
405  {
406  Info<< "No FV fields" << nl << endl;
407  }
408  }
409 
410  if (doFields)
411  {
412  Info<< "Reconstructing point fields" << nl << endl;
413 
414  const pointMesh& pMesh = pointMesh::New(mesh);
415  PtrList<pointMesh> pMeshes(procMeshes.meshes().size());
416 
417  forAll(pMeshes, proci)
418  {
419  pMeshes.set
420  (
421  proci,
422  new pointMesh(procMeshes.meshes()[proci])
423  );
424  }
425 
426  pointFieldReconstructor reconstructor
427  (
428  pMesh,
429  pMeshes,
430  procMeshes.pointProcAddressing(),
431  procMeshes.boundaryProcAddressing()
432  );
433 
434  reconstructor.reconstructAllFields(objects, selectedFields);
435 
436  if (reconstructor.nReconstructed() == 0)
437  {
438  Info<< "No point fields" << nl << endl;
439  }
440  }
441 
442 
443  // If there are any clouds, reconstruct them.
444  // The problem is that a cloud of size zero will not get written so
445  // in pass 1 we determine the cloud names and per cloud name the
446  // fields. Note that the fields are stored as IOobjectList from
447  // the first processor that has them. They are in pass2 only used
448  // for name and type (scalar, vector etc).
449 
450  if (doLagrangian)
451  {
452  HashTable<IOobjectList> allCloudObjects;
453 
454  forAll(databases, proci)
455  {
456  fileName lagrangianDir
457  (
458  fileHandler().filePath
459  (
460  databases[proci].timePath()
461  / regionDir
462  / cloud::prefix
463  )
464  );
465 
466  fileNameList cloudDirs;
467  if (!lagrangianDir.empty())
468  {
469  cloudDirs = fileHandler().readDir
470  (
471  lagrangianDir,
472  fileName::DIRECTORY
473  );
474  }
475 
476  for (const fileName& cloudDir : cloudDirs)
477  {
478  // Check if we already have cloud objects for this
479  // cloudname
480  if (!allCloudObjects.found(cloudDir))
481  {
482  // Do local scan for valid cloud objects
483  IOobjectList localObjs
484  (
485  procMeshes.meshes()[proci],
486  databases[proci].timeName(),
487  cloud::prefix/cloudDir
488  );
489 
490  if
491  (
492  localObjs.found("coordinates")
493  || localObjs.found("positions")
494  )
495  {
496  allCloudObjects.insert(cloudDir, localObjs);
497  }
498  }
499  }
500  }
501 
502 
503  if (allCloudObjects.size())
504  {
505  lagrangianReconstructor reconstructor
506  (
507  mesh,
508  procMeshes.meshes(),
509  procMeshes.faceProcAddressing(),
510  procMeshes.cellProcAddressing()
511  );
512 
513  // Pass2: reconstruct the cloud
514  forAllConstIters(allCloudObjects, iter)
515  {
516  const word cloudName = word::validate(iter.key());
517 
518  // Objects (on arbitrary processor)
519  const IOobjectList& cloudObjs = iter.val();
520 
521  Info<< "Reconstructing lagrangian fields for cloud "
522  << cloudName << nl << endl;
523 
524  reconstructor.reconstructPositions(cloudName);
525 
526  reconstructor.reconstructAllFields
527  (
528  cloudName,
529  cloudObjs,
530  selectedLagrangianFields
531  );
532  }
533  }
534  else
535  {
536  Info<< "No lagrangian fields" << nl << endl;
537  }
538  }
539 
540 
541  // If there are any FA fields, reconstruct them
542 
543  if (!doFiniteArea)
544  {
545  }
546  else if
547  (
548  objects.count<areaScalarField>()
549  || objects.count<areaVectorField>()
550  || objects.count<areaSphericalTensorField>()
551  || objects.count<areaSymmTensorField>()
552  || objects.count<areaTensorField>()
553  || objects.count<edgeScalarField>()
554  )
555  {
556  Info << "Reconstructing FA fields" << nl << endl;
557 
558  faMesh aMesh(mesh);
559 
560  processorFaMeshes procFaMeshes(procMeshes.meshes());
561 
562  faFieldReconstructor reconstructor
563  (
564  aMesh,
565  procFaMeshes.meshes(),
566  procFaMeshes.edgeProcAddressing(),
567  procFaMeshes.faceProcAddressing(),
568  procFaMeshes.boundaryProcAddressing()
569  );
570 
571  reconstructor.reconstructAllFields(objects);
572  }
573  else
574  {
575  Info << "No FA fields" << nl << endl;
576  }
577 
578  if (doReconstructSets)
579  {
580  // Scan to find all sets
581  HashTable<label> cSetNames;
582  HashTable<label> fSetNames;
583  HashTable<label> pSetNames;
584 
585  forAll(procMeshes.meshes(), proci)
586  {
587  const fvMesh& procMesh = procMeshes.meshes()[proci];
588 
589  // Note: look at sets in current time only or between
590  // mesh and current time?. For now current time. This will
591  // miss out on sets in intermediate times that have not
592  // been reconstructed.
593  IOobjectList objects
594  (
595  procMesh,
596  databases[0].timeName(), //procMesh.facesInstance()
597  polyMesh::meshSubDir/"sets"
598  );
599 
600  for (const IOobject& io : objects.csorted<cellSet>())
601  {
602  cSetNames.insert(io.name(), cSetNames.size());
603  }
604 
605  for (const IOobject& io : objects.csorted<faceSet>())
606  {
607  fSetNames.insert(io.name(), fSetNames.size());
608  }
609 
610  for (const IOobject& io : objects.csorted<pointSet>())
611  {
612  pSetNames.insert(io.name(), pSetNames.size());
613  }
614  }
615 
616  if (cSetNames.size() || fSetNames.size() || pSetNames.size())
617  {
618  // Construct all sets
619  PtrList<cellSet> cellSets(cSetNames.size());
620  PtrList<faceSet> faceSets(fSetNames.size());
621  PtrList<pointSet> pointSets(pSetNames.size());
622 
623  Info<< "Reconstructing sets:" << endl;
624  if (cSetNames.size())
625  {
626  Info<< " cellSets "
627  << cSetNames.sortedToc() << endl;
628  }
629  if (fSetNames.size())
630  {
631  Info<< " faceSets "
632  << fSetNames.sortedToc() << endl;
633  }
634  if (pSetNames.size())
635  {
636  Info<< " pointSets "
637  << pSetNames.sortedToc() << endl;
638  }
639 
640  // Load sets
641  forAll(procMeshes.meshes(), proci)
642  {
643  const fvMesh& procMesh = procMeshes.meshes()[proci];
644 
645  IOobjectList objects
646  (
647  procMesh,
648  databases[0].timeName(),
649  polyMesh::meshSubDir/"sets"
650  );
651 
652  // cellSets
653  const labelList& cellMap =
654  procMeshes.cellProcAddressing()[proci];
655 
656  for (const IOobject& io : objects.csorted<cellSet>())
657  {
658  // Load cellSet
659  const cellSet procSet(io);
660  const label seti = cSetNames[io.name()];
661  if (!cellSets.set(seti))
662  {
663  cellSets.set
664  (
665  seti,
666  new cellSet
667  (
668  mesh,
669  io.name(),
670  procSet.size()
671  )
672  );
673  }
674  cellSet& cSet = cellSets[seti];
675  cSet.instance() = runTime.timeName();
676 
677  for (const label celli : procSet)
678  {
679  cSet.insert(cellMap[celli]);
680  }
681  }
682 
683  // faceSets
684  const labelList& faceMap =
685  procMeshes.faceProcAddressing()[proci];
686 
687  for (const IOobject& io : objects.csorted<faceSet>())
688  {
689  // Load faceSet
690  const faceSet procSet(io);
691  const label seti = fSetNames[io.name()];
692  if (!faceSets.set(seti))
693  {
694  faceSets.set
695  (
696  seti,
697  new faceSet
698  (
699  mesh,
700  io.name(),
701  procSet.size()
702  )
703  );
704  }
705  faceSet& fSet = faceSets[seti];
706  fSet.instance() = runTime.timeName();
707 
708  for (const label facei : procSet)
709  {
710  fSet.insert(mag(faceMap[facei])-1);
711  }
712  }
713  // pointSets
714  const labelList& pointMap =
715  procMeshes.pointProcAddressing()[proci];
716 
717  for (const IOobject& io : objects.csorted<pointSet>())
718  {
719  // Load pointSet
720  const pointSet procSet(io);
721  const label seti = pSetNames[io.name()];
722  if (!pointSets.set(seti))
723  {
724  pointSets.set
725  (
726  seti,
727  new pointSet
728  (
729  mesh,
730  io.name(),
731  procSet.size()
732  )
733  );
734  }
735  pointSet& pSet = pointSets[seti];
736  pSet.instance() = runTime.timeName();
737 
738  for (const label pointi : procSet)
739  {
740  pSet.insert(pointMap[pointi]);
741  }
742  }
743  }
744 
745  // Write sets
746 
747  for (const auto& set : cellSets)
748  {
749  set.write();
750  }
751  for (const auto& set : faceSets)
752  {
753  set.write();
754  }
755  for (const auto& set : pointSets)
756  {
757  set.write();
758  }
759  }
760 
761 
762  // Reconstruct refinement data
763  {
764  PtrList<hexRef8Data> procData(procMeshes.meshes().size());
765 
766  forAll(procMeshes.meshes(), procI)
767  {
768  const fvMesh& procMesh = procMeshes.meshes()[procI];
769 
770  procData.set
771  (
772  procI,
773  new hexRef8Data
774  (
775  IOobject
776  (
777  "dummy",
778  procMesh.time().timeName(),
779  polyMesh::meshSubDir,
780  procMesh,
781  IOobject::READ_IF_PRESENT,
782  IOobject::NO_WRITE,
783  IOobject::NO_REGISTER
784  )
785  )
786  );
787  }
788 
789  // Combine individual parts
790 
791  const PtrList<labelIOList>& cellAddr =
792  procMeshes.cellProcAddressing();
793 
794  UPtrList<const labelList> cellMaps(cellAddr.size());
795  forAll(cellAddr, i)
796  {
797  cellMaps.set(i, &cellAddr[i]);
798  }
799 
800  const PtrList<labelIOList>& pointAddr =
801  procMeshes.pointProcAddressing();
802 
803  UPtrList<const labelList> pointMaps(pointAddr.size());
804  forAll(pointAddr, i)
805  {
806  pointMaps.set(i, &pointAddr[i]);
807  }
808 
809  UPtrList<const hexRef8Data> procRefs(procData.size());
810  forAll(procData, i)
811  {
812  procRefs.set(i, &procData[i]);
813  }
814 
815  hexRef8Data
816  (
817  IOobject
818  (
819  "dummy",
820  mesh.time().timeName(),
821  polyMesh::meshSubDir,
822  mesh,
823  IOobject::NO_READ,
824  IOobject::NO_WRITE,
825  IOobject::NO_REGISTER
826  ),
827  cellMaps,
828  pointMaps,
829  procRefs
830  ).write();
831  }
832  }
833 
834 
835  // Reconstruct refinement data
836  {
837  PtrList<hexRef8Data> procData(procMeshes.meshes().size());
838 
839  forAll(procMeshes.meshes(), procI)
840  {
841  const fvMesh& procMesh = procMeshes.meshes()[procI];
842 
843  procData.set
844  (
845  procI,
846  new hexRef8Data
847  (
848  IOobject
849  (
850  "dummy",
851  procMesh.time().timeName(),
852  polyMesh::meshSubDir,
853  procMesh,
854  IOobject::READ_IF_PRESENT,
855  IOobject::NO_WRITE,
856  IOobject::NO_REGISTER
857  )
858  )
859  );
860  }
861 
862  // Combine individual parts
863 
864  const PtrList<labelIOList>& cellAddr =
865  procMeshes.cellProcAddressing();
866 
867  UPtrList<const labelList> cellMaps(cellAddr.size());
868  forAll(cellAddr, i)
869  {
870  cellMaps.set(i, &cellAddr[i]);
871  }
872 
873  const PtrList<labelIOList>& pointAddr =
874  procMeshes.pointProcAddressing();
875 
876  UPtrList<const labelList> pointMaps(pointAddr.size());
877  forAll(pointAddr, i)
878  {
879  pointMaps.set(i, &pointAddr[i]);
880  }
881 
882  UPtrList<const hexRef8Data> procRefs(procData.size());
883  forAll(procData, i)
884  {
885  procRefs.set(i, &procData[i]);
886  }
887 
888  hexRef8Data
889  (
890  IOobject
891  (
892  "dummy",
893  mesh.time().timeName(),
894  polyMesh::meshSubDir,
895  mesh,
896  IOobject::NO_READ,
897  IOobject::NO_WRITE,
898  IOobject::NO_REGISTER
899  ),
900  cellMaps,
901  pointMaps,
902  procRefs
903  ).write();
904  }
905 
906  // If there is a "uniform" directory in the time region
907  // directory copy from the master processor
908  {
909  fileName uniformDir0
910  (
911  fileHandler().filePath
912  (
913  databases[0].timePath()/regionDir/"uniform"
914  )
915  );
916 
917  if (!uniformDir0.empty() && fileHandler().isDir(uniformDir0))
918  {
919  fileHandler().cp(uniformDir0, runTime.timePath()/regionDir);
920  }
921  }
922 
923  // For the first region of a multi-region case additionally
924  // copy the "uniform" directory in the time directory
925  if (regioni == 0 && !regionDir.empty())
926  {
927  fileName uniformDir0
928  (
929  fileHandler().filePath
930  (
931  databases[0].timePath()/"uniform"
932  )
933  );
934 
935  if (!uniformDir0.empty() && fileHandler().isDir(uniformDir0))
936  {
937  fileHandler().cp(uniformDir0, runTime.timePath());
938  }
939  }
940  }
941  }
942 
943  Info<< "\nEnd\n" << endl;
944 
945  return 0;
946 }
947 
948 
949 // ************************************************************************* //
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:598
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:531
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:232
bool allowFunctionObjects() const
The controlDict &#39;functions&#39; entry is allowed to be used.
Definition: argList.C:2088
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
StringType validate(const std::string &str, const UnaryPredicate &accept, const bool invert=false)
Return a copy of the input string with validated characters.
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:2106
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:46
GeometricField< tensor, faPatchField, areaMesh > areaTensorField
Definition: areaFieldsFwd.H:85
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:83
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:62
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:81
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:80