redistributePar.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 | Copyright (C) 2011-2016 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 Application
25  redistributePar
26 
27 Description
28  Redistributes existing decomposed mesh and fields according to the current
29  settings in the decomposeParDict file.
30 
31  Must be run on maximum number of source and destination processors.
32  Balances mesh and writes new mesh to new time directory.
33 
34  Can also work like decomposePar:
35  \verbatim
36  # Create empty processor directories (have to exist for argList)
37  mkdir processor0
38  ..
39  mkdir processorN
40 
41  # Copy undecomposed polyMesh
42  cp -r constant processor0
43 
44  # Distribute
45  mpirun -np ddd redistributePar -parallel
46  \endverbatim
47 \*---------------------------------------------------------------------------*/
48 
49 #include "fvMesh.H"
50 #include "decompositionMethod.H"
51 #include "PstreamReduceOps.H"
52 #include "fvCFD.H"
53 #include "fvMeshDistribute.H"
54 #include "mapDistributePolyMesh.H"
55 #include "IOobjectList.H"
56 #include "globalIndex.H"
57 #include "loadOrCreateMesh.H"
58 
59 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
60 
61 // Tolerance (as fraction of the bounding box). Needs to be fairly lax since
62 // usually meshes get written with limited precision (6 digits)
63 static const scalar defaultMergeTol = 1e-6;
64 
65 
66 // Get merging distance when matching face centres
67 scalar getMergeDistance
68 (
69  const argList& args,
70  const Time& runTime,
71  const boundBox& bb
72 )
73 {
74  scalar mergeTol = defaultMergeTol;
75  args.optionReadIfPresent("mergeTol", mergeTol);
76 
77  scalar writeTol =
78  Foam::pow(scalar(10.0), -scalar(IOstream::defaultPrecision()));
79 
80  Info<< "Merge tolerance : " << mergeTol << nl
81  << "Write tolerance : " << writeTol << endl;
82 
83  if (runTime.writeFormat() == IOstream::ASCII && mergeTol < writeTol)
84  {
86  << "Your current settings specify ASCII writing with "
87  << IOstream::defaultPrecision() << " digits precision." << endl
88  << "Your merging tolerance (" << mergeTol << ") is finer than this."
89  << endl
90  << "Please change your writeFormat to binary"
91  << " or increase the writePrecision" << endl
92  << "or adjust the merge tolerance (-mergeTol)."
93  << exit(FatalError);
94  }
95 
96  scalar mergeDist = mergeTol * bb.mag();
97 
98  Info<< "Overall meshes bounding box : " << bb << nl
99  << "Relative tolerance : " << mergeTol << nl
100  << "Absolute matching distance : " << mergeDist << nl
101  << endl;
102 
103  return mergeDist;
104 }
105 
106 
107 //void printMeshData(Ostream& os, const polyMesh& mesh)
108 //{
109 // os << "Number of points: " << mesh.points().size() << nl
110 // << " faces: " << mesh.faces().size() << nl
111 // << " internal faces: " << mesh.faceNeighbour().size() << nl
112 // << " cells: " << mesh.cells().size() << nl
113 // << " boundary patches: " << mesh.boundaryMesh().size() << nl
114 // << " point zones: " << mesh.pointZones().size() << nl
115 // << " face zones: " << mesh.faceZones().size() << nl
116 // << " cell zones: " << mesh.cellZones().size() << nl;
117 //}
118 
119 
120 void printMeshData(const polyMesh& mesh)
121 {
122  // Collect all data on master
123 
124  globalIndex globalCells(mesh.nCells());
125  labelListList patchNeiProcNo(Pstream::nProcs());
126  labelListList patchSize(Pstream::nProcs());
127  const labelList& pPatches = mesh.globalData().processorPatches();
128  patchNeiProcNo[Pstream::myProcNo()].setSize(pPatches.size());
129  patchSize[Pstream::myProcNo()].setSize(pPatches.size());
130  forAll(pPatches, i)
131  {
132  const processorPolyPatch& ppp = refCast<const processorPolyPatch>
133  (
134  mesh.boundaryMesh()[pPatches[i]]
135  );
136  patchNeiProcNo[Pstream::myProcNo()][i] = ppp.neighbProcNo();
137  patchSize[Pstream::myProcNo()][i] = ppp.size();
138  }
139  Pstream::gatherList(patchNeiProcNo);
140  Pstream::gatherList(patchSize);
141 
142 
143  // Print stats
144 
145  globalIndex globalBoundaryFaces(mesh.nFaces()-mesh.nInternalFaces());
146 
147  label maxProcCells = 0;
148  label totProcFaces = 0;
149  label maxProcPatches = 0;
150  label totProcPatches = 0;
151  label maxProcFaces = 0;
152 
153  for (label proci = 0; proci < Pstream::nProcs(); proci++)
154  {
155  Info<< endl
156  << "Processor " << proci << nl
157  << " Number of cells = " << globalCells.localSize(proci)
158  << endl;
159 
160  label nProcFaces = 0;
161 
162  const labelList& nei = patchNeiProcNo[proci];
163 
164  forAll(patchNeiProcNo[proci], i)
165  {
166  Info<< " Number of faces shared with processor "
167  << patchNeiProcNo[proci][i] << " = " << patchSize[proci][i]
168  << endl;
169 
170  nProcFaces += patchSize[proci][i];
171  }
172 
173  Info<< " Number of processor patches = " << nei.size() << nl
174  << " Number of processor faces = " << nProcFaces << nl
175  << " Number of boundary faces = "
176  << globalBoundaryFaces.localSize(proci) << endl;
177 
178  maxProcCells = max(maxProcCells, globalCells.localSize(proci));
179  totProcFaces += nProcFaces;
180  totProcPatches += nei.size();
181  maxProcPatches = max(maxProcPatches, nei.size());
182  maxProcFaces = max(maxProcFaces, nProcFaces);
183  }
184 
185  // Stats
186 
187  scalar avgProcCells = scalar(globalCells.size())/Pstream::nProcs();
188  scalar avgProcPatches = scalar(totProcPatches)/Pstream::nProcs();
189  scalar avgProcFaces = scalar(totProcFaces)/Pstream::nProcs();
190 
191  // In case of all faces on one processor. Just to avoid division by 0.
192  if (totProcPatches == 0)
193  {
194  avgProcPatches = 1;
195  }
196  if (totProcFaces == 0)
197  {
198  avgProcFaces = 1;
199  }
200 
201  Info<< nl
202  << "Number of processor faces = " << totProcFaces/2 << nl
203  << "Max number of cells = " << maxProcCells
204  << " (" << 100.0*(maxProcCells-avgProcCells)/avgProcCells
205  << "% above average " << avgProcCells << ")" << nl
206  << "Max number of processor patches = " << maxProcPatches
207  << " (" << 100.0*(maxProcPatches-avgProcPatches)/avgProcPatches
208  << "% above average " << avgProcPatches << ")" << nl
209  << "Max number of faces between processors = " << maxProcFaces
210  << " (" << 100.0*(maxProcFaces-avgProcFaces)/avgProcFaces
211  << "% above average " << avgProcFaces << ")" << nl
212  << endl;
213 }
214 
215 
216 // Debugging: write volScalarField with decomposition for post processing.
217 void writeDecomposition
218 (
219  const word& name,
220  const fvMesh& mesh,
221  const labelList& decomp
222 )
223 {
224  Info<< "Writing wanted cell distribution to volScalarField " << name
225  << " for postprocessing purposes." << nl << endl;
226 
227  volScalarField procCells
228  (
229  IOobject
230  (
231  name,
232  mesh.time().timeName(),
233  mesh,
234  IOobject::NO_READ,
235  IOobject::AUTO_WRITE,
236  false // do not register
237  ),
238  mesh,
239  dimensionedScalar(name, dimless, -1),
240  extrapolatedCalculatedFvPatchScalarField::typeName
241  );
242 
243  forAll(procCells, cI)
244  {
245  procCells[cI] = decomp[cI];
246  }
247  procCells.correctBoundaryConditions();
248 
249  procCells.write();
250 }
251 
252 
253 template<class GeoField>
254 void readFields
255 (
256  const boolList& haveMesh,
257  const fvMesh& mesh,
258  const autoPtr<fvMeshSubset>& subsetterPtr,
259  IOobjectList& allObjects,
260  PtrList<GeoField>& fields
261 )
262 {
263  // Get my objects of type
264  IOobjectList objects(allObjects.lookupClass(GeoField::typeName));
265 
266  // Check that we all have all objects
267  wordList objectNames = objects.toc();
268 
269  // Get master names
270  wordList masterNames(objectNames);
271  Pstream::scatter(masterNames);
272 
273  if (haveMesh[Pstream::myProcNo()] && objectNames != masterNames)
274  {
276  << "differing fields of type " << GeoField::typeName
277  << " on processors." << endl
278  << "Master has:" << masterNames << endl
279  << Pstream::myProcNo() << " has:" << objectNames
280  << abort(FatalError);
281  }
282 
283  fields.setSize(masterNames.size());
284 
285  // Have master send all fields to processors that don't have a mesh
286  if (Pstream::master())
287  {
288  forAll(masterNames, i)
289  {
290  const word& name = masterNames[i];
291  IOobject& io = *objects[name];
292  io.writeOpt() = IOobject::AUTO_WRITE;
293 
294  // Load field
295  fields.set(i, new GeoField(io, mesh));
296 
297  // Create zero sized field and send
298  if (subsetterPtr.valid())
299  {
300  tmp<GeoField> tsubfld = subsetterPtr().interpolate(fields[i]);
301 
302  // Send to all processors that don't have a mesh
303  for (label proci = 1; proci < Pstream::nProcs(); proci++)
304  {
305  if (!haveMesh[proci])
306  {
307  OPstream toProc(Pstream::blocking, proci);
308  toProc<< tsubfld();
309  }
310  }
311  }
312  }
313  }
314  else if (!haveMesh[Pstream::myProcNo()])
315  {
316  // Don't have mesh (nor fields). Receive empty field from master.
317 
318  forAll(masterNames, i)
319  {
320  const word& name = masterNames[i];
321 
322  // Receive field
323  IPstream fromMaster(Pstream::blocking, Pstream::masterNo());
324  dictionary fieldDict(fromMaster);
325 
326  fields.set
327  (
328  i,
329  new GeoField
330  (
331  IOobject
332  (
333  name,
334  mesh.time().timeName(),
335  mesh,
336  IOobject::NO_READ,
337  IOobject::AUTO_WRITE
338  ),
339  mesh,
340  fieldDict
341  )
342  );
343 
345  //fields[i].write();
346  }
347  }
348  else
349  {
350  // Have mesh so just try to load
351  forAll(masterNames, i)
352  {
353  const word& name = masterNames[i];
354  IOobject& io = *objects[name];
355  io.writeOpt() = IOobject::AUTO_WRITE;
356 
357  // Load field
358  fields.set(i, new GeoField(io, mesh));
359  }
360  }
361 }
362 
363 
364 // Debugging: compare two fields.
365 void compareFields
366 (
367  const scalar tolDim,
368  const volVectorField& a,
369  const volVectorField& b
370 )
371 {
372  forAll(a, celli)
373  {
374  if (mag(b[celli] - a[celli]) > tolDim)
375  {
377  << "Did not map volVectorField correctly:" << nl
378  << "cell:" << celli
379  << " transfer b:" << b[celli]
380  << " real cc:" << a[celli]
381  << abort(FatalError);
382  }
383  }
384  forAll(a.boundaryField(), patchi)
385  {
386  // We have real mesh cellcentre and
387  // mapped original cell centre.
388 
389  const fvPatchVectorField& aBoundary =
390  a.boundaryField()[patchi];
391 
392  const fvPatchVectorField& bBoundary =
393  b.boundaryField()[patchi];
394 
395  if (!bBoundary.coupled())
396  {
397  forAll(aBoundary, i)
398  {
399  if (mag(aBoundary[i] - bBoundary[i]) > tolDim)
400  {
402  << "Did not map volVectorField correctly:"
403  << endl
404  << "patch:" << patchi << " patchFace:" << i
405  << " cc:" << endl
406  << " real :" << aBoundary[i] << endl
407  << " mapped :" << bBoundary[i] << endl
408  << "This might be just a precision entry"
409  << " on writing the mesh." << endl;
410  //<< abort(FatalError);
411  }
412  }
413  }
414  }
415 }
416 
417 
418 
419 int main(int argc, char *argv[])
420 {
421  #include "addRegionOption.H"
422  #include "addOverwriteOption.H"
423  argList::addOption
424  (
425  "mergeTol",
426  "scalar",
427  "specify the merge distance relative to the bounding box size "
428  "(default 1e-6)"
429  );
430  // Include explicit constant options, have zero from time range
431  timeSelector::addOptions();
432 
433  #include "setRootCase.H"
434 
435  if (env("FOAM_SIGFPE"))
436  {
438  << "Detected floating point exception trapping (FOAM_SIGFPE)."
439  << " This might give" << nl
440  << " problems when mapping fields. Switch it off in case"
441  << " of problems." << endl;
442  }
443 
444 
445  // Create processor directory if non-existing
446  if (!Pstream::master() && !isDir(args.path()))
447  {
448  Pout<< "Creating case directory " << args.path() << endl;
449  mkDir(args.path());
450  }
451 
452 
453  // Make sure we do not use the master-only reading.
454  regIOobject::fileModificationChecking = regIOobject::timeStamp;
455  #include "createTime.H"
456  // Allow override of time
457  instantList times = timeSelector::selectIfPresent(runTime, args);
458  runTime.setTime(times[0], 0);
459 
460  runTime.functionObjects().off();
461 
462  word regionName = polyMesh::defaultRegion;
463  fileName meshSubDir;
464 
465  if (args.optionReadIfPresent("region", regionName))
466  {
467  meshSubDir = regionName/polyMesh::meshSubDir;
468  }
469  else
470  {
471  meshSubDir = polyMesh::meshSubDir;
472  }
473  Info<< "Using mesh subdirectory " << meshSubDir << nl << endl;
474 
475  const bool overwrite = args.optionFound("overwrite");
476 
477 
478  // Get time instance directory. Since not all processors have meshes
479  // just use the master one everywhere.
480 
481  fileName masterInstDir;
482  if (Pstream::master())
483  {
484  masterInstDir = runTime.findInstance(meshSubDir, "points");
485  }
486  Pstream::scatter(masterInstDir);
487 
488  // Check who has a mesh
489  const fileName meshPath = runTime.path()/masterInstDir/meshSubDir;
490 
491  Info<< "Found points in " << meshPath << nl << endl;
492 
493 
494  boolList haveMesh(Pstream::nProcs(), false);
495  haveMesh[Pstream::myProcNo()] = isDir(meshPath);
496  Pstream::gatherList(haveMesh);
497  Pstream::scatterList(haveMesh);
498  Info<< "Per processor mesh availability : " << haveMesh << endl;
499  const bool allHaveMesh = (findIndex(haveMesh, false) == -1);
500 
501  autoPtr<fvMesh> meshPtr = loadOrCreateMesh
502  (
503  IOobject
504  (
505  regionName,
506  masterInstDir,
507  runTime,
509  )
510  );
511 
512  fvMesh& mesh = meshPtr();
513 
514  // Print some statistics
515  Info<< "Before distribution:" << endl;
516  printMeshData(mesh);
517 
518 
519 
520  IOdictionary decompositionDict
521  (
522  IOobject
523  (
524  "decomposeParDict",
525  runTime.system(),
526  mesh,
527  IOobject::MUST_READ_IF_MODIFIED,
528  IOobject::NO_WRITE
529  )
530  );
531 
532  labelList finalDecomp;
533 
534  // Create decompositionMethod and new decomposition
535  {
536  autoPtr<decompositionMethod> decomposer
537  (
539  (
540  decompositionDict
541  )
542  );
543 
544  if (!decomposer().parallelAware())
545  {
547  << "You have selected decomposition method "
548  << decomposer().typeName
549  << " which does" << endl
550  << "not synchronise the decomposition across"
551  << " processor patches." << endl
552  << " You might want to select a decomposition method which"
553  << " is aware of this. Continuing."
554  << endl;
555  }
556 
557  finalDecomp = decomposer().decompose(mesh, mesh.cellCentres());
558  }
559 
560  // Dump decomposition to volScalarField
561  if (!overwrite)
562  {
563  writeDecomposition("decomposition", mesh, finalDecomp);
564  }
565 
566 
567  // Create 0 sized mesh to do all the generation of zero sized
568  // fields on processors that have zero sized meshes. Note that this is
569  // only nessecary on master but since polyMesh construction with
570  // Pstream::parRun does parallel comms we have to do it on all
571  // processors
572  autoPtr<fvMeshSubset> subsetterPtr;
573 
574  if (!allHaveMesh)
575  {
576  // Find last non-processor patch.
577  const polyBoundaryMesh& patches = mesh.boundaryMesh();
578 
579  label nonProci = -1;
580 
581  forAll(patches, patchi)
582  {
583  if (isA<processorPolyPatch>(patches[patchi]))
584  {
585  break;
586  }
587  nonProci++;
588  }
589 
590  if (nonProci == -1)
591  {
593  << "Cannot find non-processor patch on processor "
594  << Pstream::myProcNo() << endl
595  << " Current patches:" << patches.names() << abort(FatalError);
596  }
597 
598  // Subset 0 cells, no parallel comms. This is used to create zero-sized
599  // fields.
600  subsetterPtr.reset(new fvMeshSubset(mesh));
601  subsetterPtr().setLargeCellSubset(labelHashSet(0), nonProci, false);
602  }
603 
604 
605  // Get original objects (before incrementing time!)
606  IOobjectList objects(mesh, runTime.timeName());
607  // We don't want to map the decomposition (mapping already tested when
608  // mapping the cell centre field)
609  IOobjectList::iterator iter = objects.find("decomposition");
610  if (iter != objects.end())
611  {
612  objects.erase(iter);
613  }
614 
615 
616  // volFields
617 
618  PtrList<volScalarField> volScalarFields;
619  readFields
620  (
621  haveMesh,
622  mesh,
623  subsetterPtr,
624  objects,
625  volScalarFields
626  );
627 
628  PtrList<volVectorField> volVectorFields;
629  readFields
630  (
631  haveMesh,
632  mesh,
633  subsetterPtr,
634  objects,
635  volVectorFields
636  );
637 
638  PtrList<volSphericalTensorField> volSphereTensorFields;
639  readFields
640  (
641  haveMesh,
642  mesh,
643  subsetterPtr,
644  objects,
645  volSphereTensorFields
646  );
647 
648  PtrList<volSymmTensorField> volSymmTensorFields;
649  readFields
650  (
651  haveMesh,
652  mesh,
653  subsetterPtr,
654  objects,
655  volSymmTensorFields
656  );
657 
658  PtrList<volTensorField> volTensorFields;
659  readFields
660  (
661  haveMesh,
662  mesh,
663  subsetterPtr,
664  objects,
665  volTensorFields
666  );
667 
668 
669  // surfaceFields
670 
671  PtrList<surfaceScalarField> surfScalarFields;
672  readFields
673  (
674  haveMesh,
675  mesh,
676  subsetterPtr,
677  objects,
678  surfScalarFields
679  );
680 
681  PtrList<surfaceVectorField> surfVectorFields;
682  readFields
683  (
684  haveMesh,
685  mesh,
686  subsetterPtr,
687  objects,
688  surfVectorFields
689  );
690 
691  PtrList<surfaceSphericalTensorField> surfSphereTensorFields;
692  readFields
693  (
694  haveMesh,
695  mesh,
696  subsetterPtr,
697  objects,
698  surfSphereTensorFields
699  );
700 
701  PtrList<surfaceSymmTensorField> surfSymmTensorFields;
702  readFields
703  (
704  haveMesh,
705  mesh,
706  subsetterPtr,
707  objects,
708  surfSymmTensorFields
709  );
710 
711  PtrList<surfaceTensorField> surfTensorFields;
712  readFields
713  (
714  haveMesh,
715  mesh,
716  subsetterPtr,
717  objects,
718  surfTensorFields
719  );
720 
721 
722  // Debugging: Create additional volField that will be mapped.
723  // Used to test correctness of mapping
724  //volVectorField mapCc("mapCc", 1*mesh.C());
725 
726  // Global matching tolerance
727  const scalar tolDim = getMergeDistance
728  (
729  args,
730  runTime,
731  mesh.bounds()
732  );
733 
734  // Mesh distribution engine
735  fvMeshDistribute distributor(mesh, tolDim);
736 
737  //Pout<< "Wanted distribution:"
738  // << distributor.countCells(finalDecomp) << nl << endl;
739 
740  // Do actual sending/receiving of mesh
741  autoPtr<mapDistributePolyMesh> map = distributor.distribute(finalDecomp);
742 
744  //map().distributeFaceData(faceCc);
745 
746 
747  // Print some statistics
748  Info<< "After distribution:" << endl;
749  printMeshData(mesh);
750 
751 
752  if (!overwrite)
753  {
754  runTime++;
755  }
756  else
757  {
758  mesh.setInstance(masterInstDir);
759  }
760  Info<< "Writing redistributed mesh to " << runTime.timeName() << nl << endl;
761  mesh.write();
762 
763 
764  // Debugging: test mapped cellcentre field.
765  //compareFields(tolDim, mesh.C(), mapCc);
766 
767  // Print nice message
768  // ~~~~~~~~~~~~~~~~~~
769 
770  // Determine which processors remain so we can print nice final message.
771  labelList nFaces(Pstream::nProcs());
772  nFaces[Pstream::myProcNo()] = mesh.nFaces();
773  Pstream::gatherList(nFaces);
774  Pstream::scatterList(nFaces);
775 
776  Info<< nl
777  << "You can pick up the redecomposed mesh from the polyMesh directory"
778  << " in " << runTime.timeName() << "." << nl
779  << "If you redecomposed the mesh to less processors you can delete"
780  << nl
781  << "the processor directories with 0 sized meshes in them." << nl
782  << "Below is a sample set of commands to do this."
783  << " Take care when issuing these" << nl
784  << "commands." << nl << endl;
785 
786  forAll(nFaces, proci)
787  {
788  fileName procDir = "processor" + name(proci);
789 
790  if (nFaces[proci] == 0)
791  {
792  Info<< " rm -r " << procDir.c_str() << nl;
793  }
794  else
795  {
796  fileName timeDir = procDir/runTime.timeName()/meshSubDir;
797  fileName constDir = procDir/runTime.constant()/meshSubDir;
798 
799  Info<< " rm -r " << constDir.c_str() << nl
800  << " mv " << timeDir.c_str()
801  << ' ' << constDir.c_str() << nl;
802  }
803  }
804  Info<< endl;
805 
806 
807  Info<< "End\n" << endl;
808 
809  return 0;
810 }
811 
812 
813 // ************************************************************************* //
List< instant > instantList
List of instants.
Definition: instantList.H:42
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:57
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Inter-processor communication reduction functions.
fvPatchField< vector > fvPatchVectorField
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
static fvMesh * meshPtr
Definition: globalFoam.H:52
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
bool isDir(const fileName &)
Does the name exist as a DIRECTORY in the file system?
Definition: POSIX.C:486
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:55
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const HashSet< word > &selectedFields, LIFOStack< regIOobject * > &storedObjects)
Read the selected GeometricFields of the specified type.
Definition: ReadFields.C:244
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
List< bool > boolList
Bool container classes.
Definition: boolList.H:50
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
HashSet< label, Hash< label > > labelHashSet
A HashSet with label keys.
Definition: HashSet.H:210
const dimensionedScalar e
Elementary charge.
Definition: doubleFloat.H:78
dynamicFvMesh & mesh
Load or create (0 size) a mesh. Used in distributing meshes to a larger number of processors...
List< label > labelList
A List of labels.
Definition: labelList.H:56
errorManip< error > abort(error &err)
Definition: errorManip.H:131
prefixOSstream Pout(cout,"Pout")
Definition: IOstreams.H:53
static const char nl
Definition: Ostream.H:262
autoPtr< fvMesh > loadOrCreateMesh(const IOobject &io)
Load (if it exists) or create zero cell mesh given an IOobject:
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurence of given element and return index,.
bool mkDir(const fileName &, mode_t=0777)
Make a directory and return an error if it could not be created.
Definition: POSIX.C:295
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
List< word > wordList
A List of words.
Definition: fileName.H:54
void setSize(const label)
Reset size of List.
Definition: List.C:295
label patchi
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
#define WarningInFunction
Report a warning using Foam::Warning.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
messageStream Info
dimensioned< scalar > mag(const dimensioned< Type > &)
virtual Ostream & write(const token &)=0
Write next token to stream.
bool env(const word &)
Return true if environment variable of given name is defined.
Definition: POSIX.C:96