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-2025 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 
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 
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 "addNoOverwriteOption.H"
333  #include "addMeshOption.H"
334  #include "addRegionOption.H"
335 
336  #include "setRootCase.H"
339 
341 
342  #include "setNoOverwrite.H"
343 
344  const word oldInstance = mesh.pointsInstance();
345 
346  // Read the system dictionary and global controls
347  const dictionary dict(systemDict("createBafflesDict", args, mesh));
348  Switch internalFacesOnly = dict.lookup<Switch>("internalFacesOnly");
349  const Switch fields = dict.lookupOrDefault("fields", false);
350 
351  // Create face selections
352  PtrList<faceSelection> selectors;
353  {
354  const dictionary& selectionsDict = dict.subDict("baffles");
355 
356  label n = 0;
357  forAllConstIter(dictionary, selectionsDict, iter)
358  {
359  if (iter().isDict())
360  {
361  n++;
362  }
363  }
364  selectors.setSize(n);
365 
366  n = 0;
367  forAllConstIter(dictionary, selectionsDict, iter)
368  {
369  if (iter().isDict())
370  {
371  selectors.set
372  (
373  n++,
374  faceSelection::New(iter().keyword(), mesh, iter().dict())
375  );
376  }
377  }
378  }
379 
380  if (internalFacesOnly)
381  {
382  Info<< "Not converting faces on non-coupled patches." << nl << endl;
383  }
384 
385  // Read fields
386  IOobjectList objects(mesh, runTime.name());
387  if (fields) Info<< "Reading geometric fields" << nl << endl;
388  #include "readVolFields.H"
389  #include "readSurfaceFields.H"
390  #include "readPointFields.H"
391  if (fields) Info<< endl;
392 
393  // Creating faceZones for selectors that are not already faceZones
394  forAll(selectors, selectorI)
395  {
396  const word& name = selectors[selectorI].name();
397 
398  if (mesh.faceZones().findIndex(name) == -1)
399  {
400  const label sz = mesh.faceZones().size();
401 
402  labelList addr(0);
403  boolList flip(0);
404 
405  mesh.faceZones().setSize(sz+1);
406  mesh.faceZones().set
407  (
408  sz,
409  new faceZone(name, addr, flip, mesh.faceZones())
410  );
411  }
412  }
413 
414  // Per face, its zone ID and flip status
415  labelList faceToZoneID(mesh.nFaces(), -1);
416  boolList faceToFlip(mesh.nFaces(), false);
417  forAll(selectors, selectorI)
418  {
419  const word& name = selectors[selectorI].name();
420  label zoneID = mesh.faceZones().findIndex(name);
421 
422  selectors[selectorI].select(zoneID, faceToZoneID, faceToFlip);
423  }
424 
425  // Add faces to faceZones
426  labelList nFaces(mesh.faceZones().size(), 0);
427  forAll(faceToZoneID, facei)
428  {
429  label zoneID = faceToZoneID[facei];
430  if (zoneID != -1)
431  {
432  nFaces[zoneID]++;
433  }
434  }
435  forAll(selectors, selectorI)
436  {
437  const word& name = selectors[selectorI].name();
438  const label zoneID = mesh.faceZones().findIndex(name);
439 
440  label& n = nFaces[zoneID];
441  labelList addr(n);
442  boolList flip(n);
443 
444  n = 0;
445  forAll(faceToZoneID, facei)
446  {
447  label zone = faceToZoneID[facei];
448  if (zone == zoneID)
449  {
450  addr[n] = facei;
451  flip[n] = faceToFlip[facei];
452  n++;
453  }
454  }
455 
456  Info<< "Created zone '" << name << "' with "
457  << returnReduce(n, sumOp<label>()) << " faces" << endl;
458 
459  mesh.faceZones().set
460  (
461  zoneID,
462  new faceZone(name, addr, flip, mesh.faceZones())
463  );
464  }
465 
466  Info<< endl;
467 
468  // Keywords for the two patch entries in each selector dictionary. Note
469  // "owner" and "neighbour" are preferred; "master" and "slave" are still
470  // permitted for backwards compatibility.
471  const Pair<wordList> patchEntryKeywords
472  ({
473  wordList({"owner", "master"}),
474  wordList({"neighbour", "slave"})
475  });
476 
477  // Create the baffles. Notes:
478  //
479  // - This is done in multiple steps
480  // - Patches are created with 'calculated' patchFields
481  // - Faces are moved into these patches
482  // - The patchFields are changed to the desired type
483  // - This helps resolve ordering issues
484  // - patchFields cannot be created at the same time as patches since a
485  // coupled patchField constructor frequently needs access to the
486  // neighbouring patch
487  // - patchFields need to be created after the faces have been added to
488  // the patch so that they don't later have to be mapped from nothing
489  //
490 
491  // Add patches
492  HashSet<word> bafflePatches;
493  forAll(selectors, selectorI)
494  {
495  const word& groupName = selectors[selectorI].name();
496  const dictionary& dict =
497  selectors[selectorI].dict().optionalSubDict("patches");
498 
499  forAll(patchEntryKeywords, i)
500  {
501  dictionary patchDict =
503  (
504  patchEntryKeywords[i],
505  false,
506  false
507  ).dict();
508 
509  const word patchName(patchDict.lookup<word>("name"));
510 
511  if (bMesh.findIndex(patchName) == -1)
512  {
513  Info<< "Adding patch '" << patchName << "' to the mesh" << endl;
514 
515  // Create the patch
516  patchDict.set("nFaces", 0);
517  patchDict.set("startFace", 0);
518  autoPtr<polyPatch> ppPtr
519  (
521  (
522  patchName,
523  patchDict,
524  0,
525  bMesh
526  )
527  );
528  polyPatch& pp = ppPtr();
529 
530  // Check the validity of the patch type
531  if (isA<cyclicPolyPatch>(pp) && Pstream::parRun())
532  {
534  << "Patches of type '" << pp.type() << "' cannot be "
535  << "created in parallel. They must be created in "
536  << "serial and then decomposed." << exit(FatalError);
537  }
538  if
539  (
540  isA<processorPolyPatch>(pp)
541  || isA<nonConformalPolyPatch>(pp)
542  )
543  {
545  << "Patches of type '" << pp.type() << "' cannot be "
546  << "created by createBaffles" << exit(FatalError);
547  }
548 
549  // Add it to the group
550  if (!groupName.empty() && !pp.inGroup(groupName))
551  {
552  pp.inGroups().append(groupName);
553  }
554 
555  // Add the patch to the mesh, creating calculated boundary
556  // conditions everywhere. These will be changed later.
558  }
559  else
560  {
561  Info<< "Patch '" << patchName
562  << "' already exists in the mesh" << endl;
563 
564  patchDict.remove("name");
565  patchDict.remove("patchFields");
566 
567  if (!patchDict.empty())
568  {
570  << "Patch settings found in " << patchDict.name()
571  << " for patch '" << patchName << "', but this patch "
572  << "already exists so these settings will not be used"
573  << endl;
574  }
575  }
576 
577  bafflePatches.insert(patchName);
578  }
579  }
580 
581  // Finish adding the patches to the mesh
583 
584  Info<< endl;
585 
586  // Make sure patches and zoneFaces are synchronised across couples
589 
590  // Do the topology changes. Notes:
591  //
592  // - Loop in incrementing face order (not necessary if faceZone ordered).
593  // Preserves any existing ordering on patch faces.
594  // - Two passes, do non-flip faces first and flip faces second. This
595  // guarantees that when e.g. creating a cyclic all faces from one
596  // side come first and faces from the other side next.
597  //
598  polyTopoChange meshMod(mesh);
599  {
600  PackedBoolList modifiedFace(mesh.nFaces());
601 
602  forAll(selectors, selectorI)
603  {
604  const word& groupName = selectors[selectorI].name();
605 
606  const label fZonei = mesh.faceZones().findIndex(groupName);
607  const faceZone& fZone = mesh.faceZones()[fZonei];
608 
609  const dictionary& dict =
610  selectors[selectorI].dict().optionalSubDict("patches");
611 
612  labelPair newPatchIDs;
613  forAll(patchEntryKeywords, i)
614  {
615  const dictionary& patchDict =
617  (
618  patchEntryKeywords[i],
619  false,
620  false
621  ).dict();
622 
623  const word patchName(patchDict.lookup<word>("name"));
624 
625  newPatchIDs[i] = bMesh.findIndex(patchName);
626  }
627 
628  const label nModified =
629  createFaces
630  (
631  internalFacesOnly,
632  mesh,
633  fZone,
634  newPatchIDs[0],
635  newPatchIDs[1],
636  meshMod,
637  modifiedFace
638  );
639 
640  Info<< "Converted " << nModified
641  << " faces into boundary faces on ";
642  if (newPatchIDs[0] == newPatchIDs[1])
643  {
644  Info<< "patch '" << bMesh[newPatchIDs[0]].name() << "'";
645  }
646  else
647  {
648  Info<< "patches '" << bMesh[newPatchIDs[0]].name()
649  << "' and '" << bMesh[newPatchIDs[1]].name() << "'";
650  }
651  Info<< endl;
652  }
653  }
654 
655  Info<< endl;
656 
657  if (!overwrite)
658  {
659  runTime++;
660  }
661 
662  // Change the mesh. Change points directly
664 
665  // Update fields
666  mesh.topoChange(map);
667 
668  // Change the patch fields
669  HashSet<word> bafflePatchFields;
670  forAll(selectors, selectorI)
671  {
672  const dictionary& dict =
673  selectors[selectorI].dict().optionalSubDict("patches");
674 
675  forAll(patchEntryKeywords, i)
676  {
677  const dictionary& patchDict =
679  (
680  patchEntryKeywords[i],
681  false,
682  false
683  ).dict();
684 
685  const word patchName(patchDict.lookup<word>("name"));
686 
687  if (!bafflePatchFields.found(patchName))
688  {
689  bafflePatchFields.insert(patchName);
690 
691  dictionary patchFieldsDict =
692  patchDict.subOrEmptyDict("patchFields");
693 
694  fvMeshTools::setPatchFields
695  (
696  mesh,
697  bMesh.findIndex(patchName),
698  patchFieldsDict
699  );
700  }
701  else
702  {
703  if (patchDict.found("patchFields"))
704  {
706  << "Patch field settings found in " << patchDict.name()
707  << " for patch '" << patchName << "', but fields have "
708  << "already been set for this patch so these settings "
709  << "will not be used" << endl;
710  }
711  }
712  }
713  }
714 
715  // Remove any now zero-sized patches
716  filterPatches(mesh, bafflePatches);
717 
718  if (overwrite)
719  {
720  mesh.setInstance(oldInstance);
721  }
722 
723  Info<< "Writing mesh to " << runTime.name() << endl;
724 
725  mesh.write();
726 
727  Info<< "End\n" << endl;
728 
729  return 0;
730 }
731 
732 
733 // ************************************************************************* //
Field reading functions for post-processing utilities.
label n
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:476
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:307
void append(const T &)
Append an element at the end of the list.
Definition: ListI.H:178
A bit-packed bool list.
An ordered pair of two objects of type <Type> with first() and second() elements.
Definition: Pair.H:66
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
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:399
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:309
label localIndex(const label globalIndex) const
Map storing the local index for every global index. Used to find.
Definition: Zone.C:221
static void addNote(const string &)
Add extra notes for the usage information.
Definition: argList.C:159
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
const fileName & name() const
Return the dictionary name.
Definition: dictionary.H:111
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
const dictionary & subOrEmptyDict(const word &, const bool mustRead=false) const
Find and return a sub-dictionary.
Definition: dictionary.C:900
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:740
bool remove(const word &)
Remove an entry specified by keyword.
Definition: dictionary.C:1183
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:849
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:718
T lookupOrDefault(const word &, const T &, const bool writeDefault=writeOptionalEntries > 0) const
Find and return a T, if not found return the given default.
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:539
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.
Named list of face indices representing a sub-set of the mesh faces.
Definition: faceZone.H:66
const boolList & flipMap() const
Return face flip map.
Definition: faceZone.H:262
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.
Definition: fvMeshTools.C:104
static label addPatch(fvMesh &mesh, const polyPatch &patch)
Add patch. Inserts patch before all processor patches. Returns the.
Definition: fvMeshTools.C:33
static void addedPatches(fvMesh &mesh)
Complete adding patches.
Definition: fvMeshTools.C:73
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:96
virtual bool write(const bool write=true) const
Write mesh using IO settings from time.
Definition: fvMesh.C:1785
const word & name() const
Return reference to name.
Definition: fvMesh.H:434
virtual void topoChange(const polyTopoChangeMap &map)
Update mesh corresponding to the given map.
Definition: fvMesh.C:1369
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:1344
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1357
const fileName & pointsInstance() const
Return the current instance directory for points.
Definition: polyMesh.C:988
const faceZoneList & faceZones() const
Return face zones.
Definition: polyMesh.H:443
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:401
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1363
void setInstance(const fileName &)
Set the instance for mesh files.
Definition: polyMeshIO.C:91
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
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
#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:234
#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:242
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
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:258
messageStream Info
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition: Ostream.H:235
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
error FatalError
IOdictionary systemDict(const word &dictName, const argList &args, const objectRegistry &ob, const word &regionName=polyMesh::defaultRegion, const fileName &path=fileName::null)
Definition: systemDict.C:93
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:228
static const char nl
Definition: Ostream.H:267
labelList f(nPoints)
objects
dictionary dict
const bool overwrite
Definition: setNoOverwrite.H:1
Foam::argList args(argc, argv)
Foam::surfaceFields.