fvMeshTools.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) 2012-2016 OpenFOAM Foundation
9  Copyright (C) 2015-2019 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 \*---------------------------------------------------------------------------*/
28 
29 #include "fvMeshTools.H"
32 
33 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 
35 // Adds patch if not yet there. Returns patchID.
37 (
38  fvMesh& mesh,
39  const polyPatch& patch,
40  const dictionary& patchFieldDict,
41  const word& defaultPatchFieldType,
42  const bool validBoundary
43 )
44 {
45  polyBoundaryMesh& polyPatches =
46  const_cast<polyBoundaryMesh&>(mesh.boundaryMesh());
47 
48  label patchi = polyPatches.findPatchID(patch.name());
49  if (patchi != -1)
50  {
51  // Already there
52  return patchi;
53  }
54 
55 
56  // Append at end unless there are processor patches
57  label insertPatchi = polyPatches.size();
58  label startFacei = mesh.nFaces();
59 
60  if (!isA<processorPolyPatch>(patch))
61  {
62  forAll(polyPatches, patchi)
63  {
64  const polyPatch& pp = polyPatches[patchi];
65 
66  if (isA<processorPolyPatch>(pp))
67  {
68  insertPatchi = patchi;
69  startFacei = pp.start();
70  break;
71  }
72  }
73  }
74 
75 
76  // Below is all quite a hack. Feel free to change once there is a better
77  // mechanism to insert and reorder patches.
78 
79  // Clear local fields and e.g. polyMesh parallelInfo.
80  mesh.clearOut();
81 
82  label sz = polyPatches.size();
83 
84  fvBoundaryMesh& fvPatches = const_cast<fvBoundaryMesh&>(mesh.boundary());
85 
86  // Add polyPatch at the end
87  polyPatches.setSize(sz+1);
88  polyPatches.set
89  (
90  sz,
91  patch.clone
92  (
93  polyPatches,
94  insertPatchi, //index
95  0, //size
96  startFacei //start
97  )
98  );
99  fvPatches.setSize(sz+1);
100  fvPatches.set
101  (
102  sz,
104  (
105  polyPatches[sz], // point to newly added polyPatch
106  mesh.boundary()
107  )
108  );
109 
110  addPatchFields<volScalarField>
111  (
112  mesh,
113  patchFieldDict,
114  defaultPatchFieldType,
115  Zero
116  );
117  addPatchFields<volVectorField>
118  (
119  mesh,
120  patchFieldDict,
121  defaultPatchFieldType,
122  Zero
123  );
124  addPatchFields<volSphericalTensorField>
125  (
126  mesh,
127  patchFieldDict,
128  defaultPatchFieldType,
129  Zero
130  );
131  addPatchFields<volSymmTensorField>
132  (
133  mesh,
134  patchFieldDict,
135  defaultPatchFieldType,
136  Zero
137  );
138  addPatchFields<volTensorField>
139  (
140  mesh,
141  patchFieldDict,
142  defaultPatchFieldType,
143  Zero
144  );
145 
146  // Surface fields
147 
148  addPatchFields<surfaceScalarField>
149  (
150  mesh,
151  patchFieldDict,
152  defaultPatchFieldType,
153  Zero
154  );
155  addPatchFields<surfaceVectorField>
156  (
157  mesh,
158  patchFieldDict,
159  defaultPatchFieldType,
160  Zero
161  );
162  addPatchFields<surfaceSphericalTensorField>
163  (
164  mesh,
165  patchFieldDict,
166  defaultPatchFieldType,
167  Zero
168  );
169  addPatchFields<surfaceSymmTensorField>
170  (
171  mesh,
172  patchFieldDict,
173  defaultPatchFieldType,
174  Zero
175  );
176  addPatchFields<surfaceTensorField>
177  (
178  mesh,
179  patchFieldDict,
180  defaultPatchFieldType,
181  Zero
182  );
183 
184  // Create reordering list
185  // patches before insert position stay as is
186  labelList oldToNew(sz+1);
187  for (label i = 0; i < insertPatchi; i++)
188  {
189  oldToNew[i] = i;
190  }
191  // patches after insert position move one up
192  for (label i = insertPatchi; i < sz; i++)
193  {
194  oldToNew[i] = i+1;
195  }
196  // appended patch gets moved to insert position
197  oldToNew[sz] = insertPatchi;
198 
199  // Shuffle into place
200  polyPatches.reorder(oldToNew, validBoundary);
201  fvPatches.reorder(oldToNew);
202 
203  reorderPatchFields<volScalarField>(mesh, oldToNew);
204  reorderPatchFields<volVectorField>(mesh, oldToNew);
205  reorderPatchFields<volSphericalTensorField>(mesh, oldToNew);
206  reorderPatchFields<volSymmTensorField>(mesh, oldToNew);
207  reorderPatchFields<volTensorField>(mesh, oldToNew);
208  reorderPatchFields<surfaceScalarField>(mesh, oldToNew);
209  reorderPatchFields<surfaceVectorField>(mesh, oldToNew);
210  reorderPatchFields<surfaceSphericalTensorField>(mesh, oldToNew);
211  reorderPatchFields<surfaceSymmTensorField>(mesh, oldToNew);
212  reorderPatchFields<surfaceTensorField>(mesh, oldToNew);
213 
214  return insertPatchi;
215 }
216 
217 
218 void Foam::fvMeshTools::setPatchFields
219 (
220  fvMesh& mesh,
221  const label patchi,
222  const dictionary& patchFieldDict
223 )
224 {
225  setPatchFields<volScalarField>(mesh, patchi, patchFieldDict);
226  setPatchFields<volVectorField>(mesh, patchi, patchFieldDict);
227  setPatchFields<volSphericalTensorField>(mesh, patchi, patchFieldDict);
228  setPatchFields<volSymmTensorField>(mesh, patchi, patchFieldDict);
229  setPatchFields<volTensorField>(mesh, patchi, patchFieldDict);
230  setPatchFields<surfaceScalarField>(mesh, patchi, patchFieldDict);
231  setPatchFields<surfaceVectorField>(mesh, patchi, patchFieldDict);
232  setPatchFields<surfaceSphericalTensorField>
233  (
234  mesh,
235  patchi,
236  patchFieldDict
237  );
238  setPatchFields<surfaceSymmTensorField>(mesh, patchi, patchFieldDict);
239  setPatchFields<surfaceTensorField>(mesh, patchi, patchFieldDict);
240 }
241 
242 
244 {
245  setPatchFields<volScalarField>(mesh, patchi, Zero);
246  setPatchFields<volVectorField>(mesh, patchi, Zero);
247  setPatchFields<volSphericalTensorField>
248  (
249  mesh,
250  patchi,
251  Zero
252  );
253  setPatchFields<volSymmTensorField>
254  (
255  mesh,
256  patchi,
257  Zero
258  );
259  setPatchFields<volTensorField>(mesh, patchi, Zero);
260  setPatchFields<surfaceScalarField>(mesh, patchi, Zero);
261  setPatchFields<surfaceVectorField>(mesh, patchi, Zero);
262  setPatchFields<surfaceSphericalTensorField>
263  (
264  mesh,
265  patchi,
266  Zero
267  );
268  setPatchFields<surfaceSymmTensorField>
269  (
270  mesh,
271  patchi,
272  Zero
273  );
274  setPatchFields<surfaceTensorField>(mesh, patchi, Zero);
275 }
276 
277 
278 // Deletes last patch
279 void Foam::fvMeshTools::trimPatches(fvMesh& mesh, const label nPatches)
280 {
281  // Clear local fields and e.g. polyMesh globalMeshData.
282  mesh.clearOut();
283 
284  polyBoundaryMesh& polyPatches =
285  const_cast<polyBoundaryMesh&>(mesh.boundaryMesh());
286  fvBoundaryMesh& fvPatches = const_cast<fvBoundaryMesh&>(mesh.boundary());
287 
288  if (polyPatches.empty())
289  {
291  << "No patches in mesh"
292  << abort(FatalError);
293  }
294 
295  label nFaces = 0;
296  for (label patchi = nPatches; patchi < polyPatches.size(); patchi++)
297  {
298  nFaces += polyPatches[patchi].size();
299  }
300  reduce(nFaces, sumOp<label>());
301 
302  if (nFaces)
303  {
305  << "There are still " << nFaces
306  << " faces in " << polyPatches.size()-nPatches
307  << " patches to be deleted" << abort(FatalError);
308  }
309 
310  // Remove actual patches
311  polyPatches.setSize(nPatches);
312  fvPatches.setSize(nPatches);
313 
314  trimPatchFields<volScalarField>(mesh, nPatches);
315  trimPatchFields<volVectorField>(mesh, nPatches);
316  trimPatchFields<volSphericalTensorField>(mesh, nPatches);
317  trimPatchFields<volSymmTensorField>(mesh, nPatches);
318  trimPatchFields<volTensorField>(mesh, nPatches);
319 
320  trimPatchFields<surfaceScalarField>(mesh, nPatches);
321  trimPatchFields<surfaceVectorField>(mesh, nPatches);
322  trimPatchFields<surfaceSphericalTensorField>(mesh, nPatches);
323  trimPatchFields<surfaceSymmTensorField>(mesh, nPatches);
324  trimPatchFields<surfaceTensorField>(mesh, nPatches);
325 }
326 
327 
329 (
330  fvMesh& mesh,
331  const labelList& oldToNew,
332  const label nNewPatches,
333  const bool validBoundary
334 )
335 {
336  polyBoundaryMesh& polyPatches =
337  const_cast<polyBoundaryMesh&>(mesh.boundaryMesh());
338  fvBoundaryMesh& fvPatches = const_cast<fvBoundaryMesh&>(mesh.boundary());
339 
340  // Shuffle into place
341  polyPatches.reorder(oldToNew, validBoundary);
342  fvPatches.reorder(oldToNew);
343 
344  reorderPatchFields<volScalarField>(mesh, oldToNew);
345  reorderPatchFields<volVectorField>(mesh, oldToNew);
346  reorderPatchFields<volSphericalTensorField>(mesh, oldToNew);
347  reorderPatchFields<volSymmTensorField>(mesh, oldToNew);
348  reorderPatchFields<volTensorField>(mesh, oldToNew);
349  reorderPatchFields<surfaceScalarField>(mesh, oldToNew);
350  reorderPatchFields<surfaceVectorField>(mesh, oldToNew);
351  reorderPatchFields<surfaceSphericalTensorField>(mesh, oldToNew);
352  reorderPatchFields<surfaceSymmTensorField>(mesh, oldToNew);
353  reorderPatchFields<surfaceTensorField>(mesh, oldToNew);
354 
355  // Remove last.
356  trimPatches(mesh, nNewPatches);
357 }
358 
359 
361 (
362  fvMesh& mesh,
363  const bool validBoundary
364 )
365 {
366  const polyBoundaryMesh& pbm = mesh.boundaryMesh();
367 
368  labelList newToOld(pbm.size());
369  labelList oldToNew(pbm.size(), -1);
370  label newI = 0;
371 
372 
373  // Assumes all non-coupled boundaries are on all processors!
374  forAll(pbm, patchI)
375  {
376  const polyPatch& pp = pbm[patchI];
377 
378  if (!isA<processorPolyPatch>(pp))
379  {
380  label nFaces = pp.size();
381  if (validBoundary)
382  {
383  reduce(nFaces, sumOp<label>());
384  }
385 
386  if (nFaces > 0)
387  {
388  newToOld[newI] = patchI;
389  oldToNew[patchI] = newI++;
390  }
391  }
392  }
393 
394  // Same for processor patches (but need no reduction)
395  forAll(pbm, patchI)
396  {
397  const polyPatch& pp = pbm[patchI];
398 
399  if (isA<processorPolyPatch>(pp) && pp.size())
400  {
401  newToOld[newI] = patchI;
402  oldToNew[patchI] = newI++;
403  }
404  }
405 
406  newToOld.setSize(newI);
407 
408  // Move all deletable patches to the end
409  forAll(oldToNew, patchI)
410  {
411  if (oldToNew[patchI] == -1)
412  {
413  oldToNew[patchI] = newI++;
414  }
415  }
416 
417  reorderPatches(mesh, oldToNew, newToOld.size(), validBoundary);
418 
419  return newToOld;
420 }
421 
422 
424 (
425  const IOobject& io,
426  const bool masterOnlyReading
427 )
428 {
429  // Region name
430  // ~~~~~~~~~~~
431 
432  fileName meshSubDir;
433 
434  if (io.name() == polyMesh::defaultRegion)
435  {
436  meshSubDir = polyMesh::meshSubDir;
437  }
438  else
439  {
440  meshSubDir = io.name()/polyMesh::meshSubDir;
441  }
442 
443 
444  fileName facesInstance;
445 
446  // Patch types
447  // ~~~~~~~~~~~
448  // Read and scatter master patches (without reading master mesh!)
449 
450  PtrList<entry> patchEntries;
451  if (Pstream::master())
452  {
453  facesInstance = io.time().findInstance
454  (
455  meshSubDir,
456  "faces",
458  );
459 
460  patchEntries = polyBoundaryMeshEntries
461  (
462  IOobject
463  (
464  "boundary",
465  facesInstance,
466  meshSubDir,
467  io.db(),
470  false
471  )
472  );
473 
474  // Send patches
475  for
476  (
477  int slave=Pstream::firstSlave();
478  slave<=Pstream::lastSlave();
479  slave++
480  )
481  {
483  toSlave << patchEntries;
484  }
485  }
486  else
487  {
488  // Receive patches
489  IPstream fromMaster
490  (
493  );
494  fromMaster >> patchEntries;
495  }
496 
497 
498  Pstream::scatter(facesInstance);
499 
500  // Dummy meshes
501  // ~~~~~~~~~~~~
502 
503  // Check who has a mesh
504  const fileName meshDir = io.time().path()/facesInstance/meshSubDir;
505  bool haveMesh = isDir(meshDir);
506  if (masterOnlyReading && !Pstream::master())
507  {
508  haveMesh = false;
509  }
510 
511 
512  // Set up to read-if-present. Note: does not search for mesh so set
513  // instance explicitly
514  IOobject meshIO(io);
515  meshIO.instance() = facesInstance;
516  if (masterOnlyReading && !Pstream::master())
517  {
518  meshIO.readOpt() = IOobject::NO_READ;
519  }
520  else
521  {
523  }
524 
525 
526  // Read mesh
527  // ~~~~~~~~~
528  // Now all processors read a mesh and use supplied points,faces etc
529  // if there is none.
530  // Note: fvSolution, fvSchemes are also using the supplied Ioobject so
531  // on slave will be NO_READ, on master READ_IF_PRESENT. This will
532  // conflict with e.g. timeStampMaster reading so switch off.
533 
534  const regIOobject::fileCheckTypes oldCheckType =
537 
538  auto meshPtr = autoPtr<fvMesh>::New(meshIO, Zero);
539  fvMesh& mesh = *meshPtr;
540 
542 
543 
544 
545  // Add patches
546  // ~~~~~~~~~~~
547 
548 
549  bool havePatches = mesh.boundary().size();
550  reduce(havePatches, andOp<bool>());
551 
552  if (!havePatches)
553  {
554  // There are one or more processors without full complement of
555  // patches
556 
557  DynamicList<polyPatch*> newPatches;
558 
559  if (haveMesh) //Pstream::master())
560  {
561  forAll(mesh.boundary(), patchI)
562  {
563  newPatches.append
564  (
565  mesh.boundaryMesh()[patchI].clone(mesh.boundaryMesh()).ptr()
566  );
567  }
568  }
569  else
570  {
571  forAll(patchEntries, patchI)
572  {
573  const entry& e = patchEntries[patchI];
574  const word type(e.dict().get<word>("type"));
575 
576  if
577  (
578  type == processorPolyPatch::typeName
579  || type == processorCyclicPolyPatch::typeName
580  )
581  {}
582  else
583  {
584  dictionary patchDict(e.dict());
585  patchDict.set("nFaces", 0);
586  patchDict.set("startFace", 0);
587 
588  newPatches.append
589  (
591  (
592  patchEntries[patchI].keyword(),
593  patchDict,
594  newPatches.size(),
596  ).ptr()
597  );
598  }
599  }
600  }
602  mesh.addFvPatches(newPatches);
603  }
604 
605  //Pout<< "patches:" << endl;
606  //forAll(mesh.boundary(), patchI)
607  //{
608  // Pout<< " type" << mesh.boundary()[patchI].type()
609  // << " size:" << mesh.boundary()[patchI].size()
610  // << " start:" << mesh.boundary()[patchI].start() << endl;
611  //}
612  //Pout<< endl;
613 
614 
615  // Determine zones
616  // ~~~~~~~~~~~~~~~
617 
618  wordList pointZoneNames(mesh.pointZones().names());
619  Pstream::scatter(pointZoneNames);
620  wordList faceZoneNames(mesh.faceZones().names());
621  Pstream::scatter(faceZoneNames);
622  wordList cellZoneNames(mesh.cellZones().names());
623  Pstream::scatter(cellZoneNames);
624 
625  if (!haveMesh)
626  {
627  // Add the zones. Make sure to remove the old dummy ones first
628  mesh.pointZones().clear();
629  mesh.faceZones().clear();
630  mesh.cellZones().clear();
631 
632  List<pointZone*> pz(pointZoneNames.size());
633  forAll(pointZoneNames, i)
634  {
635  pz[i] = new pointZone
636  (
637  pointZoneNames[i],
638  i,
639  mesh.pointZones()
640  );
641  }
642  List<faceZone*> fz(faceZoneNames.size());
643  forAll(faceZoneNames, i)
644  {
645  fz[i] = new faceZone
646  (
647  faceZoneNames[i],
648  i,
649  mesh.faceZones()
650  );
651  }
652  List<cellZone*> cz(cellZoneNames.size());
653  forAll(cellZoneNames, i)
654  {
655  cz[i] = new cellZone
656  (
657  cellZoneNames[i],
658  i,
659  mesh.cellZones()
660  );
661  }
662 
663  if (pz.size() && fz.size() && cz.size())
664  {
665  mesh.addZones(pz, fz, cz);
666  }
667  }
668 
669  return meshPtr;
670 }
671 
672 
674 (
675  const objectRegistry& mesh,
676  const word& regionName,
677  const bool verbose
678 )
679 {
680  // Create dummy system/fv*
681  {
682  IOobject io
683  (
684  "fvSchemes",
685  mesh.time().system(),
686  regionName,
687  mesh,
690  false
691  );
692 
693  if (!io.typeHeaderOk<IOdictionary>(false))
694  {
695  if (verbose)
696  {
697  Info<< "Writing dummy " << regionName/io.name() << endl;
698  }
699  dictionary dummyDict;
700  dictionary divDict;
701  dummyDict.add("divSchemes", divDict);
702  dictionary gradDict;
703  dummyDict.add("gradSchemes", gradDict);
704  dictionary laplDict;
705  dummyDict.add("laplacianSchemes", laplDict);
706 
707  IOdictionary(io, dummyDict).regIOobject::write();
708  }
709  }
710  {
711  IOobject io
712  (
713  "fvSolution",
714  mesh.time().system(),
715  regionName,
716  mesh,
719  false
720  );
721 
722  if (!io.typeHeaderOk<IOdictionary>(false))
723  {
724  if (verbose)
725  {
726  Info<< "Writing dummy " << regionName/io.name() << endl;
727  }
728  dictionary dummyDict;
729  IOdictionary(io, dummyDict).regIOobject::write();
730  }
731  }
732 }
733 
734 
735 // ************************************************************************* //
Foam::entry
A keyword and a list of tokens is an 'entry'.
Definition: entry.H:67
Foam::IOobject::NO_WRITE
Definition: IOobject.H:130
Foam::autoPtr::New
static autoPtr< T > New(Args &&... args)
Construct autoPtr of T with forwarding arguments.
Foam::IOdictionary
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:54
Foam::ZoneMesh::clear
void clear()
Clear the zones.
Definition: ZoneMesh.C:635
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
meshPtr
Foam::autoPtr< Foam::fvMesh > meshPtr(nullptr)
Foam::IOobject::fileCheckTypes
fileCheckTypes
Enumeration defining the file checking options.
Definition: IOobject.H:134
Foam::UPstream::masterNo
static constexpr int masterNo() noexcept
Process index of the master.
Definition: UPstream.H:432
Foam::fvMeshTools::zeroPatchFields
static void zeroPatchFields(fvMesh &mesh, const label patchi)
Change patchField to zero on registered fields.
Definition: fvMeshTools.C:243
Foam::IOobject::name
const word & name() const
Return name.
Definition: IOobjectI.H:46
Foam::IOobject::timeStamp
Definition: IOobject.H:136
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::fileName
A class for handling file names.
Definition: fileName.H:69
Foam::polyBoundaryMeshEntries
Foam::polyBoundaryMeshEntries.
Definition: polyBoundaryMeshEntries.H:50
Foam::fvMeshTools::createDummyFvMeshFiles
static void createDummyFvMeshFiles(const objectRegistry &parent, const word &regionName, const bool verbose=false)
Create additional fv* files.
Definition: fvMeshTools.C:674
Foam::polyBoundaryMesh
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
Definition: polyBoundaryMesh.H:62
nPatches
label nPatches
Definition: readKivaGrid.H:396
Foam::polyMesh::defaultRegion
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:312
Foam::Zero
static constexpr const zero Zero
Global zero.
Definition: zero.H:128
Foam::DynamicList
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:57
Foam::OPstream
Output inter-processor communications stream.
Definition: OPstream.H:52
Foam::IOobject::typeHeaderOk
bool typeHeaderOk(const bool checkType=true, const bool search=true, const bool verbose=true)
Read header (uses typeFilePath to find file) and check its info.
Definition: IOobjectTemplates.C:39
Foam::IOobject::instance
const fileName & instance() const
Definition: IOobjectI.H:167
Foam::pointZone
A subset of mesh points.
Definition: pointZone.H:65
Foam::fvMeshTools::removeEmptyPatches
static labelList removeEmptyPatches(fvMesh &, const bool validBoundary)
Remove zero sized patches. All but processor patches are.
Definition: fvMeshTools.C:361
Foam::polyMesh::meshSubDir
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
Definition: polyMesh.H:315
Foam::polyBoundaryMesh::reorder
void reorder(const labelUList &oldToNew, const bool validBoundary)
Reorders patches. Ordering does not have to be done in.
Definition: polyBoundaryMesh.C:1192
Foam::IOobject::fileModificationChecking
static fileCheckTypes fileModificationChecking
Type of file modification checking.
Definition: IOobject.H:207
Foam::polyMesh::cellZones
const cellZoneMesh & cellZones() const
Return cell zone mesh.
Definition: polyMesh.H:483
Foam::Pstream::scatter
static void scatter(const List< commsStruct > &comms, T &Value, const int tag, const label comm)
Scatter data. Distribute without modification. Reverse of gather.
Definition: gatherScatter.C:150
Foam::fvMesh::removeFvBoundary
void removeFvBoundary()
Remove boundary patches. Warning: fvPatchFields hold ref to.
Definition: fvMesh.C:522
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:435
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:337
Foam::dictionary::set
entry * set(entry *entryPtr)
Assign a new entry, overwriting any existing entry.
Definition: dictionary.C:842
polyBoundaryMeshEntries.H
Foam::IOobject::time
const Time & time() const
Return time.
Definition: IOobject.C:438
processorCyclicPolyPatch.H
Foam::sumOp
Definition: ops.H:213
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
Foam::cellZone
A subset of mesh cells.
Definition: cellZone.H:62
Foam::polyMesh::faceZones
const faceZoneMesh & faceZones() const
Return face zone mesh.
Definition: polyMesh.H:477
Foam::IOobject::db
const objectRegistry & db() const
Return the local objectRegistry.
Definition: IOobject.C:432
Foam::fvBoundaryMesh
Foam::fvBoundaryMesh.
Definition: fvBoundaryMesh.H:57
Foam::objectRegistry
Registry of regIOobjects.
Definition: objectRegistry.H:60
regionName
Foam::word regionName
Definition: createNamedDynamicFvMesh.H:1
Foam::reduce
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
Definition: PstreamReduceOps.H:51
Foam::label
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:62
Foam::faceZone
A subset of mesh faces organised as a primitive patch.
Definition: faceZone.H:65
Foam::UPstream::lastSlave
static int lastSlave(const label communicator=0)
Process index of last slave.
Definition: UPstream.H:467
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
Foam::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
Foam::polyMesh::pointZones
const pointZoneMesh & pointZones() const
Return point zone mesh.
Definition: polyMesh.H:471
Foam::fvMeshTools::newMesh
static autoPtr< fvMesh > newMesh(const IOobject &io, const bool masterOnlyReading)
Read mesh or create dummy mesh (0 cells, >0 patches). Works in two.
Definition: fvMeshTools.C:424
Foam::DynamicList::append
DynamicList< T, SizeMin > & append(const T &val)
Append an element to the end of this list.
Definition: DynamicListI.H:472
Foam::IOobject::READ_IF_PRESENT
Definition: IOobject.H:122
Foam::PtrList::setSize
void setSize(const label newLen)
Same as resize()
Definition: PtrListI.H:108
Foam::andOp
Definition: ops.H:233
Foam::PtrList< entry >
Foam::fvMeshTools::reorderPatches
static void reorderPatches(fvMesh &, const labelList &oldToNew, const label nPatches, const bool validBoundary)
Reorder and remove trailing patches. If validBoundary call is parallel.
Definition: fvMeshTools.C:329
Foam::polyBoundaryMesh::findPatchID
label findPatchID(const word &patchName, const bool allowNotFound=true) const
Find patch index given a name, return -1 if not found.
Definition: polyBoundaryMesh.C:766
Foam::UPstream::commsTypes::scheduled
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
Foam::fvMesh::addFvPatches
void addFvPatches(PtrList< polyPatch > &plist, const bool validBoundary=true)
Add boundary patches. Constructor helper.
Definition: fvMesh.C:492
Foam::PtrList::clone
PtrList< T > clone(Args &&... args) const
Make a copy by cloning each of the list elements.
Foam::polyPatch::New
static autoPtr< polyPatch > New(const word &patchType, const word &name, const label size, const label start, const label index, const polyBoundaryMesh &bm)
Return a pointer to a new patch created on freestore from.
Definition: polyPatchNew.C:35
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:84
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:137
Foam::PtrList::set
const T * set(const label i) const
Return const pointer to element (if set) or nullptr.
Definition: PtrListI.H:143
Foam::polyPatch::start
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:311
Foam::Time::findInstance
word findInstance(const fileName &dir, const word &name=word::null, const IOobject::readOption rOpt=IOobject::MUST_READ, const word &stopInstance=word::null) const
Definition: Time.C:781
Foam::ZoneMesh::names
wordList names() const
A list of the zone names.
Definition: ZoneMesh.C:340
Foam::UPstream::master
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:438
Foam::UPstream::firstSlave
static constexpr int firstSlave() noexcept
Process index of first slave.
Definition: UPstream.H:461
Foam::New
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
Definition: DimensionedFieldReuseFunctions.H:105
Foam::fvMesh::boundary
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:589
Foam::autoPtr< Foam::fvMesh >
Foam::TimePaths::system
const word & system() const
Return system name.
Definition: TimePathsI.H:94
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:355
Foam::Time::path
fileName path() const
Return path.
Definition: Time.H:303
fvMeshTools.H
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::List< label >
Foam::IOobject::readOpt
readOption readOpt() const
The read option.
Definition: IOobjectI.H:141
Foam::fvMeshTools::addPatch
static label addPatch(fvMesh &mesh, const polyPatch &patch, const dictionary &patchFieldDict, const word &defaultPatchFieldType, const bool validBoundary)
Add patch. Inserts patch before all processor patches.
Definition: fvMeshTools.C:37
Foam::type
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:590
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:246
Foam::dictionary::add
entry * add(entry *entryPtr, bool mergeEntry=false)
Add a new entry.
Definition: dictionary.C:703
Foam::IPstream
Input inter-processor communications stream.
Definition: IPstream.H:52
Foam::List::setSize
void setSize(const label newSize)
Alias for resize(const label)
Definition: ListI.H:146
Foam::IOobject::NO_READ
Definition: IOobject.H:123
Foam::polyMesh::addZones
void addZones(const List< pointZone * > &pz, const List< faceZone * > &fz, const List< cellZone * > &cz)
Add mesh zones.
Definition: polyMesh.C:968
Foam::fvMesh::clearOut
void clearOut()
Clear all geometry and addressing.
Definition: fvMesh.C:234
Foam::isDir
bool isDir(const fileName &name, const bool followLink=true)
Does the name exist as a DIRECTORY in the file system?
Definition: MSwindows.C:643
Foam::IOobject::MUST_READ
Definition: IOobject.H:120