createBaffles.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2016-2022 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 Application
28  createBaffles
29 
30 Group
31  grpMeshManipulationUtilities
32 
33 Description
34  Makes internal faces into boundary faces. Does not duplicate points, unlike
35  mergeOrSplitBaffles.
36 
37  Note: if any coupled patch face is selected for baffling the opposite
38  member has to be selected for baffling as well.
39 
40  - if the patch already exists will not override it nor its fields
41  - if the patch does not exist it will be created together with 'calculated'
42  patchfields unless the field is mentioned in the patchFields section.
43  - any 0-sized patches (since faces have been moved out) will get removed
44 
45 \*---------------------------------------------------------------------------*/
46 
47 #include "argList.H"
48 #include "Time.H"
49 #include "polyTopoChange.H"
50 #include "polyModifyFace.H"
51 #include "polyAddFace.H"
52 #include "ReadFields.H"
53 #include "volFields.H"
54 #include "surfaceFields.H"
55 #include "fvMeshMapper.H"
56 #include "faceSelection.H"
57 
58 #include "fvMeshTools.H"
59 #include "topoSet.H"
60 #include "processorPolyPatch.H"
61 #include "processorMeshes.H"
62 
63 using namespace Foam;
64 
65 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
66 
67 label addPatch
68 (
69  fvMesh& mesh,
70  const word& patchName,
71  const word& groupName,
72  const dictionary& patchDict
73 )
74 {
75  const polyBoundaryMesh& pbm = mesh.boundaryMesh();
76 
77  if (pbm.findPatchID(patchName) == -1)
78  {
79  autoPtr<polyPatch> ppPtr
80  (
82  (
83  patchName,
84  patchDict,
85  0,
86  pbm
87  )
88  );
89  auto& pp = *ppPtr;
90 
91  if (!groupName.empty())
92  {
93  pp.inGroups().appendUniq(groupName);
94  }
95 
96 
97  // Add patch, create calculated everywhere
99  (
100  mesh,
101  pp,
102  dictionary(), // do not set specialised patchFields
104  true // parallel sync'ed addition
105  );
106  }
107  else
108  {
109  Info<< "Patch '" << patchName
110  << "' already exists. Only "
111  << "moving patch faces - type will remain the same"
112  << endl;
113  }
114 
115  return pbm.findPatchID(patchName);
116 }
117 
118 
119 // Filter out the empty patches.
120 void filterPatches(fvMesh& mesh, const wordHashSet& addedPatchNames)
121 {
122  // Remove any zero-sized ones. Assumes
123  // - processor patches are already only there if needed
124  // - all other patches are available on all processors
125  // - but coupled ones might still be needed, even if zero-size
126  // (e.g. processorCyclic)
127  // See also logic in createPatch.
128  const polyBoundaryMesh& pbm = mesh.boundaryMesh();
129 
130  labelList oldToNew(pbm.size(), -1);
131  label newPatchi = 0;
132  forAll(pbm, patchi)
133  {
134  const polyPatch& pp = pbm[patchi];
135 
136  if (!isA<processorPolyPatch>(pp))
137  {
138  if
139  (
140  isA<coupledPolyPatch>(pp)
141  || returnReduceOr(pp.size())
142  || addedPatchNames.found(pp.name())
143  )
144  {
145  // Coupled (and unknown size) or uncoupled and used
146  oldToNew[patchi] = newPatchi++;
147  }
148  }
149  }
150 
151  forAll(pbm, patchi)
152  {
153  const polyPatch& pp = pbm[patchi];
154 
155  if (isA<processorPolyPatch>(pp))
156  {
157  oldToNew[patchi] = newPatchi++;
158  }
159  }
160 
161 
162  const label nKeepPatches = newPatchi;
163 
164  // Shuffle unused ones to end
165  if (nKeepPatches != pbm.size())
166  {
167  Info<< endl
168  << "Removing zero-sized patches:" << endl << incrIndent;
169 
170  forAll(oldToNew, patchi)
171  {
172  if (oldToNew[patchi] == -1)
173  {
174  Info<< indent << pbm[patchi].name()
175  << " type " << pbm[patchi].type()
176  << " at position " << patchi << endl;
177  oldToNew[patchi] = newPatchi++;
178  }
179  }
180  Info<< decrIndent;
181 
182  fvMeshTools::reorderPatches(mesh, oldToNew, nKeepPatches, true);
183  Info<< endl;
184  }
185 }
186 
187 
188 void modifyOrAddFace
189 (
190  polyTopoChange& meshMod,
191  const face& f,
192  const label facei,
193  const label own,
194  const bool flipFaceFlux,
195  const label newPatchi,
196  const label zoneID,
197  const bool zoneFlip,
198 
199  bitSet& modifiedFace
200 )
201 {
202  if (modifiedFace.set(facei))
203  {
204  // First usage of face. Modify.
205  meshMod.setAction
206  (
208  (
209  f, // modified face
210  facei, // label of face
211  own, // owner
212  -1, // neighbour
213  flipFaceFlux, // face flip
214  newPatchi, // patch for face
215  false, // remove from zone
216  zoneID, // zone for face
217  zoneFlip // face flip in zone
218  )
219  );
220  }
221  else
222  {
223  // Second or more usage of face. Add.
224  meshMod.setAction
225  (
227  (
228  f, // modified face
229  own, // owner
230  -1, // neighbour
231  -1, // master point
232  -1, // master edge
233  facei, // master face
234  flipFaceFlux, // face flip
235  newPatchi, // patch for face
236  zoneID, // zone for face
237  zoneFlip // face flip in zone
238  )
239  );
240  }
241 }
242 
243 
244 // Create faces for fZone faces. Usually newMasterPatches, newSlavePatches
245 // only size one but can be more for duplicate baffle sets
246 void createFaces
247 (
248  const bool internalFacesOnly,
249  const fvMesh& mesh,
250  const faceZone& fZone,
251  const labelList& newMasterPatches,
252  const labelList& newSlavePatches,
253  polyTopoChange& meshMod,
254  bitSet& modifiedFace,
255  label& nModified
256 )
257 {
258  const polyBoundaryMesh& pbm = mesh.boundaryMesh();
259 
260  forAll(newMasterPatches, i)
261  {
262  // Pass 1. Do selected side of zone
263  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
264 
265  for (label facei = 0; facei < mesh.nInternalFaces(); facei++)
266  {
267  label zoneFacei = fZone.whichFace(facei);
268 
269  if (zoneFacei != -1)
270  {
271  if (!fZone.flipMap()[zoneFacei])
272  {
273  // Use owner side of face
274  modifyOrAddFace
275  (
276  meshMod,
277  mesh.faces()[facei], // modified face
278  facei, // label of face
279  mesh.faceOwner()[facei],// owner
280  false, // face flip
281  newMasterPatches[i], // patch for face
282  fZone.index(), // zone for face
283  false, // face flip in zone
284  modifiedFace // modify or add status
285  );
286  }
287  else
288  {
289  // Use neighbour side of face.
290  // To keep faceZone pointing out of original neighbour
291  // we don't need to set faceFlip since that cell
292  // now becomes the owner
293  modifyOrAddFace
294  (
295  meshMod,
296  mesh.faces()[facei].reverseFace(), // modified face
297  facei, // label of face
298  mesh.faceNeighbour()[facei],// owner
299  true, // face flip
300  newMasterPatches[i], // patch for face
301  fZone.index(), // zone for face
302  false, // face flip in zone
303  modifiedFace // modify or add status
304  );
305  }
306 
307  nModified++;
308  }
309  }
310 
311 
312  // Pass 2. Do other side of zone
313  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
314 
315  for (label facei = 0; facei < mesh.nInternalFaces(); facei++)
316  {
317  label zoneFacei = fZone.whichFace(facei);
318 
319  if (zoneFacei != -1)
320  {
321  if (!fZone.flipMap()[zoneFacei])
322  {
323  // Use neighbour side of face
324  modifyOrAddFace
325  (
326  meshMod,
327  mesh.faces()[facei].reverseFace(), // modified face
328  facei, // label of face
329  mesh.faceNeighbour()[facei], // owner
330  true, // face flip
331  newSlavePatches[i], // patch for face
332  fZone.index(), // zone for face
333  true, // face flip in zone
334  modifiedFace // modify or add
335  );
336  }
337  else
338  {
339  // Use owner side of face
340  modifyOrAddFace
341  (
342  meshMod,
343  mesh.faces()[facei], // modified face
344  facei, // label of face
345  mesh.faceOwner()[facei],// owner
346  false, // face flip
347  newSlavePatches[i], // patch for face
348  fZone.index(), // zone for face
349  true, // face flip in zone
350  modifiedFace // modify or add status
351  );
352  }
353  }
354  }
355 
356 
357  // Modify any boundary faces
358  // ~~~~~~~~~~~~~~~~~~~~~~~~~
359 
360  // Normal boundary:
361  // - move to new patch. Might already be back-to-back baffle
362  // you want to add cyclic to. Do warn though.
363  //
364  // Processor boundary:
365  // - do not move to cyclic
366  // - add normal patches though.
367 
368  // For warning once per patch.
369  labelHashSet patchWarned;
370 
371  forAll(pbm, patchi)
372  {
373  const polyPatch& pp = pbm[patchi];
374 
375  const label newMasterPatchi = newMasterPatches[i];
376  const label newSlavePatchi = newSlavePatches[i];
377 
378  if
379  (
380  pp.coupled()
381  && (
382  pbm[newMasterPatchi].coupled()
383  || pbm[newSlavePatchi].coupled()
384  )
385  )
386  {
387  // Do not allow coupled faces to be moved to different
388  // coupled patches.
389  }
390  else if (pp.coupled() || !internalFacesOnly)
391  {
392  forAll(pp, i)
393  {
394  label facei = pp.start()+i;
395 
396  label zoneFacei = fZone.whichFace(facei);
397 
398  if (zoneFacei != -1)
399  {
400  if (patchWarned.insert(patchi))
401  {
403  << "Found boundary face (in patch "
404  << pp.name()
405  << ") in faceZone " << fZone.name()
406  << " to convert to baffle patches "
407  << pbm[newMasterPatchi].name() << "/"
408  << pbm[newSlavePatchi].name()
409  << endl
410  << " Set internalFacesOnly to true in the"
411  << " createBaffles control dictionary if you"
412  << " don't wish to convert boundary faces."
413  << endl;
414  }
415 
416  modifyOrAddFace
417  (
418  meshMod,
419  mesh.faces()[facei], // modified face
420  facei, // label of face
421  mesh.faceOwner()[facei], // owner
422  false, // face flip
423  fZone.flipMap()[zoneFacei]
424  ? newSlavePatchi
425  : newMasterPatchi, // patch for face
426  fZone.index(), // zone for face
427  fZone.flipMap()[zoneFacei], // face flip in zone
428  modifiedFace // modify or add
429  );
430 
431  nModified++;
432  }
433  }
434  }
435  }
436  }
437 }
438 
439 
440 int main(int argc, char *argv[])
441 {
443  (
444  "Makes internal faces into boundary faces.\n"
445  "Does not duplicate points."
446  );
447 
448  argList::addOption("dict", "file", "Alternative createBafflesDict");
449  #include "addOverwriteOption.H"
450  #include "addRegionOption.H"
451 
452  argList::noFunctionObjects(); // Never use function objects
453 
454  #include "setRootCase.H"
455  #include "createTime.H"
456  #include "createNamedMesh.H"
457 
458 
459  const bool overwrite = args.found("overwrite");
460 
461  const word oldInstance = mesh.pointsInstance();
462 
463  const word dictName("createBafflesDict");
464  #include "setSystemMeshDictionaryIO.H"
465 
466  bool internalFacesOnly(false);
467 
468  bool noFields(false);
469 
470  PtrList<faceSelection> selectors;
471  {
472  Info<< "Reading baffle criteria from " << dictIO.name() << nl << endl;
474 
475  internalFacesOnly = dict.get<bool>("internalFacesOnly");
476  noFields = dict.getOrDefault("noFields", false);
477 
478  const dictionary& selectionsDict = dict.subDict("baffles");
479 
480  selectors.resize(selectionsDict.size());
481 
482  label nselect = 0;
483  for (const entry& dEntry : selectionsDict)
484  {
485  if (dEntry.isDict())
486  {
487  selectors.set
488  (
489  nselect,
490  faceSelection::New(dEntry.keyword(), mesh, dEntry.dict())
491  );
492 
493  ++nselect;
494  }
495  }
496 
497  selectors.resize(nselect);
498  }
499 
500 
501  if (internalFacesOnly)
502  {
503  Info<< "Not converting faces on non-coupled patches." << nl << endl;
504  }
505 
506 
507  // Read objects in time directory
508  IOobjectList objects(mesh, runTime.timeName());
509 
510  // Read vol fields.
511  Info<< "Reading geometric fields" << nl << endl;
513  if (!noFields) ReadFields(mesh, objects, vsFlds);
514 
516  if (!noFields) ReadFields(mesh, objects, vvFlds);
517 
519  if (!noFields) ReadFields(mesh, objects, vstFlds);
520 
521  PtrList<volSymmTensorField> vsymtFlds;
522  if (!noFields) ReadFields(mesh, objects, vsymtFlds);
523 
525  if (!noFields) ReadFields(mesh, objects, vtFlds);
526 
527  // Read surface fields.
528 
530  if (!noFields) ReadFields(mesh, objects, ssFlds);
531 
533  if (!noFields) ReadFields(mesh, objects, svFlds);
534 
536  if (!noFields) ReadFields(mesh, objects, sstFlds);
537 
539  if (!noFields) ReadFields(mesh, objects, ssymtFlds);
540 
542  if (!noFields) ReadFields(mesh, objects, stFlds);
543 
544 
545 
546 
547  // Creating (if necessary) faceZones
548  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
549 
550  forAll(selectors, selectorI)
551  {
552  const word& name = selectors[selectorI].name();
553 
554  if (mesh.faceZones().findZoneID(name) == -1)
555  {
557  const label zoneID = mesh.faceZones().size();
558 
559  mesh.faceZones().setSize(zoneID+1);
560  mesh.faceZones().set
561  (
562  zoneID,
563  new faceZone(name, labelList(), false, zoneID, mesh.faceZones())
564  );
565  }
566  }
567 
568 
569  // Select faces
570  // ~~~~~~~~~~~~
571 
572  //- Per face zoneID it is in and flip status.
573  labelList faceToZoneID(mesh.nFaces(), -1);
574  boolList faceToFlip(mesh.nFaces(), false);
575  forAll(selectors, selectorI)
576  {
577  const word& name = selectors[selectorI].name();
578  label zoneID = mesh.faceZones().findZoneID(name);
579 
580  selectors[selectorI].select(zoneID, faceToZoneID, faceToFlip);
581  }
582 
583  // Add faces to faceZones
584  labelList nFaces(mesh.faceZones().size(), Zero);
585  forAll(faceToZoneID, facei)
586  {
587  label zoneID = faceToZoneID[facei];
588  if (zoneID != -1)
589  {
590  nFaces[zoneID]++;
591  }
592  }
593 
594  forAll(selectors, selectorI)
595  {
596  const word& name = selectors[selectorI].name();
597  label zoneID = mesh.faceZones().findZoneID(name);
598 
599  label& n = nFaces[zoneID];
600  labelList addr(n);
601  boolList flip(n);
602  n = 0;
603  forAll(faceToZoneID, facei)
604  {
605  label zone = faceToZoneID[facei];
606  if (zone == zoneID)
607  {
608  addr[n] = facei;
609  flip[n] = faceToFlip[facei];
610  ++n;
611  }
612  }
613 
614  Info<< "Created zone " << name
615  << " at index " << zoneID
616  << " with " << returnReduce(n, sumOp<label>()) << " faces" << endl;
617 
618  mesh.faceZones().set
619  (
620  zoneID,
621  new faceZone
622  (
623  name,
624  std::move(addr),
625  std::move(flip),
626  zoneID,
627  mesh.faceZones()
628  )
629  );
630  }
631 
632 
633 
634  // Count patches to add
635  // ~~~~~~~~~~~~~~~~~~~~
636  wordHashSet bafflePatches;
637  {
638  forAll(selectors, selectorI)
639  {
640  const dictionary& dict = selectors[selectorI].dict();
641 
642  if (dict.found("patches"))
643  {
644  for (const entry& dEntry : dict.subDict("patches"))
645  {
646  const word patchName(dEntry.dict().get<word>("name"));
647 
648  bafflePatches.insert(patchName);
649  }
650  }
651  else
652  {
653  const word masterName = selectors[selectorI].name() + "_master";
654  bafflePatches.insert(masterName);
655 
656  const word slaveName = selectors[selectorI].name() + "_slave";
657  bafflePatches.insert(slaveName);
658  }
659  }
660  }
661 
662 
663  // Create baffles
664  // ~~~~~~~~~~~~~~
665  // Is done in multiple steps
666  // - create patches with 'calculated' patchFields
667  // - move faces into these patches
668  // - change the patchFields to the wanted type
669  // This order is done so e.g. fixedJump works:
670  // - you cannot create patchfields at the same time as patches since
671  // they do an evaluate upon construction
672  // - you want to create the patchField only after you have faces
673  // so you don't get the 'create-from-nothing' mapping problem.
674 
675 
676  // Pass 1: add patches
677  // ~~~~~~~~~~~~~~~~~~~
678 
679  // wordHashSet addedPatches;
680  {
681  const polyBoundaryMesh& pbm = mesh.boundaryMesh();
682  forAll(selectors, selectorI)
683  {
684  const dictionary& dict = selectors[selectorI].dict();
685  const word& groupName = selectors[selectorI].name();
686 
687  if (dict.found("patches"))
688  {
689  for (const entry& dEntry : dict.subDict("patches"))
690  {
691  const dictionary& dict = dEntry.dict();
692 
693  const word patchName(dict.get<word>("name"));
694 
695  if (pbm.findPatchID(patchName) == -1)
696  {
697  dictionary patchDict = dict;
698  patchDict.set("nFaces", 0);
699  patchDict.set("startFace", 0);
700 
701  // Note: do not set coupleGroup if constructed from
702  // baffles so you have freedom specifying it
703  // yourself.
704  //patchDict.set("coupleGroup", groupName);
705 
706  addPatch(mesh, patchName, groupName, patchDict);
707  }
708  else
709  {
710  Info<< "Patch '" << patchName
711  << "' already exists. Only "
712  << "moving patch faces - type will remain the same"
713  << endl;
714  }
715  }
716  }
717  else
718  {
719  const dictionary& patchSource = dict.subDict("patchPairs");
720  const word masterName = groupName + "_master";
721  const word slaveName = groupName + "_slave";
722 
723  word groupNameMaster = groupName;
724  word groupNameSlave = groupName;
725 
726 
727  dictionary patchDictMaster(patchSource);
728  patchDictMaster.set("nFaces", 0);
729  patchDictMaster.set("startFace", 0);
730  patchDictMaster.set("coupleGroup", groupName);
731 
732  dictionary patchDictSlave(patchDictMaster);
733 
734  // Note: This is added for the particular case where we want
735  // master and slave in different groupNames
736  // (ie 3D thermal baffles)
737 
738  const bool sameGroup =
739  patchSource.getOrDefault("sameGroup", true);
740 
741  if (!sameGroup)
742  {
743  groupNameMaster = groupName + "Group_master";
744  groupNameSlave = groupName + "Group_slave";
745  patchDictMaster.set("coupleGroup", groupNameMaster);
746  patchDictSlave.set("coupleGroup", groupNameSlave);
747  }
748 
749  addPatch(mesh, masterName, groupNameMaster, patchDictMaster);
750  addPatch(mesh, slaveName, groupNameSlave, patchDictSlave);
751  }
752  }
753  }
754 
755 
756  // Make sure patches and zoneFaces are synchronised across couples
759 
760 
761 
762  // Mesh change container
763  polyTopoChange meshMod(mesh);
764 
765  const polyBoundaryMesh& pbm = mesh.boundaryMesh();
766 
767 
768  // Do the actual changes. Note:
769  // - loop in incrementing face order (not necessary if faceZone ordered).
770  // Preserves any existing ordering on patch faces.
771  // - two passes, do non-flip faces first and flip faces second. This
772  // guarantees that when e.g. creating a cyclic all faces from one
773  // side come first and faces from the other side next.
774 
775  // Whether first use of face (modify) or consecutive (add)
776  bitSet modifiedFace(mesh.nFaces());
777  label nModified = 0;
778 
779  forAll(selectors, selectorI)
780  {
781  const word& name = selectors[selectorI].name();
782  label zoneID = mesh.faceZones().findZoneID(name);
783  const faceZone& fZone = mesh.faceZones()[zoneID];
784 
785  const dictionary& dict = selectors[selectorI].dict();
786 
787  DynamicList<label> newMasterPatches;
788  DynamicList<label> newSlavePatches;
789 
790  if (dict.found("patches"))
791  {
792  bool master = true;
793 
794  for (const entry& dEntry : dict.subDict("patches"))
795  {
796  const word patchName(dEntry.dict().get<word>("name"));
797 
798  const label patchi = pbm.findPatchID(patchName);
799 
800  if (master)
801  {
802  newMasterPatches.append(patchi);
803  }
804  else
805  {
806  newSlavePatches.append(patchi);
807  }
808  master = !master;
809  }
810  }
811  else
812  {
813  const word masterName = selectors[selectorI].name() + "_master";
814  newMasterPatches.append(pbm.findPatchID(masterName));
815 
816  const word slaveName = selectors[selectorI].name() + "_slave";
817  newSlavePatches.append(pbm.findPatchID(slaveName));
818  }
819 
820 
821  createFaces
822  (
823  internalFacesOnly,
824  mesh,
825  fZone,
826  newMasterPatches,
827  newSlavePatches,
828  meshMod,
829  modifiedFace,
830  nModified
831  );
832  }
833 
834 
835  Info<< "Converted " << returnReduce(nModified, sumOp<label>())
836  << " faces into boundary faces in patches "
837  << bafflePatches.sortedToc() << nl << endl;
838 
839  if (!overwrite)
840  {
841  ++runTime;
842  }
843 
844  // Change the mesh. Change points directly (no inflation).
845  autoPtr<mapPolyMesh> map = meshMod.changeMesh(mesh, false);
846 
847  // Update fields
848  mesh.updateMesh(map());
849 
850 
851  // Correct boundary faces mapped-out-of-nothing.
852  // This is just a hack to correct the value field.
853  {
854  fvMeshMapper mapper(mesh, map());
855  bool hasWarned = false;
856 
857  for (const word& patchName : bafflePatches)
858  {
859  label patchi = mesh.boundaryMesh().findPatchID(patchName);
860 
861  const fvPatchMapper& pm = mapper.boundaryMap()[patchi];
862 
863  if (pm.sizeBeforeMapping() == 0)
864  {
865  if (!hasWarned)
866  {
867  hasWarned = true;
869  << "Setting field on boundary faces to zero." << endl
870  << "You might have to edit these fields." << endl;
871  }
872 
874  }
875  }
876  }
877 
878 
879  // Pass 2: change patchFields
880  // ~~~~~~~~~~~~~~~~~~~~~~~~~~
881 
882  {
883  const polyBoundaryMesh& pbm = mesh.boundaryMesh();
884  forAll(selectors, selectorI)
885  {
886  const dictionary& dict = selectors[selectorI].dict();
887  if (dict.found("patches"))
888  {
889  for (const entry& dEntry : dict.subDict("patches"))
890  {
891  const dictionary& dict = dEntry.dict();
892 
893  const word patchName(dict.get<word>("name"));
894 
895  label patchi = pbm.findPatchID(patchName);
896 
897  if (dEntry.dict().found("patchFields"))
898  {
899  const dictionary& patchFieldsDict =
900  dEntry.dict().subDict
901  (
902  "patchFields"
903  );
904 
905  fvMeshTools::setPatchFields
906  (
907  mesh,
908  patchi,
909  patchFieldsDict
910  );
911  }
912  }
913  }
914  else
915  {
916  const dictionary& patchSource = dict.subDict("patchPairs");
917 
918  const bool sameGroup =
919  patchSource.getOrDefault("sameGroup", true);
920 
921  const word& groupName = selectors[selectorI].name();
922 
923  if (patchSource.found("patchFields"))
924  {
925  dictionary patchFieldsDict = patchSource.subDict
926  (
927  "patchFields"
928  );
929 
930  if (sameGroup)
931  {
932  // Add coupleGroup to all entries
933  for (entry& dEntry : patchFieldsDict)
934  {
935  if (dEntry.isDict())
936  {
937  dictionary& dict = dEntry.dict();
938  dict.set("coupleGroup", groupName);
939  }
940  }
941 
942  const labelList& patchIDs =
943  pbm.groupPatchIDs()[groupName];
944 
945  forAll(patchIDs, i)
946  {
947  fvMeshTools::setPatchFields
948  (
949  mesh,
950  patchIDs[i],
951  patchFieldsDict
952  );
953  }
954  }
955  else
956  {
957  const word masterPatchName(groupName + "_master");
958  const word slavePatchName(groupName + "_slave");
959 
960  label patchiMaster = pbm.findPatchID(masterPatchName);
961  label patchiSlave = pbm.findPatchID(slavePatchName);
962 
963  fvMeshTools::setPatchFields
964  (
965  mesh,
966  patchiMaster,
967  patchFieldsDict
968  );
969 
970  fvMeshTools::setPatchFields
971  (
972  mesh,
973  patchiSlave,
974  patchFieldsDict
975  );
976  }
977  }
978  }
979  }
980  }
981 
982 
983  // Move mesh (since morphing might not do this)
984  if (map().hasMotionPoints())
985  {
986  mesh.movePoints(map().preMotionPoints());
987  }
988 
989 
990  // Remove any now zero-sized patches
991  filterPatches(mesh, bafflePatches);
992 
993 
994  if (overwrite)
995  {
996  mesh.setInstance(oldInstance);
997  }
998 
999  Info<< "Writing mesh to " << runTime.timeName() << endl;
1000 
1001  mesh.write();
1004 
1005  Info<< "End\n" << endl;
1006 
1007  return 0;
1008 }
1009 
1010 
1011 // ************************************************************************* //
label findPatchID(const word &patchName, const bool allowNotFound=true) const
Find patch index given a name, return -1 if not found.
Foam::surfaceFields.
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:584
bool checkParallelSync(const bool report=false) const
Check whether all procs have all zones and in same order.
Definition: ZoneMesh.C:746
virtual label sizeBeforeMapping() const
Return size of field before mapping.
dictionary dict
static void noFunctionObjects(bool addWithOption=false)
Remove &#39;-noFunctionObjects&#39; option and ignore any occurrences.
Definition: argList.C:514
static void reorderPatches(fvMesh &, const labelList &oldToNew, const label nPatches, const bool validBoundary)
Reorder and remove trailing patches.
Definition: fvMeshTools.C:315
static void addNote(const string &note)
Add extra notes for the usage information.
Definition: argList.C:453
static void removeFiles(const polyMesh &mesh)
Helper: remove all procAddressing files from mesh instance.
void clearAddressing()
Clear addressing.
Definition: ZoneMesh.C:703
void set(const bitSet &bitset)
Set specified bits from another bitset.
Definition: bitSetI.H:583
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:449
wordList ReadFields(const typename GeoMesh::Mesh &mesh, const IOobjectList &objects, PtrList< GeometricField< Type, PatchField, GeoMesh >> &fields, const bool syncPar=true, const bool readOldTime=false)
Read Geometric fields of templated type.
List of IOobjects with searching and retrieving facilities. Implemented as a HashTable, so the various sorted methods should be used if traversing in parallel.
Definition: IOobjectList.H:55
Class describing modification of a face.
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:68
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:120
bool found(const Key &key) const
True if hashed key is found in table.
Definition: HashTableI.H:93
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1110
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:150
virtual bool write(const bool valid=true) const
Write mesh using IO settings from time.
Definition: fvMesh.C:1072
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:49
const wordList & inGroups() const noexcept
The (optional) groups that the patch belongs to.
const word dictName("faMeshDefinition")
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:487
label appendUniq(const T &val)
Append an element if not already in the list.
Definition: List.H:526
virtual void movePoints(const pointField &)
Move points, returns volumes swept by faces in motion.
Definition: fvMesh.C:888
Required Variables.
label whichFace(const label globalCellID) const
Helper function to re-direct to zone::localID(...)
Definition: faceZone.C:365
Field reading functions for post-processing utilities.
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T. FatalIOError if not found, or if the number of tokens is incorrect.
label nFaces() const noexcept
Number of mesh faces.
autoPtr< mapPolyMesh > changeMesh(polyMesh &mesh, const labelUList &patchMap, const bool inflate, const bool syncParallel=true, const bool orderCells=false, const bool orderPoints=false)
Inplace changes mesh without change of patches.
bool insert(const Key &key)
Insert a new entry, not overwriting existing entries.
Definition: HashSet.H:227
T returnReduce(const T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
Class holds all the necessary information for mapping fields associated with fvMesh.
Definition: fvMeshMapper.H:53
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary.
Definition: dictionary.C:453
A face addition data class. A face can be inflated either from a point or from another face and can e...
Definition: polyAddFace.H:47
bool found(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find an entry (const access) with the given keyword.
Definition: dictionaryI.H:100
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:50
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:413
static void zeroPatchFields(fvMesh &mesh, const label patchi)
Change patchField to zero on registered fields.
Definition: fvMeshTools.C:221
static void removeFiles(const polyMesh &)
Helper: remove all sets files from mesh instance.
Definition: topoSet.C:618
const fileName & pointsInstance() const
Return the current instance directory for points.
Definition: polyMesh.C:847
virtual void updateMesh(const mapPolyMesh &mpm)
Update mesh corresponding to the given map.
Definition: fvMesh.C:964
dynamicFvMesh & mesh
virtual bool coupled() const
Return true if this patch is geometrically coupled (i.e. faces and.
Definition: polyPatch.H:464
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:52
label findZoneID(const word &zoneName) const
Find zone index by name, return -1 if not found.
Definition: ZoneMesh.C:513
Base class for mesh zones.
Definition: zone.H:59
A class for handling words, derived from Foam::string.
Definition: word.H:63
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:36
label size() const noexcept
The number of elements in the list.
Definition: UPtrListI.H:99
This boundary condition is not designed to be evaluated; it is assmued that the value is assigned via...
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1104
static void addOption(const word &optName, const string &param="", const string &usage="", bool advanced=false)
Add an option to validOptions with usage information.
Definition: argList.C:376
label nInternalFaces() const noexcept
Number of internal faces.
void append(const T &val)
Copy append an element to the end of this list.
Definition: DynamicList.H:558
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1091
const HashTable< labelList > & groupPatchIDs() const
The patch indices per patch group.
const T * set(const label i) const
Return const pointer to element (can be nullptr), or nullptr for out-of-range access (ie...
Definition: PtrList.H:163
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO...
const word & name() const noexcept
The patch name.
label index() const noexcept
The index of this zone in the zone list.
bool checkParallelSync(const bool report=false) const
Check whether all procs have all patches and in same order.
static word timeName(const scalar t, const int precision=precision_)
Return time name of given scalar time formatted with the given precision.
Definition: Time.C:760
const faceZoneMesh & faceZones() const noexcept
Return face zone mesh.
Definition: polyMesh.H:646
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition: Ostream.H:467
labelList f(nPoints)
void setSize(const label newLen)
Same as resize()
Definition: PtrList.H:183
void resize(const label newLen)
Adjust size of PtrList.
Definition: PtrList.C:95
static autoPtr< faceSelection > New(const word &name, const fvMesh &mesh, const dictionary &dict)
Return a reference to the selected faceSelection.
void setInstance(const fileName &instance, const IOobjectOption::writeOption wOpt=IOobject::AUTO_WRITE)
Set the instance for mesh files.
Definition: polyMeshIO.C:29
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition: bitSet.H:59
#define WarningInFunction
Report a warning using Foam::Warning.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers...
Definition: List.H:55
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:79
List< Key > sortedToc() const
The table of contents (the keys) in sorted order.
Definition: HashTable.C:130
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:441
const word & name() const noexcept
The zone name.
const boolList & flipMap() const noexcept
Return face flip map.
Definition: faceZone.H:325
messageStream Info
Information stream (stdout output on master, null elsewhere)
label n
IOobject dictIO
static autoPtr< polyPatch > New(const word &patchType, const word &name, const label size, const label start, const label index, const polyBoundaryMesh &bm)
Return pointer to a new patch created on freestore from components.
Definition: polyPatchNew.C:28
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
entry * set(entry *entryPtr)
Assign a new entry, overwriting any existing entry.
Definition: dictionary.C:777
A subset of mesh faces organised as a primitive patch.
Definition: faceZone.H:60
Mapping class for a fvPatchField.
Definition: fvPatchMapper.H:53
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T, or return the given default value. FatalIOError if it is found and the number of...
List< label > labelList
A List of labels.
Definition: List.H:62
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:69
Foam::argList args(argc, argv)
bool returnReduceOr(const bool value, const label comm=UPstream::worldComm)
Perform logical (or) MPI Allreduce on a copy. Uses UPstream::reduceOr.
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition: Ostream.H:458
bool found(const word &optName) const
Return true if the named option is found.
Definition: argListI.H:171
Namespace for OpenFOAM.
A keyword and a list of tokens is an &#39;entry&#39;.
Definition: entry.H:63
label setAction(const topoAction &action)
For compatibility with polyTopoChange: set topological action.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:157