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-2023 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
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
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
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[newNeighbourPatchi]
313  : bMesh[newOwnerPatchi];
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"
365  #include "createNamedMesh.H"
366 
367  const polyBoundaryMesh& bMesh = mesh.boundaryMesh();
368 
369  const bool overwrite = args.optionFound("overwrite");
370 
371  const word oldInstance = mesh.pointsInstance();
372 
373  // Read the system dictionary and global controls
374  const dictionary dict(systemDict("createBafflesDict", args, mesh));
375  Switch internalFacesOnly = dict.lookup<Switch>("internalFacesOnly");
376  const Switch fields = dict.lookupOrDefault("fields", false);
377 
378  // Create face selections
379  PtrList<faceSelection> selectors;
380  {
381  const dictionary& selectionsDict = dict.subDict("baffles");
382 
383  label n = 0;
384  forAllConstIter(dictionary, selectionsDict, iter)
385  {
386  if (iter().isDict())
387  {
388  n++;
389  }
390  }
391  selectors.setSize(n);
392 
393  n = 0;
394  forAllConstIter(dictionary, selectionsDict, iter)
395  {
396  if (iter().isDict())
397  {
398  selectors.set
399  (
400  n++,
401  faceSelection::New(iter().keyword(), mesh, iter().dict())
402  );
403  }
404  }
405  }
406 
407  if (internalFacesOnly)
408  {
409  Info<< "Not converting faces on non-coupled patches." << nl << endl;
410  }
411 
412  // Read fields
413  IOobjectList objects(mesh, runTime.name());
414  if (fields) Info<< "Reading geometric fields" << nl << endl;
415  #include "readVolFields.H"
416  #include "readSurfaceFields.H"
417  #include "readPointFields.H"
418  if (fields) Info<< endl;
419 
420  // Creating faceZones for selectors that are not already faceZones
421  forAll(selectors, selectorI)
422  {
423  const word& name = selectors[selectorI].name();
424 
425  if (mesh.faceZones().findZoneID(name) == -1)
426  {
427  mesh.faceZones().clearAddressing();
428 
429  const label sz = mesh.faceZones().size();
430 
431  labelList addr(0);
432  boolList flip(0);
433 
434  mesh.faceZones().setSize(sz+1);
435  mesh.faceZones().set
436  (
437  sz,
438  new faceZone(name, addr, flip, sz, mesh.faceZones())
439  );
440  }
441  }
442 
443  // Per face, its zone ID and flip status
444  labelList faceToZoneID(mesh.nFaces(), -1);
445  boolList faceToFlip(mesh.nFaces(), false);
446  forAll(selectors, selectorI)
447  {
448  const word& name = selectors[selectorI].name();
449  label zoneID = mesh.faceZones().findZoneID(name);
450 
451  selectors[selectorI].select(zoneID, faceToZoneID, faceToFlip);
452  }
453 
454  // Add faces to faceZones
455  labelList nFaces(mesh.faceZones().size(), 0);
456  forAll(faceToZoneID, facei)
457  {
458  label zoneID = faceToZoneID[facei];
459  if (zoneID != -1)
460  {
461  nFaces[zoneID]++;
462  }
463  }
464  forAll(selectors, selectorI)
465  {
466  const word& name = selectors[selectorI].name();
467  const label zoneID = mesh.faceZones().findZoneID(name);
468 
469  label& n = nFaces[zoneID];
470  labelList addr(n);
471  boolList flip(n);
472 
473  n = 0;
474  forAll(faceToZoneID, facei)
475  {
476  label zone = faceToZoneID[facei];
477  if (zone == zoneID)
478  {
479  addr[n] = facei;
480  flip[n] = faceToFlip[facei];
481  n++;
482  }
483  }
484 
485  Info<< "Created zone '" << name << "' with "
486  << returnReduce(n, sumOp<label>()) << " faces" << endl;
487 
488  mesh.faceZones().set
489  (
490  zoneID,
491  new faceZone(name, addr, flip, zoneID, mesh.faceZones())
492  );
493  }
494 
495  Info<< endl;
496 
497  // Keywords for the two patch entries in each selector dictionary. Note
498  // "owner" and "neighbour" are preferred; "master" and "slave" are still
499  // permitted for backwards compatibility.
500  const FixedList<wordList, 2> patchEntryKeywords
501  ({
502  wordList({"owner", "master"}),
503  wordList({"neighbour", "slave"})
504  });
505 
506  // Create the baffles. Notes:
507  //
508  // - This is done in multiple steps
509  // - Patches are created with 'calculated' patchFields
510  // - Faces are moved into these patches
511  // - The patchFields are changed to the desired type
512  // - This helps resolve ordering issues
513  // - patchFields cannot be created at the same time as patches since a
514  // coupled patchField constructor frequently needs access to the
515  // neighbouring patch
516  // - patchFields need to be created after the faces have been added to
517  // the patch so that they don't later have to be mapped from nothing
518  //
519 
520  // Add patches
521  HashSet<word> bafflePatches;
522  forAll(selectors, selectorI)
523  {
524  const word& groupName = selectors[selectorI].name();
525  const dictionary& dict =
526  selectors[selectorI].dict().optionalSubDict("patches");
527 
528  forAll(patchEntryKeywords, i)
529  {
530  dictionary patchDict =
532  (
533  patchEntryKeywords[i],
534  false,
535  false
536  ).dict();
537 
538  const word patchName(patchDict.lookup<word>("name"));
539 
540  if (bMesh.findPatchID(patchName) == -1)
541  {
542  Info<< "Adding patch '" << patchName << "' to the mesh" << endl;
543 
544  // Create the patch
545  patchDict.set("nFaces", 0);
546  patchDict.set("startFace", 0);
547  autoPtr<polyPatch> ppPtr
548  (
550  (
551  patchName,
552  patchDict,
553  0,
554  bMesh
555  )
556  );
557  polyPatch& pp = ppPtr();
558 
559  // Add it to the group
560  if (!groupName.empty() && !pp.inGroup(groupName))
561  {
562  pp.inGroups().append(groupName);
563  }
564 
565  // Add the patch to the mesh, creating calculated boundary
566  // conditions everywhere. These will be changed later.
568  (
569  mesh,
570  pp,
571  dictionary(), // do not set specialised patchFields
573  true // parallel sync'ed addition
574  );
575  }
576  else
577  {
578  Info<< "Patch '" << patchName
579  << "' already exists in the mesh" << endl;
580 
581  patchDict.remove("name");
582  patchDict.remove("patchFields");
583 
584  if (!patchDict.empty())
585  {
587  << "Patch settings found in " << patchDict.name()
588  << " for patch '" << patchName << "', but this patch "
589  << "already exists so these settings will not be used"
590  << endl;
591  }
592  }
593 
594  bafflePatches.insert(patchName);
595  }
596  }
597 
598  Info<< endl;
599 
600  // Make sure patches and zoneFaces are synchronised across couples
601  mesh.boundaryMesh().checkParallelSync(true);
602  mesh.faceZones().checkParallelSync(true);
603 
604  // Do the topology changes. Notes:
605  //
606  // - Loop in incrementing face order (not necessary if faceZone ordered).
607  // Preserves any existing ordering on patch faces.
608  // - Two passes, do non-flip faces first and flip faces second. This
609  // guarantees that when e.g. creating a cyclic all faces from one
610  // side come first and faces from the other side next.
611  //
612  polyTopoChange meshMod(mesh);
613  {
614  PackedBoolList modifiedFace(mesh.nFaces());
615 
616  forAll(selectors, selectorI)
617  {
618  const word& groupName = selectors[selectorI].name();
619 
620  const faceZone& fZone = mesh.faceZones()[groupName];
621 
622  const dictionary& dict =
623  selectors[selectorI].dict().optionalSubDict("patches");
624 
625  FixedList<label, 2> newPatchIDs;
626  forAll(patchEntryKeywords, i)
627  {
628  const dictionary& patchDict =
630  (
631  patchEntryKeywords[i],
632  false,
633  false
634  ).dict();
635 
636  const word patchName(patchDict.lookup<word>("name"));
637 
638  newPatchIDs[i] = bMesh.findPatchID(patchName);
639  }
640 
641  const label nModified =
642  createFaces
643  (
644  internalFacesOnly,
645  mesh,
646  fZone,
647  newPatchIDs[0],
648  newPatchIDs[1],
649  meshMod,
650  modifiedFace
651  );
652 
653  Info<< "Converted " << nModified
654  << " faces into boundary faces on ";
655  if (newPatchIDs[0] == newPatchIDs[1])
656  {
657  Info<< "patch '" << bMesh[newPatchIDs[0]].name() << "'";
658  }
659  else
660  {
661  Info<< "patches '" << bMesh[newPatchIDs[0]].name()
662  << "' and '" << bMesh[newPatchIDs[1]].name() << "'";
663  }
664  Info<< endl;
665  }
666  }
667 
668  Info<< endl;
669 
670  if (!overwrite)
671  {
672  runTime++;
673  }
674 
675  // Change the mesh. Change points directly (no inflation).
676  autoPtr<polyTopoChangeMap> map = meshMod.changeMesh(mesh, false);
677 
678  // Update fields
679  mesh.topoChange(map);
680 
681  // Change the patch fields
682  HashSet<word> bafflePatchFields;
683  forAll(selectors, selectorI)
684  {
685  const dictionary& dict =
686  selectors[selectorI].dict().optionalSubDict("patches");
687 
688  forAll(patchEntryKeywords, i)
689  {
690  const dictionary& patchDict =
692  (
693  patchEntryKeywords[i],
694  false,
695  false
696  ).dict();
697 
698  const word patchName(patchDict.lookup<word>("name"));
699 
700  if (!bafflePatchFields.found(patchName))
701  {
702  bafflePatchFields.insert(patchName);
703 
704  dictionary patchFieldsDict =
705  patchDict.subOrEmptyDict("patchFields");
706 
707  fvMeshTools::setPatchFields
708  (
709  mesh,
710  bMesh.findPatchID(patchName),
711  patchFieldsDict
712  );
713  }
714  else
715  {
716  if (patchDict.found("patchFields"))
717  {
719  << "Patch field settings found in " << patchDict.name()
720  << " for patch '" << patchName << "', but fields have "
721  << "already been set for this patch so these settings "
722  << "will not be used" << endl;
723  }
724  }
725  }
726  }
727 
728  // Move mesh (since morphing might not do this)
729  if (map().hasMotionPoints())
730  {
731  mesh.setPoints(map().preMotionPoints());
732  }
733 
734  // Remove any now zero-sized patches
735  filterPatches(mesh, bafflePatches);
736 
737  if (overwrite)
738  {
739  mesh.setInstance(oldInstance);
740  }
741 
742  Info<< "Writing mesh to " << runTime.name() << endl;
743 
744  mesh.write();
745 
746  Info<< "End\n" << endl;
747 
748  return 0;
749 }
750 
751 
752 // ************************************************************************* //
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:111
bool found(const Key &) const
Return true if hashedEntry is found in table.
Definition: HashTable.C:113
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
label findZoneID(const word &zoneName) const
Find zone index given a name.
Definition: MeshZones.C:341
void clearAddressing()
Clear addressing.
Definition: MeshZones.C:387
bool checkParallelSync(const bool report=false) const
Check whether all procs have all zones and in same order. Return.
Definition: MeshZones.C:428
A bit-packed bool list.
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:65
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 UPtrList.
Definition: UPtrListI.H:29
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:109
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
dictionary subOrEmptyDict(const word &, const bool mustRead=false) const
Find and return a sub-dictionary as a copy, or.
Definition: dictionary.C:1049
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:860
void set(entry *)
Assign a new entry, overwrite any existing entry.
Definition: dictionary.C:1307
bool remove(const word &)
Remove an entry specified by keyword.
Definition: dictionary.C:1332
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:998
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:838
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:659
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:68
const boolList & flipMap() const
Return face flip map.
Definition: faceZone.H:249
label whichFace(const label globalCellID) const
Helper function to re-direct to zone::localID(...)
Definition: faceZone.C:307
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:101
virtual bool write(const bool write=true) const
Write mesh using IO settings from time.
Definition: fvMesh.C:1736
virtual void topoChange(const polyTopoChangeMap &map)
Update mesh corresponding to the given map.
Definition: fvMesh.C:1236
virtual void setPoints(const pointField &)
Reset the points.
Definition: fvMesh.C:1133
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.
A face addition data class. A face can be inflated either from a point or from another face and can e...
Definition: polyAddFace.H:54
Foam::polyBoundaryMesh.
label findPatchID(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.
const meshFaceZones & faceZones() const
Return face zones.
Definition: polyMesh.H:447
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1373
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1386
const fileName & pointsInstance() const
Return the current instance directory for points.
Definition: polyMesh.C:1019
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:405
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1392
void setInstance(const fileName &)
Set the instance for mesh files.
Definition: polyMeshIO.C:78
Class describing modification of a face.
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 inflate, const bool syncParallel=true, const bool orderCells=false, const bool orderPoints=false)
Inplace changes mesh without change of patches.
label setAction(const topoAction &action)
For compatibility with polyTopoChange: set topological action.
label nInternalFaces() const
label nFaces() const
A class for handling words, derived from string.
Definition: word.H:62
Base class for zones.
Definition: zone.H:60
label index() const
Return the index of this zone in zone list.
Definition: zone.H:158
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
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:230
#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:235
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:251
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:228
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:221
static const char nl
Definition: Ostream.H:260
labelList f(nPoints)
objects
dictionary dict
Foam::argList args(argc, argv)
Foam::surfaceFields.