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 | Copyright (C) 2012-2015 OpenFOAM Foundation
6  \\/ M anipulation | Copyright (C) 2015 OpenCFD Ltd.
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 \*---------------------------------------------------------------------------*/
25 
26 #include "fvMeshTools.H"
29 #include "fvMeshTools.H"
30 
31 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32 
33 // Adds patch if not yet there. Returns patchID.
35 (
36  fvMesh& mesh,
37  const polyPatch& patch,
38  const dictionary& patchFieldDict,
39  const word& defaultPatchFieldType,
40  const bool validBoundary
41 )
42 {
43  polyBoundaryMesh& polyPatches =
44  const_cast<polyBoundaryMesh&>(mesh.boundaryMesh());
45 
46  label patchI = polyPatches.findPatchID(patch.name());
47  if (patchI != -1)
48  {
49  // Already there
50  return patchI;
51  }
52 
53 
54  // Append at end unless there are processor patches
55  label insertPatchI = polyPatches.size();
56  label startFaceI = mesh.nFaces();
57 
58  if (!isA<processorPolyPatch>(patch))
59  {
60  forAll(polyPatches, patchI)
61  {
62  const polyPatch& pp = polyPatches[patchI];
63 
64  if (isA<processorPolyPatch>(pp))
65  {
66  insertPatchI = patchI;
67  startFaceI = pp.start();
68  break;
69  }
70  }
71  }
72 
73 
74  // Below is all quite a hack. Feel free to change once there is a better
75  // mechanism to insert and reorder patches.
76 
77  // Clear local fields and e.g. polyMesh parallelInfo.
78  mesh.clearOut();
79 
80  label sz = polyPatches.size();
81 
82  fvBoundaryMesh& fvPatches = const_cast<fvBoundaryMesh&>(mesh.boundary());
83 
84  // Add polyPatch at the end
85  polyPatches.setSize(sz+1);
86  polyPatches.set
87  (
88  sz,
89  patch.clone
90  (
91  polyPatches,
92  insertPatchI, //index
93  0, //size
94  startFaceI //start
95  )
96  );
97  fvPatches.setSize(sz+1);
98  fvPatches.set
99  (
100  sz,
102  (
103  polyPatches[sz], // point to newly added polyPatch
104  mesh.boundary()
105  )
106  );
107 
108  addPatchFields<volScalarField>
109  (
110  mesh,
111  patchFieldDict,
112  defaultPatchFieldType,
114  );
115  addPatchFields<volVectorField>
116  (
117  mesh,
118  patchFieldDict,
119  defaultPatchFieldType,
121  );
122  addPatchFields<volSphericalTensorField>
123  (
124  mesh,
125  patchFieldDict,
126  defaultPatchFieldType,
128  );
129  addPatchFields<volSymmTensorField>
130  (
131  mesh,
132  patchFieldDict,
133  defaultPatchFieldType,
135  );
136  addPatchFields<volTensorField>
137  (
138  mesh,
139  patchFieldDict,
140  defaultPatchFieldType,
142  );
143 
144  // Surface fields
145 
146  addPatchFields<surfaceScalarField>
147  (
148  mesh,
149  patchFieldDict,
150  defaultPatchFieldType,
152  );
153  addPatchFields<surfaceVectorField>
154  (
155  mesh,
156  patchFieldDict,
157  defaultPatchFieldType,
159  );
160  addPatchFields<surfaceSphericalTensorField>
161  (
162  mesh,
163  patchFieldDict,
164  defaultPatchFieldType,
166  );
167  addPatchFields<surfaceSymmTensorField>
168  (
169  mesh,
170  patchFieldDict,
171  defaultPatchFieldType,
173  );
174  addPatchFields<surfaceTensorField>
175  (
176  mesh,
177  patchFieldDict,
178  defaultPatchFieldType,
180  );
181 
182  // Create reordering list
183  // patches before insert position stay as is
184  labelList oldToNew(sz+1);
185  for (label i = 0; i < insertPatchI; i++)
186  {
187  oldToNew[i] = i;
188  }
189  // patches after insert position move one up
190  for (label i = insertPatchI; i < sz; i++)
191  {
192  oldToNew[i] = i+1;
193  }
194  // appended patch gets moved to insert position
195  oldToNew[sz] = insertPatchI;
196 
197  // Shuffle into place
198  polyPatches.reorder(oldToNew, validBoundary);
199  fvPatches.reorder(oldToNew);
200 
201  reorderPatchFields<volScalarField>(mesh, oldToNew);
202  reorderPatchFields<volVectorField>(mesh, oldToNew);
203  reorderPatchFields<volSphericalTensorField>(mesh, oldToNew);
204  reorderPatchFields<volSymmTensorField>(mesh, oldToNew);
205  reorderPatchFields<volTensorField>(mesh, oldToNew);
206  reorderPatchFields<surfaceScalarField>(mesh, oldToNew);
207  reorderPatchFields<surfaceVectorField>(mesh, oldToNew);
208  reorderPatchFields<surfaceSphericalTensorField>(mesh, oldToNew);
209  reorderPatchFields<surfaceSymmTensorField>(mesh, oldToNew);
210  reorderPatchFields<surfaceTensorField>(mesh, oldToNew);
211 
212  return insertPatchI;
213 }
214 
215 
217 (
218  fvMesh& mesh,
219  const label patchI,
220  const dictionary& patchFieldDict
221 )
222 {
223  setPatchFields<volScalarField>(mesh, patchI, patchFieldDict);
224  setPatchFields<volVectorField>(mesh, patchI, patchFieldDict);
225  setPatchFields<volSphericalTensorField>(mesh, patchI, patchFieldDict);
226  setPatchFields<volSymmTensorField>(mesh, patchI, patchFieldDict);
227  setPatchFields<volTensorField>(mesh, patchI, patchFieldDict);
228  setPatchFields<surfaceScalarField>(mesh, patchI, patchFieldDict);
229  setPatchFields<surfaceVectorField>(mesh, patchI, patchFieldDict);
230  setPatchFields<surfaceSphericalTensorField>
231  (
232  mesh,
233  patchI,
234  patchFieldDict
235  );
236  setPatchFields<surfaceSymmTensorField>(mesh, patchI, patchFieldDict);
237  setPatchFields<surfaceTensorField>(mesh, patchI, patchFieldDict);
238 }
239 
240 
242 {
243  setPatchFields<volScalarField>(mesh, patchI, pTraits<scalar>::zero);
244  setPatchFields<volVectorField>(mesh, patchI, pTraits<vector>::zero);
245  setPatchFields<volSphericalTensorField>
246  (
247  mesh,
248  patchI,
250  );
251  setPatchFields<volSymmTensorField>
252  (
253  mesh,
254  patchI,
256  );
257  setPatchFields<volTensorField>(mesh, patchI, pTraits<tensor>::zero);
258  setPatchFields<surfaceScalarField>(mesh, patchI, pTraits<scalar>::zero);
259  setPatchFields<surfaceVectorField>(mesh, patchI, pTraits<vector>::zero);
260  setPatchFields<surfaceSphericalTensorField>
261  (
262  mesh,
263  patchI,
265  );
266  setPatchFields<surfaceSymmTensorField>
267  (
268  mesh,
269  patchI,
271  );
272  setPatchFields<surfaceTensorField>(mesh, patchI, pTraits<tensor>::zero);
273 }
274 
275 
276 // Deletes last patch
278 {
279  // Clear local fields and e.g. polyMesh globalMeshData.
280  mesh.clearOut();
281 
282  polyBoundaryMesh& polyPatches =
283  const_cast<polyBoundaryMesh&>(mesh.boundaryMesh());
284  fvBoundaryMesh& fvPatches = const_cast<fvBoundaryMesh&>(mesh.boundary());
285 
286  if (polyPatches.empty())
287  {
289  << "No patches in mesh"
290  << abort(FatalError);
291  }
292 
293  label nFaces = 0;
294  for (label patchI = nPatches; patchI < polyPatches.size(); patchI++)
295  {
296  nFaces += polyPatches[patchI].size();
297  }
298  reduce(nFaces, sumOp<label>());
299 
300  if (nFaces)
301  {
303  << "There are still " << nFaces
304  << " faces in " << polyPatches.size()-nPatches
305  << " patches to be deleted" << abort(FatalError);
306  }
307 
308  // Remove actual patches
309  polyPatches.setSize(nPatches);
310  fvPatches.setSize(nPatches);
311 
312  trimPatchFields<volScalarField>(mesh, nPatches);
313  trimPatchFields<volVectorField>(mesh, nPatches);
314  trimPatchFields<volSphericalTensorField>(mesh, nPatches);
315  trimPatchFields<volSymmTensorField>(mesh, nPatches);
316  trimPatchFields<volTensorField>(mesh, nPatches);
317 
318  trimPatchFields<surfaceScalarField>(mesh, nPatches);
319  trimPatchFields<surfaceVectorField>(mesh, nPatches);
320  trimPatchFields<surfaceSphericalTensorField>(mesh, nPatches);
321  trimPatchFields<surfaceSymmTensorField>(mesh, nPatches);
322  trimPatchFields<surfaceTensorField>(mesh, nPatches);
323 }
324 
325 
327 (
328  fvMesh& mesh,
329  const labelList& oldToNew,
330  const label nNewPatches,
331  const bool validBoundary
332 )
333 {
334  polyBoundaryMesh& polyPatches =
335  const_cast<polyBoundaryMesh&>(mesh.boundaryMesh());
336  fvBoundaryMesh& fvPatches = const_cast<fvBoundaryMesh&>(mesh.boundary());
337 
338  // Shuffle into place
339  polyPatches.reorder(oldToNew, validBoundary);
340  fvPatches.reorder(oldToNew);
341 
342  reorderPatchFields<volScalarField>(mesh, oldToNew);
343  reorderPatchFields<volVectorField>(mesh, oldToNew);
344  reorderPatchFields<volSphericalTensorField>(mesh, oldToNew);
345  reorderPatchFields<volSymmTensorField>(mesh, oldToNew);
346  reorderPatchFields<volTensorField>(mesh, oldToNew);
347  reorderPatchFields<surfaceScalarField>(mesh, oldToNew);
348  reorderPatchFields<surfaceVectorField>(mesh, oldToNew);
349  reorderPatchFields<surfaceSphericalTensorField>(mesh, oldToNew);
350  reorderPatchFields<surfaceSymmTensorField>(mesh, oldToNew);
351  reorderPatchFields<surfaceTensorField>(mesh, oldToNew);
352 
353  // Remove last.
354  trimPatches(mesh, nNewPatches);
355 }
356 
357 
359 (
360  fvMesh& mesh,
361  const bool validBoundary
362 )
363 {
364  const polyBoundaryMesh& pbm = mesh.boundaryMesh();
365 
366  labelList newToOld(pbm.size());
367  labelList oldToNew(pbm.size(), -1);
368  label newI = 0;
369 
370 
371  // Assumes all non-coupled boundaries are on all processors!
372  forAll(pbm, patchI)
373  {
374  const polyPatch& pp = pbm[patchI];
375 
376  if (!isA<processorPolyPatch>(pp))
377  {
378  label nFaces = pp.size();
379  if (validBoundary)
380  {
381  reduce(nFaces, sumOp<label>());
382  }
383 
384  if (nFaces > 0)
385  {
386  newToOld[newI] = patchI;
387  oldToNew[patchI] = newI++;
388  }
389  }
390  }
391 
392  // Same for processor patches (but need no reduction)
393  forAll(pbm, patchI)
394  {
395  const polyPatch& pp = pbm[patchI];
396 
397  if (isA<processorPolyPatch>(pp) && pp.size())
398  {
399  newToOld[newI] = patchI;
400  oldToNew[patchI] = newI++;
401  }
402  }
403 
404  newToOld.setSize(newI);
405 
406  // Move all deleteable patches to the end
407  forAll(oldToNew, patchI)
408  {
409  if (oldToNew[patchI] == -1)
410  {
411  oldToNew[patchI] = newI++;
412  }
413  }
414 
415  reorderPatches(mesh, oldToNew, newToOld.size(), validBoundary);
416 
417  return newToOld;
418 }
419 
420 
422 (
423  const IOobject& io,
424  const bool masterOnlyReading
425 )
426 {
427  // Region name
428  // ~~~~~~~~~~~
429 
430  fileName meshSubDir;
431 
432  if (io.name() == polyMesh::defaultRegion)
433  {
434  meshSubDir = polyMesh::meshSubDir;
435  }
436  else
437  {
438  meshSubDir = io.name()/polyMesh::meshSubDir;
439  }
440 
441 
442  fileName facesInstance;
443 
444  // Patch types
445  // ~~~~~~~~~~~
446  // Read and scatter master patches (without reading master mesh!)
447 
448  PtrList<entry> patchEntries;
449  if (Pstream::master())
450  {
451  facesInstance = io.time().findInstance
452  (
453  meshSubDir,
454  "faces",
456  );
457 
458  patchEntries = polyBoundaryMeshEntries
459  (
460  IOobject
461  (
462  "boundary",
463  facesInstance,
464  meshSubDir,
465  io.db(),
468  false
469  )
470  );
471 
472  // Send patches
473  for
474  (
475  int slave=Pstream::firstSlave();
476  slave<=Pstream::lastSlave();
477  slave++
478  )
479  {
480  OPstream toSlave(Pstream::scheduled, slave);
481  toSlave << patchEntries;
482  }
483  }
484  else
485  {
486  // Receive patches
488  fromMaster >> patchEntries;
489  }
490 
491 
492  Pstream::scatter(facesInstance);
493 
494  // Dummy meshes
495  // ~~~~~~~~~~~~
496 
497  // Check who has a mesh
498  const fileName meshDir = io.time().path()/facesInstance/meshSubDir;
499  bool haveMesh = isDir(meshDir);
500  if (masterOnlyReading && !Pstream::master())
501  {
502  haveMesh = false;
503  }
504 
505 
506  // Set up to read-if-present. Note: does not search for mesh so set
507  // instance explicitly
508  IOobject meshIO(io);
509  meshIO.instance() = facesInstance;
510  if (masterOnlyReading && !Pstream::master())
511  {
512  meshIO.readOpt() = IOobject::NO_READ;
513  }
514  else
515  {
517  }
518 
519 
520  // Read mesh
521  // ~~~~~~~~~
522  // Now all processors read a mesh and use supplied points,faces etc
523  // if there is none.
524  // Note: fvSolution, fvSchemes are also using the supplied Ioobject so
525  // on slave will be NO_READ, on master READ_IF_PRESENT. This will
526  // conflict with e.g. timeStampMaster reading so switch off.
527 
528  const regIOobject::fileCheckTypes oldCheckType =
531 
533  (
534  new fvMesh
535  (
536  meshIO,
537  xferCopy(pointField()),
538  xferCopy(faceList()),
539  xferCopy(labelList()),
541  )
542  );
543  fvMesh& mesh = meshPtr();
544 
546 
547 
548 
549  // Add patches
550  // ~~~~~~~~~~~
551 
552 
553  bool havePatches = mesh.boundary().size();
554  reduce(havePatches, andOp<bool>());
555 
556  if (!havePatches)
557  {
558  // There are one or more processors without full complement of
559  // patches
560 
561  DynamicList<polyPatch*> newPatches;
562 
563  if (haveMesh) //Pstream::master())
564  {
565  forAll(mesh.boundary(), patchI)
566  {
567  newPatches.append
568  (
569  mesh.boundaryMesh()[patchI].clone(mesh.boundaryMesh()).ptr()
570  );
571  }
572  }
573  else
574  {
575  forAll(patchEntries, patchI)
576  {
577  const entry& e = patchEntries[patchI];
578  const word type(e.dict().lookup("type"));
579 
580  if
581  (
582  type == processorPolyPatch::typeName
583  || type == processorCyclicPolyPatch::typeName
584  )
585  {}
586  else
587  {
588  dictionary patchDict(e.dict());
589  patchDict.set("nFaces", 0);
590  patchDict.set("startFace", 0);
591 
592  newPatches.append
593  (
595  (
596  patchEntries[patchI].keyword(),
597  patchDict,
598  newPatches.size(),
600  ).ptr()
601  );
602  }
603  }
604  }
606  mesh.addFvPatches(newPatches);
607  }
608 
609  //Pout<< "patches:" << endl;
610  //forAll(mesh.boundary(), patchI)
611  //{
612  // Pout<< " type" << mesh.boundary()[patchI].type()
613  // << " size:" << mesh.boundary()[patchI].size()
614  // << " start:" << mesh.boundary()[patchI].start() << endl;
615  //}
616  //Pout<< endl;
617 
618 
619  // Determine zones
620  // ~~~~~~~~~~~~~~~
621 
622  wordList pointZoneNames(mesh.pointZones().names());
623  Pstream::scatter(pointZoneNames);
624  wordList faceZoneNames(mesh.faceZones().names());
625  Pstream::scatter(faceZoneNames);
626  wordList cellZoneNames(mesh.cellZones().names());
627  Pstream::scatter(cellZoneNames);
628 
629  if (!haveMesh)
630  {
631  // Add the zones. Make sure to remove the old dummy ones first
632  mesh.pointZones().clear();
633  mesh.faceZones().clear();
634  mesh.cellZones().clear();
635 
636  List<pointZone*> pz(pointZoneNames.size());
637  forAll(pointZoneNames, i)
638  {
639  pz[i] = new pointZone
640  (
641  pointZoneNames[i],
642  labelList(0),
643  i,
644  mesh.pointZones()
645  );
646  }
647  List<faceZone*> fz(faceZoneNames.size());
648  forAll(faceZoneNames, i)
649  {
650  fz[i] = new faceZone
651  (
652  faceZoneNames[i],
653  labelList(0),
654  boolList(0),
655  i,
656  mesh.faceZones()
657  );
658  }
659  List<cellZone*> cz(cellZoneNames.size());
660  forAll(cellZoneNames, i)
661  {
662  cz[i] = new cellZone
663  (
664  cellZoneNames[i],
665  labelList(0),
666  i,
667  mesh.cellZones()
668  );
669  }
670 
671  if (pz.size() && fz.size() && cz.size())
672  {
673  mesh.addZones(pz, fz, cz);
674  }
675  }
676 
677  return meshPtr;
678 }
679 
680 
681 // ************************************************************************* //
Foam::entry
A keyword and a list of tokens is an 'entry'.
Definition: entry.H:65
Foam::pointField
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
Foam::ZoneMesh::clear
void clear()
Clear the zones.
Definition: ZoneMesh.C:408
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Foam::polyBoundaryMesh::findPatchID
label findPatchID(const word &patchName) const
Find patch index given a name.
Definition: polyBoundaryMesh.C:678
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::compressible::New
autoPtr< BasicCompressibleTurbulenceModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleTurbulenceModel::transportModel &transport, const word &propertiesName)
Definition: turbulentFluidThermoModel.C:36
Foam::fileName
A class for handling file names.
Definition: fileName.H:69
Foam::labelList
List< label > labelList
A List of labels.
Definition: labelList.H:56
Foam::polyBoundaryMeshEntries
Foam::polyBoundaryMeshEntries.
Definition: polyBoundaryMeshEntries.H:48
Foam::polyBoundaryMesh
Foam::polyBoundaryMesh.
Definition: polyBoundaryMesh.H:60
nPatches
label nPatches
Definition: readKivaGrid.H:402
Foam::polyMesh::defaultRegion
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:306
Foam::UPstream::scheduled
@ scheduled
Definition: UPstream.H:67
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::DynamicList
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:56
Foam::OPstream
Output inter-processor communications stream.
Definition: OPstream.H:50
Foam::fvMesh::addFvPatches
void addFvPatches(const List< polyPatch * > &, const bool validBoundary=true)
Add boundary patches. Constructor helper.
Definition: fvMesh.C:462
Foam::pointZone
A subset of mesh points. The labels of points in the zone can be obtained from the addressing() list.
Definition: pointZone.H:62
Foam::fvMeshTools::removeEmptyPatches
static labelList removeEmptyPatches(fvMesh &, const bool validBoundary)
Remove zero sized patches. All but processor patches are.
Definition: fvMeshTools.C:359
Foam::polyMesh::meshSubDir
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
Definition: polyMesh.H:309
Foam::polyMesh::cellZones
const cellZoneMesh & cellZones() const
Return cell zone mesh.
Definition: polyMesh.H:469
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:147
Foam::fvMeshTools::zeroPatchFields
static void zeroPatchFields(fvMesh &mesh, const label patchI)
Change patchField to zero on registered fields.
Definition: fvMeshTools.C:241
Foam::fvMesh::removeFvBoundary
void removeFvBoundary()
Remove boundary patches. Warning: fvPatchFields hold ref to.
Definition: fvMesh.C:480
Foam::fvMeshTools::trimPatches
static void trimPatches(fvMesh &, const label nPatches)
Definition: fvMeshTools.C:277
Foam::IOobject::MUST_READ
@ MUST_READ
Definition: IOobject.H:108
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:421
Foam::polyBoundaryMesh::reorder
void reorder(const labelUList &, const bool validBoundary)
Reorders patches. Ordering does not have to be done in.
Definition: polyBoundaryMesh.C:1100
polyBoundaryMeshEntries.H
Foam::IOobject::instance
const fileName & instance() const
Definition: IOobject.H:350
Foam::IOobject::time
const Time & time() const
Return time.
Definition: IOobject.C:245
processorCyclicPolyPatch.H
Foam::IOobject::NO_WRITE
@ NO_WRITE
Definition: IOobject.H:118
Foam::cellZone
A subset of mesh cells.
Definition: cellZone.H:61
Foam::polyMesh::faceZones
const faceZoneMesh & faceZones() const
Return face zone mesh.
Definition: polyMesh.H:463
Foam::IOobject::db
const objectRegistry & db() const
Return the local objectRegistry.
Definition: IOobject.C:239
Foam::fvBoundaryMesh
Foam::fvBoundaryMesh.
Definition: fvBoundaryMesh.H:52
Foam::regIOobject::timeStamp
@ timeStamp
Definition: regIOobject.H:70
Foam::reduce
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
Definition: PstreamReduceOps.H:43
Foam::PtrList::set
bool set(const label) const
Is element set.
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:59
Foam::IOobject::readOpt
readOption readOpt() const
Definition: IOobject.H:317
Foam::faceZone
A subset of mesh faces organised as a primitive patch.
Definition: faceZone.H:64
Foam::UPstream::lastSlave
static int lastSlave(const label communicator=0)
Process index of last slave.
Definition: UPstream.H:428
Foam::IOobject::NO_READ
@ NO_READ
Definition: IOobject.H:111
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:457
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:422
Foam::PtrList::empty
bool empty() const
Return true if the PtrList is empty (ie, size() is zero).
Definition: PtrListI.H:39
Foam::IOobject::name
const word & name() const
Return name.
Definition: IOobject.H:273
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:327
Foam::UPstream::masterNo
static int masterNo()
Process index of the master.
Definition: UPstream.H:393
Foam::regIOobject::fileCheckTypes
fileCheckTypes
Types of communications.
Definition: regIOobject.H:68
meshPtr
static fvMesh * meshPtr
Definition: globalFoam.H:52
Foam::PtrList::reorder
void reorder(const labelUList &)
Reorders elements. Ordering does not have to be done in.
Definition: PtrList.C:208
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
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:32
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
Foam::UPstream::firstSlave
static int firstSlave()
Process index of first slave.
Definition: UPstream.H:422
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:131
Foam::isDir
bool isDir(const fileName &)
Does the name exist as a DIRECTORY in the file system?
Definition: POSIX.C:615
Foam::e
const double e
Elementary charge.
Definition: doubleFloat.H:94
Foam::boolList
List< bool > boolList
Bool container classes.
Definition: boolList.H:50
Foam::polyPatch::start
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:312
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
Return the location of "dir" containing the file "name".
Definition: findInstance.C:38
Foam::ZoneMesh::names
wordList names() const
Return a list of zone names.
Definition: ZoneMesh.C:263
Foam::DynamicList::append
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
Foam::polyPatch::clone
virtual autoPtr< polyPatch > clone(const polyBoundaryMesh &bm) const
Construct and return a clone, resetting the boundary mesh.
Definition: polyPatch.H:231
Foam::UPstream::master
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:399
Foam::regIOobject::fileModificationChecking
static fileCheckTypes fileModificationChecking
Definition: regIOobject.H:128
Foam::fvMesh::boundary
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:550
Foam::List::setSize
void setSize(const label)
Reset size of List.
Foam::autoPtr< Foam::fvMesh >
Foam::andOp
Definition: ops.H:177
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
Foam::sumOp
Definition: ops.H:162
Foam::Time::path
fileName path() const
Return path.
Definition: Time.H:281
fvMeshTools.H
Foam::faceList
List< face > faceList
Definition: faceListFwd.H:43
Foam::xferCopy
Xfer< T > xferCopy(const T &)
Construct by copying the contents of the arg.
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
Foam::pTraits
Traits class for primitives.
Definition: pTraits.H:50
Foam::IOobject::clone
Foam::autoPtr< IOobject > clone() const
Clone.
Definition: IOobject.H:252
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:35
Foam::fvMeshTools::setPatchFields
static void setPatchFields(fvMesh &mesh, const label patchI, const dictionary &patchFieldDict)
Set patchFields according to dictionary.
Definition: fvMeshToolsTemplates.C:89
Foam::PtrList::size
label size() const
Return the number of elements in the PtrList.
Definition: PtrListI.H:32
Foam::PtrList::setSize
void setSize(const label)
Reset size of PtrList. If extending the PtrList, new entries are.
Definition: PtrList.C:142
Foam::IPstream
Input inter-processor communications stream.
Definition: IPstream.H:50
Foam::List::size
void size(const label)
Override size to be inconsistent with allocated storage.
Foam::type
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:588
Foam::IOobject::READ_IF_PRESENT
@ READ_IF_PRESENT
Definition: IOobject.H:110
Foam::patchIdentifier::name
const word & name() const
Return name.
Definition: patchIdentifier.H:109
Foam::polyMesh::addZones
void addZones(const List< pointZone * > &pz, const List< faceZone * > &fz, const List< cellZone * > &cz)
Add mesh zones.
Definition: polyMesh.C:922
Foam::fvMesh::clearOut
void clearOut()
Clear all geometry and addressing.
Definition: fvMesh.C:231
Foam::dictionary::set
void set(entry *)
Assign a new entry, overwrite any existing entry.
Definition: dictionary.C:856