redistributePar.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2019 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 typename GeoField::Mesh& 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::commsTypes::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
324  (
325  Pstream::commsTypes::blocking,
326  Pstream::masterNo()
327  );
328  dictionary fieldDict(fromMaster);
329 
330  fields.set
331  (
332  i,
333  new GeoField
334  (
335  IOobject
336  (
337  name,
338  mesh.thisDb().time().timeName(),
339  mesh.thisDb(),
340  IOobject::NO_READ,
341  IOobject::AUTO_WRITE
342  ),
343  mesh,
344  fieldDict
345  )
346  );
347 
349  // fields[i].write();
350  }
351  }
352  else
353  {
354  // Have mesh so just try to load
355  forAll(masterNames, i)
356  {
357  const word& name = masterNames[i];
358  IOobject& io = *objects[name];
359  io.writeOpt() = IOobject::AUTO_WRITE;
360 
361  // Load field
362  fields.set(i, new GeoField(io, mesh));
363  }
364  }
365 }
366 
367 
368 // Debugging: compare two fields.
369 void compareFields
370 (
371  const scalar tolDim,
372  const volVectorField& a,
373  const volVectorField& b
374 )
375 {
376  forAll(a, celli)
377  {
378  if (mag(b[celli] - a[celli]) > tolDim)
379  {
381  << "Did not map volVectorField correctly:" << nl
382  << "cell:" << celli
383  << " transfer b:" << b[celli]
384  << " real cc:" << a[celli]
385  << abort(FatalError);
386  }
387  }
388  forAll(a.boundaryField(), patchi)
389  {
390  // We have real mesh cellcentre and
391  // mapped original cell centre.
392 
393  const fvPatchVectorField& aBoundary =
394  a.boundaryField()[patchi];
395 
396  const fvPatchVectorField& bBoundary =
397  b.boundaryField()[patchi];
398 
399  if (!bBoundary.coupled())
400  {
401  forAll(aBoundary, i)
402  {
403  if (mag(aBoundary[i] - bBoundary[i]) > tolDim)
404  {
406  << "Did not map volVectorField correctly:"
407  << endl
408  << "patch:" << patchi << " patchFace:" << i
409  << " cc:" << endl
410  << " real :" << aBoundary[i] << endl
411  << " mapped :" << bBoundary[i] << endl
412  << "This might be just a precision entry"
413  << " on writing the mesh." << endl;
414  //<< abort(FatalError);
415  }
416  }
417  }
418  }
419 }
420 
421 
422 
423 int main(int argc, char *argv[])
424 {
425  #include "addRegionOption.H"
426  #include "addOverwriteOption.H"
427  argList::addOption
428  (
429  "mergeTol",
430  "scalar",
431  "specify the merge distance relative to the bounding box size "
432  "(default 1e-6)"
433  );
434  // Include explicit constant options, have zero from time range
435  timeSelector::addOptions();
436 
437  #include "setRootCase.H"
438 
439  if (env("FOAM_SIGFPE"))
440  {
442  << "Detected floating point exception trapping (FOAM_SIGFPE)."
443  << " This might give" << nl
444  << " problems when mapping fields. Switch it off in case"
445  << " of problems." << endl;
446  }
447 
448 
449  // Create processor directory if non-existing
450  if (!Pstream::master() && !isDir(args.path()))
451  {
452  Pout<< "Creating case directory " << args.path() << endl;
453  mkDir(args.path());
454  }
455 
456 
457  // Make sure we do not use the master-only reading.
458  regIOobject::fileModificationChecking = regIOobject::timeStamp;
459  #include "createTime.H"
460  // Allow override of time
461  instantList times = timeSelector::selectIfPresent(runTime, args);
462  runTime.setTime(times[0], 0);
463 
464  runTime.functionObjects().off();
465 
466  word regionName = polyMesh::defaultRegion;
467  fileName meshSubDir;
468 
469  if (args.optionReadIfPresent("region", regionName))
470  {
471  meshSubDir = regionName/polyMesh::meshSubDir;
472  }
473  else
474  {
475  meshSubDir = polyMesh::meshSubDir;
476  }
477  Info<< "Using mesh subdirectory " << meshSubDir << nl << endl;
478 
479  const bool overwrite = args.optionFound("overwrite");
480 
481 
482  // Get time instance directory. Since not all processors have meshes
483  // just use the master one everywhere.
484 
485  fileName masterInstDir;
486  if (Pstream::master())
487  {
488  masterInstDir = runTime.findInstance(meshSubDir, "points");
489  }
490  Pstream::scatter(masterInstDir);
491 
492  // Check who has a mesh
493  const fileName meshPath = runTime.path()/masterInstDir/meshSubDir;
494 
495  Info<< "Found points in " << meshPath << nl << endl;
496 
497 
498  boolList haveMesh(Pstream::nProcs(), false);
499  haveMesh[Pstream::myProcNo()] = isDir(meshPath);
500  Pstream::gatherList(haveMesh);
501  Pstream::scatterList(haveMesh);
502  Info<< "Per processor mesh availability : " << haveMesh << endl;
503  const bool allHaveMesh = (findIndex(haveMesh, false) == -1);
504 
505  autoPtr<fvMesh> meshPtr = loadOrCreateMesh
506  (
507  IOobject
508  (
509  regionName,
510  masterInstDir,
511  runTime,
513  )
514  );
515 
516  fvMesh& mesh = meshPtr();
517 
518  // Print some statistics
519  Info<< "Before distribution:" << endl;
520  printMeshData(mesh);
521 
522 
523 
524  IOdictionary decompositionDict
525  (
526  IOobject
527  (
528  "decomposeParDict",
529  runTime.system(),
530  mesh,
531  IOobject::MUST_READ_IF_MODIFIED,
532  IOobject::NO_WRITE
533  )
534  );
535 
536  labelList finalDecomp;
537 
538  // Create decompositionMethod and new decomposition
539  {
540  autoPtr<decompositionMethod> decomposer
541  (
543  (
544  decompositionDict
545  )
546  );
547 
548  if (!decomposer().parallelAware())
549  {
551  << "You have selected decomposition method "
552  << decomposer().typeName
553  << " which does" << endl
554  << "not synchronise the decomposition across"
555  << " processor patches." << endl
556  << " You might want to select a decomposition method which"
557  << " is aware of this. Continuing."
558  << endl;
559  }
560 
561  finalDecomp = decomposer().decompose(mesh, mesh.cellCentres());
562  }
563 
564  // Dump decomposition to volScalarField
565  if (!overwrite)
566  {
567  writeDecomposition("decomposition", mesh, finalDecomp);
568  }
569 
570 
571  // Create 0 sized mesh to do all the generation of zero sized
572  // fields on processors that have zero sized meshes. Note that this is
573  // only necessary on master but since polyMesh construction with
574  // Pstream::parRun does parallel comms we have to do it on all
575  // processors
576  autoPtr<fvMeshSubset> subsetterPtr;
577 
578  if (!allHaveMesh)
579  {
580  // Find last non-processor patch.
581  const polyBoundaryMesh& patches = mesh.boundaryMesh();
582 
583  label nonProci = -1;
584 
585  forAll(patches, patchi)
586  {
587  if (isA<processorPolyPatch>(patches[patchi]))
588  {
589  break;
590  }
591  nonProci++;
592  }
593 
594  if (nonProci == -1)
595  {
597  << "Cannot find non-processor patch on processor "
598  << Pstream::myProcNo() << endl
599  << " Current patches:" << patches.names() << abort(FatalError);
600  }
601 
602  // Subset 0 cells, no parallel comms. This is used to create zero-sized
603  // fields.
604  subsetterPtr.reset(new fvMeshSubset(mesh));
605  subsetterPtr().setLargeCellSubset(labelHashSet(0), nonProci, false);
606  }
607 
608 
609  // Get original objects (before incrementing time!)
610  IOobjectList objects(mesh, runTime.timeName());
611  // We don't want to map the decomposition (mapping already tested when
612  // mapping the cell centre field)
613  IOobjectList::iterator iter = objects.find("decomposition");
614  if (iter != objects.end())
615  {
616  objects.erase(iter);
617  }
618 
619 
620  // volFields
621 
622  PtrList<volScalarField> volScalarFields;
623  readFields
624  (
625  haveMesh,
626  mesh,
627  subsetterPtr,
628  objects,
629  volScalarFields
630  );
631 
632  PtrList<volVectorField> volVectorFields;
633  readFields
634  (
635  haveMesh,
636  mesh,
637  subsetterPtr,
638  objects,
639  volVectorFields
640  );
641 
642  PtrList<volSphericalTensorField> volSphereTensorFields;
643  readFields
644  (
645  haveMesh,
646  mesh,
647  subsetterPtr,
648  objects,
649  volSphereTensorFields
650  );
651 
652  PtrList<volSymmTensorField> volSymmTensorFields;
653  readFields
654  (
655  haveMesh,
656  mesh,
657  subsetterPtr,
658  objects,
659  volSymmTensorFields
660  );
661 
662  PtrList<volTensorField> volTensorFields;
663  readFields
664  (
665  haveMesh,
666  mesh,
667  subsetterPtr,
668  objects,
669  volTensorFields
670  );
671 
672 
673  // surfaceFields
674 
675  PtrList<surfaceScalarField> surfScalarFields;
676  readFields
677  (
678  haveMesh,
679  mesh,
680  subsetterPtr,
681  objects,
682  surfScalarFields
683  );
684 
685  PtrList<surfaceVectorField> surfVectorFields;
686  readFields
687  (
688  haveMesh,
689  mesh,
690  subsetterPtr,
691  objects,
692  surfVectorFields
693  );
694 
695  PtrList<surfaceSphericalTensorField> surfSphereTensorFields;
696  readFields
697  (
698  haveMesh,
699  mesh,
700  subsetterPtr,
701  objects,
702  surfSphereTensorFields
703  );
704 
705  PtrList<surfaceSymmTensorField> surfSymmTensorFields;
706  readFields
707  (
708  haveMesh,
709  mesh,
710  subsetterPtr,
711  objects,
712  surfSymmTensorFields
713  );
714 
715  PtrList<surfaceTensorField> surfTensorFields;
716  readFields
717  (
718  haveMesh,
719  mesh,
720  subsetterPtr,
721  objects,
722  surfTensorFields
723  );
724 
725 
726  // pointFields
727 
728  pointMesh& pMesh =
729  const_cast<pointMesh&>(pointMesh::New(mesh));
730  PtrList<pointScalarField> pointScalarFields;
731  readFields
732  (
733  haveMesh,
734  pMesh,
735  subsetterPtr,
736  objects,
737  pointScalarFields
738  );
739 
740  PtrList<pointVectorField> pointVectorFields;
741  readFields
742  (
743  haveMesh,
744  pMesh,
745  subsetterPtr,
746  objects,
747  pointVectorFields
748  );
749 
750  PtrList<pointSphericalTensorField> pointSphereTensorFields;
751  readFields
752  (
753  haveMesh,
754  pMesh,
755  subsetterPtr,
756  objects,
757  pointSphereTensorFields
758  );
759 
760  PtrList<pointSymmTensorField> pointSymmTensorFields;
761  readFields
762  (
763  haveMesh,
764  pMesh,
765  subsetterPtr,
766  objects,
767  pointSymmTensorFields
768  );
769 
770  PtrList<pointTensorField> pointTensorFields;
771  readFields
772  (
773  haveMesh,
774  pMesh,
775  subsetterPtr,
776  objects,
777  pointTensorFields
778  );
779 
780  // Debugging: Create additional volField that will be mapped.
781  // Used to test correctness of mapping
782  // volVectorField mapCc("mapCc", 1*mesh.C());
783 
784  // Global matching tolerance
785  const scalar tolDim = getMergeDistance
786  (
787  args,
788  runTime,
789  mesh.bounds()
790  );
791 
792  // Mesh distribution engine
793  fvMeshDistribute distributor(mesh, tolDim);
794 
795  // Pout<< "Wanted distribution:"
796  // << distributor.countCells(finalDecomp) << nl << endl;
797 
798  // Do actual sending/receiving of mesh
799  autoPtr<mapDistributePolyMesh> map = distributor.distribute(finalDecomp);
800 
802  // map().distributeFaceData(faceCc);
803 
804 
805  // Print some statistics
806  Info<< "After distribution:" << endl;
807  printMeshData(mesh);
808 
809 
810  if (!overwrite)
811  {
812  runTime++;
813  }
814  else
815  {
816  mesh.setInstance(masterInstDir);
817  }
818  Info<< "Writing redistributed mesh to " << runTime.timeName() << nl << endl;
819  mesh.write();
820 
821 
822  // Debugging: test mapped cellcentre field.
823  // compareFields(tolDim, mesh.C(), mapCc);
824 
825  // Print nice message
826  // ~~~~~~~~~~~~~~~~~~
827 
828  // Determine which processors remain so we can print nice final message.
829  labelList nFaces(Pstream::nProcs());
830  nFaces[Pstream::myProcNo()] = mesh.nFaces();
831  Pstream::gatherList(nFaces);
832  Pstream::scatterList(nFaces);
833 
834  Info<< nl
835  << "You can pick up the redecomposed mesh from the polyMesh directory"
836  << " in " << runTime.timeName() << "." << nl
837  << "If you redecomposed the mesh to less processors you can delete"
838  << nl
839  << "the processor directories with 0 sized meshes in them." << nl
840  << "Below is a sample set of commands to do this."
841  << " Take care when issuing these" << nl
842  << "commands." << nl << endl;
843 
844  forAll(nFaces, proci)
845  {
846  fileName procDir = "processor" + name(proci);
847 
848  if (nFaces[proci] == 0)
849  {
850  Info<< " rm -r " << procDir.c_str() << nl;
851  }
852  else
853  {
854  fileName timeDir = procDir/runTime.timeName()/meshSubDir;
855  fileName constDir = procDir/runTime.constant()/meshSubDir;
856 
857  Info<< " rm -r " << constDir.c_str() << nl
858  << " mv " << timeDir.c_str()
859  << ' ' << constDir.c_str() << nl;
860  }
861  }
862  Info<< endl;
863 
864 
865  Info<< "End\n" << endl;
866 
867  return 0;
868 }
869 
870 
871 // ************************************************************************* //
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:434
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
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
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:256
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:55
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:208
dynamicFvMesh & mesh
volScalarField & e
Elementary charge.
Definition: createFields.H:11
bool isDir(const fileName &, const bool followLink=true)
Does the name exist as a directory in the file system?
Definition: POSIX.C:539
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
static const char nl
Definition: Ostream.H:265
objects
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 occurrence 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:290
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:281
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.
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
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:91