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