All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
createBaffles.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-2022 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  createBaffles
26 
27 Description
28  Makes faces into boundary faces. Does not duplicate points.
29 
30  Notes:
31  - If any coupled patch face is selected for baffling the opposite member
32  has to be selected for baffling as well.
33  - If fields are being modified then boundary conditions must be specified
34  within the patch's patchFields sub-dictionary.
35  - Any patches left with no faces will be removed, except for patches of
36  coupled, internal, of non-conformal type
37 
38 \*---------------------------------------------------------------------------*/
39 
40 #include "argList.H"
41 #include "Time.H"
42 #include "polyTopoChange.H"
43 #include "polyModifyFace.H"
44 #include "polyAddFace.H"
45 #include "ReadFields.H"
46 #include "volFields.H"
47 #include "surfaceFields.H"
48 #include "pointFields.H"
49 #include "fvMeshMapper.H"
50 #include "faceSelection.H"
51 #include "searchableSurface.H"
52 #include "fvMeshTools.H"
53 #include "systemDict.H"
54 #include "processorPolyPatch.H"
55 #include "internalPolyPatch.H"
56 #include "nonConformalPolyPatch.H"
57 
58 using namespace Foam;
59 
60 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
61 
62 void filterPatches(fvMesh& mesh, const HashSet<word>& bafflePatches)
63 {
64  // Remove any zero-sized patches, except for constraint types, and any
65  // specified in the system dictionary
66 
67  const polyBoundaryMesh& bMesh = mesh.boundaryMesh();
68 
69  label newPatchi = 0;
70 
71  // List of patch's new indices
72  labelList oldToNew(bMesh.size(), -1);
73 
74  // Add all the kept non-processor patches to the list first
75  forAll(bMesh, patchi)
76  {
77  const polyPatch& pp = bMesh[patchi];
78 
79  if (!isA<processorPolyPatch>(pp))
80  {
81  if
82  (
83  bafflePatches.found(pp.name())
84  || fvPatch::constraintType(pp.type())
85  || returnReduce(pp.size(), sumOp<label>())
86  )
87  {
88  oldToNew[patchi] = newPatchi++;
89  }
90  }
91  }
92 
93  // Now add all the processor patches
94  forAll(bMesh, patchi)
95  {
96  const polyPatch& pp = bMesh[patchi];
97 
98  if (isA<processorPolyPatch>(pp))
99  {
100  oldToNew[patchi] = newPatchi++;
101  }
102  }
103 
104  // Note how many patches are to be kept
105  const label nKeepPatches = newPatchi;
106 
107  // If there are any to remove ...
108  if (nKeepPatches != bMesh.size())
109  {
110  Info<< endl << "Removing zero-sized patches:" << endl << incrIndent;
111 
112  // Add all the removed patches to the list
113  forAll(oldToNew, patchi)
114  {
115  if (oldToNew[patchi] == -1)
116  {
117  Info<< indent << bMesh[patchi].name()
118  << " type " << bMesh[patchi].type()
119  << " at position " << patchi << endl;
120 
121  oldToNew[patchi] = newPatchi++;
122  }
123  }
124 
125  Info<< decrIndent;
126 
127  // Permute the mesh, keeping only the ones not to be removed
128  fvMeshTools::reorderPatches(mesh, oldToNew, nKeepPatches, true);
129 
130  Info<< endl;
131  }
132 }
133 
134 
135 void modifyOrAddFace
136 (
137  polyTopoChange& meshMod,
138  const face& f,
139  const label facei,
140  const label own,
141  const bool flipFaceFlux,
142  const label newPatchi,
143  const label zoneID,
144  const bool zoneFlip,
145  PackedBoolList& modifiedFace
146 )
147 {
148  if (!modifiedFace[facei])
149  {
150  // First usage of face. Modify.
151  meshMod.setAction
152  (
154  (
155  f, // modified face
156  facei, // label of face
157  own, // owner
158  -1, // neighbour
159  flipFaceFlux, // face flip
160  newPatchi, // patch for face
161  false, // remove from zone
162  zoneID, // zone for face
163  zoneFlip // face flip in zone
164  )
165  );
166 
167  modifiedFace[facei] = 1;
168  }
169  else
170  {
171  // Second usage of face. Add.
172  meshMod.setAction
173  (
175  (
176  f, // modified face
177  own, // owner
178  -1, // neighbour
179  -1, // master point
180  -1, // master edge
181  facei, // master face
182  flipFaceFlux, // face flip
183  newPatchi, // patch for face
184  zoneID, // zone for face
185  zoneFlip // face flip in zone
186  )
187  );
188  }
189 }
190 
191 
192 label createFaces
193 (
194  const bool internalFacesOnly,
195  const fvMesh& mesh,
196  const faceZone& fZone,
197  const label newOwnerPatchi,
198  const label newNeighbourPatchi,
199  polyTopoChange& meshMod,
200  PackedBoolList& modifiedFace
201 )
202 {
203  label nModified = 0;
204 
205  const polyBoundaryMesh& bMesh = mesh.boundaryMesh();
206 
207  // Pass 1. Do selected side of zone
208  for (label facei = 0; facei < mesh.nInternalFaces(); facei++)
209  {
210  const label zoneFacei = fZone.whichFace(facei);
211 
212  if (zoneFacei != -1)
213  {
214  if (!fZone.flipMap()[zoneFacei])
215  {
216  // Use owner side of face
217  modifyOrAddFace
218  (
219  meshMod,
220  mesh.faces()[facei], // modified face
221  facei, // label of face
222  mesh.faceOwner()[facei],// owner
223  false, // face flip
224  newOwnerPatchi, // patch for face
225  fZone.index(), // zone for face
226  false, // face flip in zone
227  modifiedFace // modify or add status
228  );
229  }
230  else
231  {
232  // Use neighbour side of face.
233  // To keep faceZone pointing out of original neighbour
234  // we don't need to set faceFlip since that cell
235  // now becomes the owner
236  modifyOrAddFace
237  (
238  meshMod,
239  mesh.faces()[facei].reverseFace(), // modified face
240  facei, // label of face
241  mesh.faceNeighbour()[facei],// owner
242  true, // face flip
243  newOwnerPatchi, // patch for face
244  fZone.index(), // zone for face
245  false, // face flip in zone
246  modifiedFace // modify or add status
247  );
248  }
249 
250  nModified++;
251  }
252  }
253 
254  // Pass 2. Do other side of zone
255  for (label facei = 0; facei < mesh.nInternalFaces(); facei++)
256  {
257  label zoneFacei = fZone.whichFace(facei);
258 
259  if (zoneFacei != -1)
260  {
261  if (!fZone.flipMap()[zoneFacei])
262  {
263  // Use neighbour side of face
264  modifyOrAddFace
265  (
266  meshMod,
267  mesh.faces()[facei].reverseFace(), // modified face
268  facei, // label of face
269  mesh.faceNeighbour()[facei], // owner
270  true, // face flip
271  newNeighbourPatchi, // patch for face
272  fZone.index(), // zone for face
273  true, // face flip in zone
274  modifiedFace // modify or add
275  );
276  }
277  else
278  {
279  // Use owner side of face
280  modifyOrAddFace
281  (
282  meshMod,
283  mesh.faces()[facei], // modified face
284  facei, // label of face
285  mesh.faceOwner()[facei],// owner
286  false, // face flip
287  newNeighbourPatchi, // patch for face
288  fZone.index(), // zone for face
289  true, // face flip in zone
290  modifiedFace // modify or add status
291  );
292  }
293  }
294  }
295 
296  // Modify any boundary faces
297  forAll(bMesh, patchi)
298  {
299  const polyPatch& pp = bMesh[patchi];
300 
301  if (pp.coupled() || !internalFacesOnly)
302  {
303  forAll(pp, i)
304  {
305  const label facei = pp.start() + i;
306  const label zoneFacei = fZone.whichFace(facei);
307 
308  if (zoneFacei != -1)
309  {
310  const polyPatch& newPp =
311  fZone.flipMap()[zoneFacei]
312  ? bMesh[newOwnerPatchi]
313  : bMesh[newNeighbourPatchi];
314 
315  // We cannot move coupled faces to different coupled
316  // faces. Generate an error if this is attempted.
317  if (pp.coupled() && newPp.coupled())
318  {
320  << "Face on coupled patch \"" << pp.name()
321  << "\" selected for conversion to coupled "
322  << "patch \"" << newPp.name() << "\". "
323  << "Conversion from coupled patch to coupled "
324  << "patch is not allowed."
325  << exit(FatalError);
326  }
327 
328  modifyOrAddFace
329  (
330  meshMod,
331  mesh.faces()[facei], // modified face
332  facei, // label of face
333  mesh.faceOwner()[facei], // owner
334  false, // face flip
335  newPp.index(), // patch for face
336  fZone.index(), // zone for face
337  fZone.flipMap()[zoneFacei], // face flip in zone
338  modifiedFace // modify or add
339  );
340 
341  nModified++;
342  }
343  }
344  }
345  }
346 
347  return returnReduce(nModified, sumOp<label>());
348 }
349 
350 
351 int main(int argc, char *argv[])
352 {
354  (
355  "Makes internal faces into boundary faces.\n"
356  "Does not duplicate points."
357  );
358  #include "addDictOption.H"
359  #include "addOverwriteOption.H"
360  #include "addDictOption.H"
361  #include "addRegionOption.H"
362 
363  #include "setRootCase.H"
364  #include "createTime.H"
365  runTime.functionObjects().off();
366  #include "createNamedMesh.H"
367 
368  const polyBoundaryMesh& bMesh = mesh.boundaryMesh();
369 
370  const bool overwrite = args.optionFound("overwrite");
371 
372  const word oldInstance = mesh.pointsInstance();
373 
374  // Read the system dictionary and global controls
375  const dictionary dict(systemDict("createBafflesDict", args, mesh));
376  Switch internalFacesOnly = dict.lookup<Switch>("internalFacesOnly");
377  const Switch fields = dict.lookupOrDefault("fields", false);
378 
379  // Create face selections
380  PtrList<faceSelection> selectors;
381  {
382  const dictionary& selectionsDict = dict.subDict("baffles");
383 
384  label n = 0;
385  forAllConstIter(dictionary, selectionsDict, iter)
386  {
387  if (iter().isDict())
388  {
389  n++;
390  }
391  }
392  selectors.setSize(n);
393 
394  n = 0;
395  forAllConstIter(dictionary, selectionsDict, iter)
396  {
397  if (iter().isDict())
398  {
399  selectors.set
400  (
401  n++,
402  faceSelection::New(iter().keyword(), mesh, iter().dict())
403  );
404  }
405  }
406  }
407 
408  if (internalFacesOnly)
409  {
410  Info<< "Not converting faces on non-coupled patches." << nl << endl;
411  }
412 
413  // Read fields
414  IOobjectList objects(mesh, runTime.timeName());
415  if (fields) Info<< "Reading geometric fields" << nl << endl;
416  #include "readVolFields.H"
417  #include "readSurfaceFields.H"
418  #include "readPointFields.H"
419  if (fields) Info<< endl;
420 
421  // Creating faceZones for selectors that are not already faceZones
422  forAll(selectors, selectorI)
423  {
424  const word& name = selectors[selectorI].name();
425 
426  if (mesh.faceZones().findZoneID(name) == -1)
427  {
428  mesh.faceZones().clearAddressing();
429 
430  const label sz = mesh.faceZones().size();
431 
432  labelList addr(0);
433  boolList flip(0);
434 
435  mesh.faceZones().setSize(sz+1);
436  mesh.faceZones().set
437  (
438  sz,
439  new faceZone(name, addr, flip, sz, mesh.faceZones())
440  );
441  }
442  }
443 
444  // Per face, its zone ID and flip status
445  labelList faceToZoneID(mesh.nFaces(), -1);
446  boolList faceToFlip(mesh.nFaces(), false);
447  forAll(selectors, selectorI)
448  {
449  const word& name = selectors[selectorI].name();
450  label zoneID = mesh.faceZones().findZoneID(name);
451 
452  selectors[selectorI].select(zoneID, faceToZoneID, faceToFlip);
453  }
454 
455  // Add faces to faceZones
456  labelList nFaces(mesh.faceZones().size(), 0);
457  forAll(faceToZoneID, facei)
458  {
459  label zoneID = faceToZoneID[facei];
460  if (zoneID != -1)
461  {
462  nFaces[zoneID]++;
463  }
464  }
465  forAll(selectors, selectorI)
466  {
467  const word& name = selectors[selectorI].name();
468  const label zoneID = mesh.faceZones().findZoneID(name);
469 
470  label& n = nFaces[zoneID];
471  labelList addr(n);
472  boolList flip(n);
473 
474  n = 0;
475  forAll(faceToZoneID, facei)
476  {
477  label zone = faceToZoneID[facei];
478  if (zone == zoneID)
479  {
480  addr[n] = facei;
481  flip[n] = faceToFlip[facei];
482  n++;
483  }
484  }
485 
486  Info<< "Created zone '" << name << "' with "
487  << returnReduce(n, sumOp<label>()) << " faces" << endl;
488 
489  mesh.faceZones().set
490  (
491  zoneID,
492  new faceZone(name, addr, flip, zoneID, mesh.faceZones())
493  );
494  }
495 
496  Info<< endl;
497 
498  // Keywords for the two patch entries in each selector dictionary. Note
499  // "owner" and "neighbour" are preferred; "master" and "slave" are still
500  // permitted for backwards compatibility.
501  const FixedList<wordList, 2> patchEntryKeywords
502  ({
503  wordList({"owner", "master"}),
504  wordList({"neighbour", "slave"})
505  });
506 
507  // Create the baffles. Notes:
508  //
509  // - This is done in multiple steps
510  // - Patches are created with 'calculated' patchFields
511  // - Faces are moved into these patches
512  // - The patchFields are changed to the desired type
513  // - This helps resolve ordering issues
514  // - patchFields cannot be created at the same time as patches since a
515  // coupled patchField constructor frequently needs access to the
516  // neighbouring patch
517  // - patchFields need to be created after the faces have been added to
518  // the patch so that they don't later have to be mapped from nothing
519  //
520 
521  // Add patches
522  HashSet<word> bafflePatches;
523  forAll(selectors, selectorI)
524  {
525  const word& groupName = selectors[selectorI].name();
526  const dictionary& dict =
527  selectors[selectorI].dict().optionalSubDict("patches");
528 
529  forAll(patchEntryKeywords, i)
530  {
531  dictionary patchDict =
533  (
534  patchEntryKeywords[i],
535  false,
536  false
537  ).dict();
538 
539  const word patchName(patchDict.lookup<word>("name"));
540 
541  if (bMesh.findPatchID(patchName) == -1)
542  {
543  Info<< "Adding patch '" << patchName << "' to the mesh" << endl;
544 
545  // Create the patch
546  patchDict.set("nFaces", 0);
547  patchDict.set("startFace", 0);
548  autoPtr<polyPatch> ppPtr
549  (
551  (
552  patchName,
553  patchDict,
554  0,
555  bMesh
556  )
557  );
558  polyPatch& pp = ppPtr();
559 
560  // Add it to the group
561  if (!groupName.empty() && !pp.inGroup(groupName))
562  {
563  pp.inGroups().append(groupName);
564  }
565 
566  // Add the patch to the mesh, creating calculated boundary
567  // conditions everywhere. These will be changed later.
569  (
570  mesh,
571  pp,
572  dictionary(), // do not set specialised patchFields
574  true // parallel sync'ed addition
575  );
576  }
577  else
578  {
579  Info<< "Patch '" << patchName
580  << "' already exists in the mesh" << endl;
581 
582  patchDict.remove("name");
583  patchDict.remove("patchFields");
584 
585  if (!patchDict.empty())
586  {
588  << "Patch settings found in " << patchDict.name()
589  << " for patch '" << patchName << "', but this patch "
590  << "already exists so these settings will not be used"
591  << endl;
592  }
593  }
594 
595  bafflePatches.insert(patchName);
596  }
597  }
598 
599  Info<< endl;
600 
601  // Make sure patches and zoneFaces are synchronised across couples
602  mesh.boundaryMesh().checkParallelSync(true);
603  mesh.faceZones().checkParallelSync(true);
604 
605  // Do the topology changes. Notes:
606  //
607  // - Loop in incrementing face order (not necessary if faceZone ordered).
608  // Preserves any existing ordering on patch faces.
609  // - Two passes, do non-flip faces first and flip faces second. This
610  // guarantees that when e.g. creating a cyclic all faces from one
611  // side come first and faces from the other side next.
612  //
613  polyTopoChange meshMod(mesh);
614  {
615  PackedBoolList modifiedFace(mesh.nFaces());
616 
617  forAll(selectors, selectorI)
618  {
619  const word& groupName = selectors[selectorI].name();
620 
621  const faceZone& fZone = mesh.faceZones()[groupName];
622 
623  const dictionary& dict =
624  selectors[selectorI].dict().optionalSubDict("patches");
625 
626  FixedList<label, 2> newPatchIDs;
627  forAll(patchEntryKeywords, i)
628  {
629  const dictionary& patchDict =
631  (
632  patchEntryKeywords[i],
633  false,
634  false
635  ).dict();
636 
637  const word patchName(patchDict.lookup<word>("name"));
638 
639  newPatchIDs[i] = bMesh.findPatchID(patchName);
640  }
641 
642  const label nModified =
643  createFaces
644  (
645  internalFacesOnly,
646  mesh,
647  fZone,
648  newPatchIDs[0],
649  newPatchIDs[1],
650  meshMod,
651  modifiedFace
652  );
653 
654  Info<< "Converted " << nModified
655  << " faces into boundary faces on ";
656  if (newPatchIDs[0] == newPatchIDs[1])
657  {
658  Info<< "patch '" << bMesh[newPatchIDs[0]].name() << "'";
659  }
660  else
661  {
662  Info<< "patches '" << bMesh[newPatchIDs[0]].name()
663  << "' and '" << bMesh[newPatchIDs[1]].name() << "'";
664  }
665  Info<< endl;
666  }
667  }
668 
669  Info<< endl;
670 
671  if (!overwrite)
672  {
673  runTime++;
674  }
675 
676  // Change the mesh. Change points directly (no inflation).
677  autoPtr<polyTopoChangeMap> map = meshMod.changeMesh(mesh, false);
678 
679  // Update fields
680  mesh.topoChange(map);
681 
682  // Change the patch fields
683  HashSet<word> bafflePatchFields;
684  forAll(selectors, selectorI)
685  {
686  const dictionary& dict =
687  selectors[selectorI].dict().optionalSubDict("patches");
688 
689  forAll(patchEntryKeywords, i)
690  {
691  const dictionary& patchDict =
693  (
694  patchEntryKeywords[i],
695  false,
696  false
697  ).dict();
698 
699  const word patchName(patchDict.lookup<word>("name"));
700 
701  if (!bafflePatchFields.found(patchName))
702  {
703  bafflePatchFields.insert(patchName);
704 
705  dictionary patchFieldsDict =
706  patchDict.subOrEmptyDict("patchFields");
707 
708  fvMeshTools::setPatchFields
709  (
710  mesh,
711  bMesh.findPatchID(patchName),
712  patchFieldsDict
713  );
714  }
715  else
716  {
717  if (patchDict.found("patchFields"))
718  {
720  << "Patch field settings found in " << patchDict.name()
721  << " for patch '" << patchName << "', but fields have "
722  << "already been set for this patch so these settings "
723  << "will not be used" << endl;
724  }
725  }
726  }
727  }
728 
729  // Move mesh (since morphing might not do this)
730  if (map().hasMotionPoints())
731  {
732  mesh.setPoints(map().preMotionPoints());
733  }
734 
735  // Remove any now zero-sized patches
736  filterPatches(mesh, bafflePatches);
737 
738  if (overwrite)
739  {
740  mesh.setInstance(oldInstance);
741  }
742 
743  Info<< "Writing mesh to " << runTime.timeName() << endl;
744 
745  mesh.write();
746 
747  Info<< "End\n" << endl;
748 
749  return 0;
750 }
751 
752 
753 // ************************************************************************* //
Foam::surfaceFields.
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:453
A HashTable with keys but without contents.
Definition: HashSet.H:59
dictionary dict
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:141
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:663
label findZoneID(const word &zoneName) const
Find zone index given a name.
Definition: MeshZones.C:341
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
bool remove(const word &)
Remove an entry specified by keyword.
Definition: dictionary.C:1316
const entry & lookupEntryBackwardsCompatible(const wordList &, bool recursive, bool patternMatch) const
Find and return an entry data stream if present, trying a list.
Definition: dictionary.C:842
autoPtr< polyTopoChangeMap > changeMesh(polyMesh &mesh, const bool inflate, const bool syncParallel=true, const bool orderCells=false, const bool orderPoints=false)
Inplace changes mesh without change of patches.
const word & name() const
Return name.
const word & name() const
Return name.
Definition: IOobject.H:315
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:221
bool set(const label) const
Is element set.
Definition: PtrListI.H:65
List of IOobjects with searching and retrieving facilities.
Definition: IOobjectList.H:50
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Class describing modification of a face.
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
A 1D vector of objects of type <T> with a fixed size <Size>.
Definition: FixedList.H:54
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
label nInternalFaces() const
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1243
label nFaces() const
bool inGroup(const word &) const
Test if in group.
const boolList & flipMap() const
Return face flip map.
Definition: faceZone.H:249
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
bool checkParallelSync(const bool report=false) const
Check whether all procs have all zones and in same order. Return.
Definition: MeshZones.C:428
bool optionFound(const word &opt) const
Return true if the named option is found.
Definition: argListI.H:114
A simple wrapper around bool so that it can be read as a word: true/false, on/off, yes/no, y/n, t/f, or none/any.
Definition: Switch.H:60
virtual const dictionary & dict() const =0
Return dictionary if this entry is a dictionary.
label findPatchID(const word &patchName) const
Find patch index given a name.
label whichFace(const label globalCellID) const
Helper function to re-direct to zone::localID(...)
Definition: faceZone.C:307
Field reading functions for post-processing utilities.
const wordList & inGroups() const
Return the optional groups patch belongs to.
virtual bool write(const bool write=true) const
Write mesh using IO settings from time.
Definition: fvMesh.C:1658
bool insert(const Key &key)
Insert a new entry.
Definition: HashSet.H:111
IOdictionary systemDict(const word &dictName, const argList &args, const objectRegistry &ob, const word &regionName=polyMesh::defaultRegion)
Definition: systemDict.C:92
static bool constraintType(const word &pt)
Return true if the given type is a constraint type.
Definition: fvPatch.C:61
virtual void topoChange(const polyTopoChangeMap &map)
Update mesh corresponding to the given map.
Definition: fvMesh.C:1151
A face addition data class. A face can be inflated either from a point or from another face and can e...
Definition: polyAddFace.H:51
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:1002
label index() const
Return the index of this zone in zone list.
Definition: zone.H:158
const fileName & pointsInstance() const
Return the current instance directory for points.
Definition: polyMesh.C:882
virtual bool coupled() const
Return true if this patch is geometrically coupled (i.e. faces and.
Definition: polyPatch.H:316
const fileName & name() const
Return the dictionary name.
Definition: dictionary.H:109
bool found(const Key &) const
Return true if hashedEntry is found in table.
Definition: HashTable.C:113
A class for handling words, derived from string.
Definition: word.H:59
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:34
void append(const T &)
Append an element at the end of the list.
Definition: ListI.H:178
This boundary condition is not designed to be evaluated; it is assumed that the value is assigned via...
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1237
virtual void setPoints(const pointField &)
Reset the points.
Definition: fvMesh.C:1048
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1224
void setSize(const label)
Reset size of PtrList. If extending the PtrList, new entries are.
Definition: PtrList.C:131
Foam::polyBoundaryMesh.
bool checkParallelSync(const bool report=false) const
Check whether all procs have all patches and in same order. Return.
static const char nl
Definition: Ostream.H:260
objects
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition: Ostream.H:235
void setInstance(const fileName &)
Set the instance for mesh files.
Definition: polyMeshIO.C:77
static autoPtr< faceSelection > New(const word &name, const fvMesh &mesh, const dictionary &dict)
Return a reference to the selected faceSelection.
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
List< word > wordList
A List of words.
Definition: fileName.H:54
A bit-packed bool list.
label patchi
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
#define WarningInFunction
Report a warning using Foam::Warning.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:70
const meshFaceZones & faceZones() const
Return face zones.
Definition: polyMesh.H:495
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:95
label index() const
Return the index of this patch in the boundaryMesh.
Direct mesh changes based on v1.3 polyTopoChange syntax.
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:306
void set(entry *)
Assign a new entry, overwrite any existing entry.
Definition: dictionary.C:1291
messageStream Info
label n
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
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
A subset of mesh faces organised as a primitive patch.
Definition: faceZone.H:65
static void addNote(const string &)
Add extra notes for the usage information.
Definition: argList.C:159
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
Foam::argList args(argc, argv)
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition: Ostream.H:228
void clearAddressing()
Clear addressing.
Definition: MeshZones.C:387
dictionary subOrEmptyDict(const word &, const bool mustRead=false) const
Find and return a sub-dictionary as a copy, or.
Definition: dictionary.C:1033
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:864
label setAction(const topoAction &action)
For compatibility with polyTopoChange: set topological action.