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-2023 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 {
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.
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.contains(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 labelUList& newMasterPatches,
252  const labelUList& newSlavePatches,
253  polyTopoChange& meshMod,
254  bitSet& modifiedFace,
255  label& nModified
256 )
257 {
259 
260  forAll(newMasterPatches, i)
261  {
262  const label newMasterPatchi = newMasterPatches[i];
263  const label newSlavePatchi = newSlavePatches[i];
264 
265  // Pass 1. Do selected side of zone
266  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
267 
268  for (label facei = 0; facei < mesh.nInternalFaces(); ++facei)
269  {
270  const label zoneFacei = fZone.whichFace(facei);
271 
272  if (zoneFacei != -1)
273  {
274  if (!fZone.flipMap()[zoneFacei])
275  {
276  // Use owner side of face
277  modifyOrAddFace
278  (
279  meshMod,
280  mesh.faces()[facei], // modified face
281  facei, // label of face
282  mesh.faceOwner()[facei],// owner
283  false, // face flip
284  newMasterPatchi, // patch for face
285  fZone.index(), // zone for face
286  false, // face flip in zone
287  modifiedFace // modify or add status
288  );
289  }
290  else
291  {
292  // Use neighbour side of face.
293  // To keep faceZone pointing out of original neighbour
294  // we don't need to set faceFlip since that cell
295  // now becomes the owner
296  modifyOrAddFace
297  (
298  meshMod,
299  mesh.faces()[facei].reverseFace(), // modified face
300  facei, // label of face
301  mesh.faceNeighbour()[facei],// owner
302  true, // face flip
303  newMasterPatchi, // patch for face
304  fZone.index(), // zone for face
305  false, // face flip in zone
306  modifiedFace // modify or add status
307  );
308  }
309 
310  ++nModified;
311  }
312  }
313 
314 
315  // Pass 2. Do other side of zone
316  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
317 
318  for (label facei = 0; facei < mesh.nInternalFaces(); ++facei)
319  {
320  const label zoneFacei = fZone.whichFace(facei);
321 
322  if (zoneFacei != -1)
323  {
324  if (!fZone.flipMap()[zoneFacei])
325  {
326  // Use neighbour side of face
327  modifyOrAddFace
328  (
329  meshMod,
330  mesh.faces()[facei].reverseFace(), // modified face
331  facei, // label of face
332  mesh.faceNeighbour()[facei], // owner
333  true, // face flip
334  newSlavePatchi, // patch for face
335  fZone.index(), // zone for face
336  true, // face flip in zone
337  modifiedFace // modify or add
338  );
339  }
340  else
341  {
342  // Use owner side of face
343  modifyOrAddFace
344  (
345  meshMod,
346  mesh.faces()[facei], // modified face
347  facei, // label of face
348  mesh.faceOwner()[facei],// owner
349  false, // face flip
350  newSlavePatchi, // patch for face
351  fZone.index(), // zone for face
352  true, // face flip in zone
353  modifiedFace // modify or add status
354  );
355  }
356  }
357  }
358 
359 
360  // Modify any boundary faces
361  // ~~~~~~~~~~~~~~~~~~~~~~~~~
362 
363  // Normal boundary:
364  // - move to new patch. Might already be back-to-back baffle
365  // you want to add cyclic to. Do warn though.
366  //
367  // Processor boundary:
368  // - do not move to cyclic
369  // - add normal patches though.
370 
371  // For warning once per patch.
372  labelHashSet patchWarned;
373 
374  forAll(pbm, patchi)
375  {
376  const polyPatch& pp = pbm[patchi];
377  const bool isCoupled = pp.coupled();
378 
379  if
380  (
381  isCoupled
382  && (
383  pbm[newMasterPatchi].coupled()
384  || pbm[newSlavePatchi].coupled()
385  )
386  )
387  {
388  // Do not allow coupled faces to be moved to different
389  // coupled patches.
390  }
391  else if (isCoupled || !internalFacesOnly)
392  {
393  forAll(pp, i)
394  {
395  label facei = pp.start()+i;
396 
397  label zoneFacei = fZone.whichFace(facei);
398 
399  if (zoneFacei != -1)
400  {
401  if (!isCoupled && patchWarned.insert(patchi))
402  {
404  << "Found boundary face (patch " << pp.name()
405  << ") in faceZone " << fZone.name()
406  << " to convert to baffle patches "
407  << pbm[newMasterPatchi].name() << "/"
408  << pbm[newSlavePatchi].name() << nl
409  << " Set internalFacesOnly to true in the"
410  << " createBaffles control dictionary if you"
411  << " don't wish to convert boundary faces."
412  << endl;
413  }
414 
415  modifyOrAddFace
416  (
417  meshMod,
418  mesh.faces()[facei], // modified face
419  facei, // label of face
420  mesh.faceOwner()[facei], // owner
421  false, // face flip
422  fZone.flipMap()[zoneFacei]
423  ? newSlavePatchi
424  : newMasterPatchi, // patch for face
425  fZone.index(), // zone for face
426  fZone.flipMap()[zoneFacei], // face flip in zone
427  modifiedFace // modify or add
428  );
429 
430  ++nModified;
431  }
432  }
433  }
434  }
435  }
436 }
437 
438 
439 int main(int argc, char *argv[])
440 {
442  (
443  "Makes internal faces into boundary faces.\n"
444  "Does not duplicate points."
445  );
446 
447  argList::addOption("dict", "file", "Alternative createBafflesDict");
448  #include "addOverwriteOption.H"
449  #include "addRegionOption.H"
450 
451  argList::noFunctionObjects(); // Never use function objects
452 
453  #include "setRootCase.H"
454  #include "createTime.H"
455  #include "createNamedMesh.H"
456 
457 
458  const bool overwrite = args.found("overwrite");
459 
460  const word oldInstance = mesh.pointsInstance();
461 
462  const word dictName("createBafflesDict");
463  #include "setSystemMeshDictionaryIO.H"
464 
465  bool internalFacesOnly(false);
466 
467  bool noFields(false);
468 
469  PtrList<faceSelection> selectors;
470  {
471  Info<< "Reading baffle criteria from " << dictIO.name() << nl << endl;
473 
474  internalFacesOnly = dict.get<bool>("internalFacesOnly");
475  noFields = dict.getOrDefault("noFields", false);
476 
477  const dictionary& selectionsDict = dict.subDict("baffles");
478 
479  selectors.resize(selectionsDict.size());
480 
481  label nselect = 0;
482  for (const entry& dEntry : selectionsDict)
483  {
484  if (dEntry.isDict())
485  {
486  selectors.set
487  (
488  nselect,
489  faceSelection::New(dEntry.keyword(), mesh, dEntry.dict())
490  );
491 
492  ++nselect;
493  }
494  }
495 
496  selectors.resize(nselect);
497  }
498 
499 
500  if (internalFacesOnly)
501  {
502  Info<< "Not converting faces on non-coupled patches." << nl << endl;
503  }
504 
505 
506  // Read objects in time directory
507  IOobjectList objects(mesh, runTime.timeName());
508 
509  // Read vol fields.
510  Info<< "Reading geometric fields" << nl << endl;
512  if (!noFields) ReadFields(mesh, objects, vsFlds);
513 
515  if (!noFields) ReadFields(mesh, objects, vvFlds);
516 
518  if (!noFields) ReadFields(mesh, objects, vstFlds);
519 
520  PtrList<volSymmTensorField> vsymtFlds;
521  if (!noFields) ReadFields(mesh, objects, vsymtFlds);
522 
524  if (!noFields) ReadFields(mesh, objects, vtFlds);
525 
526  // Read surface fields.
527 
529  if (!noFields) ReadFields(mesh, objects, ssFlds);
530 
532  if (!noFields) ReadFields(mesh, objects, svFlds);
533 
535  if (!noFields) ReadFields(mesh, objects, sstFlds);
536 
538  if (!noFields) ReadFields(mesh, objects, ssymtFlds);
539 
541  if (!noFields) ReadFields(mesh, objects, stFlds);
542 
543 
544 
545 
546  // Creating (if necessary) faceZones
547  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
548 
549  for (const faceSelection& selector : selectors)
550  {
551  const word& name = selector.name();
552 
553  if (mesh.faceZones().findZoneID(name) == -1)
554  {
556  const label zoneID = mesh.faceZones().size();
557 
559  }
560  }
561 
562 
563  // Select faces
564  // ~~~~~~~~~~~~
565 
566  //- Per face zoneID it is in and flip status.
567  labelList faceToZoneID(mesh.nFaces(), -1);
568  boolList faceToFlip(mesh.nFaces(), false);
569  for (faceSelection& selector : selectors)
570  {
571  const word& name = selector.name();
572  const label zoneID = mesh.faceZones().findZoneID(name);
573 
574  selector.select(zoneID, faceToZoneID, faceToFlip);
575  }
576 
577  // Add faces to faceZones
578  labelList nFaces(mesh.faceZones().size(), Zero);
579  for (const label zoneID : faceToZoneID)
580  {
581  if (zoneID != -1)
582  {
583  ++nFaces[zoneID];
584  }
585  }
586 
587  for (const faceSelection& selector : selectors)
588  {
589  const word& name = selector.name();
590  const label zoneID = mesh.faceZones().findZoneID(name);
591 
592  label& n = nFaces[zoneID];
593  labelList addr(n);
594  boolList flip(n);
595  n = 0;
596  forAll(faceToZoneID, facei)
597  {
598  label zone = faceToZoneID[facei];
599  if (zone == zoneID)
600  {
601  addr[n] = facei;
602  flip[n] = faceToFlip[facei];
603  ++n;
604  }
605  }
606 
607  Info<< "Created zone " << name
608  << " at index " << zoneID
609  << " with " << returnReduce(n, sumOp<label>()) << " faces" << endl;
610 
611  mesh.faceZones().set
612  (
613  zoneID,
614  new faceZone
615  (
616  name,
617  std::move(addr),
618  std::move(flip),
619  zoneID,
620  mesh.faceZones()
621  )
622  );
623  }
624 
625 
626 
627  // Count patches to add
628  // ~~~~~~~~~~~~~~~~~~~~
629  wordHashSet bafflePatches;
630  {
631  for (const faceSelection& selector : selectors)
632  {
633  const dictionary& dict = selector.dict();
634  const word& groupName = selector.name();
635  const dictionary* patchesDictPtr = dict.findDict("patches");
636 
637  if (patchesDictPtr)
638  {
639  for (const entry& dEntry : *patchesDictPtr)
640  {
641  const word patchName(dEntry.dict().get<word>("name"));
642 
643  bafflePatches.insert(patchName);
644  }
645  }
646  else
647  {
648  const word masterName(groupName + "_master");
649  bafflePatches.insert(masterName);
650 
651  const word slaveName(groupName + "_slave");
652  bafflePatches.insert(slaveName);
653  }
654  }
655  }
656 
657 
658  // Create baffles
659  // ~~~~~~~~~~~~~~
660  // Is done in multiple steps
661  // - create patches with 'calculated' patchFields
662  // - move faces into these patches
663  // - change the patchFields to the wanted type
664  // This order is done so e.g. fixedJump works:
665  // - you cannot create patchfields at the same time as patches since
666  // they do an evaluate upon construction
667  // - you want to create the patchField only after you have faces
668  // so you don't get the 'create-from-nothing' mapping problem.
669 
670 
671  // Pass 1: add patches
672  // ~~~~~~~~~~~~~~~~~~~
673 
674  // wordHashSet addedPatches;
675  {
677  for (const faceSelection& selector : selectors)
678  {
679  const dictionary& dict = selector.dict();
680  const word& groupName = selector.name();
681  const dictionary* patchesDictPtr = dict.findDict("patches");
682 
683  if (patchesDictPtr)
684  {
685  for (const entry& dEntry : *patchesDictPtr)
686  {
687  const dictionary& dict = dEntry.dict();
688 
689  const word patchName(dict.get<word>("name"));
690 
691  if (pbm.findPatchID(patchName) == -1)
692  {
693  dictionary patchDict = dict;
694  patchDict.set("nFaces", 0);
695  patchDict.set("startFace", 0);
696 
697  // Note: do not set coupleGroup if constructed from
698  // baffles so you have freedom specifying it
699  // yourself.
700  //patchDict.set("coupleGroup", groupName);
701 
702  addPatch(mesh, patchName, groupName, patchDict);
703  }
704  else
705  {
706  Info<< "Patch '" << patchName
707  << "' already exists. Only "
708  << "moving patch faces - type will remain the same"
709  << endl;
710  }
711  }
712  }
713  else
714  {
715  const dictionary& patchSource = dict.subDict("patchPairs");
716  const word masterName(groupName + "_master");
717  const word slaveName(groupName + "_slave");
718 
719  word groupNameMaster = groupName;
720  word groupNameSlave = groupName;
721 
722 
723  dictionary patchDictMaster(patchSource);
724  patchDictMaster.set("nFaces", 0);
725  patchDictMaster.set("startFace", 0);
726  patchDictMaster.set("coupleGroup", groupName);
727 
728  dictionary patchDictSlave(patchDictMaster);
729 
730  // Note: This is added for the particular case where we want
731  // master and slave in different groupNames
732  // (ie 3D thermal baffles)
733 
734  const bool sameGroup =
735  patchSource.getOrDefault("sameGroup", true);
736 
737  if (!sameGroup)
738  {
739  groupNameMaster = groupName + "Group_master";
740  groupNameSlave = groupName + "Group_slave";
741  patchDictMaster.set("coupleGroup", groupNameMaster);
742  patchDictSlave.set("coupleGroup", groupNameSlave);
743  }
744 
745  addPatch(mesh, masterName, groupNameMaster, patchDictMaster);
746  addPatch(mesh, slaveName, groupNameSlave, patchDictSlave);
747  }
748  }
749  }
750 
751 
752  // Make sure patches and zoneFaces are synchronised across couples
755 
756 
757 
758  // Mesh change container
759  polyTopoChange meshMod(mesh);
760 
762 
763 
764  // Do the actual changes. Note:
765  // - loop in incrementing face order (not necessary if faceZone ordered).
766  // Preserves any existing ordering on patch faces.
767  // - two passes, do non-flip faces first and flip faces second. This
768  // guarantees that when e.g. creating a cyclic all faces from one
769  // side come first and faces from the other side next.
770 
771  // Whether first use of face (modify) or consecutive (add)
772  bitSet modifiedFace(mesh.nFaces());
773  label nModified = 0;
774 
775  for (const faceSelection& selector : selectors)
776  {
777  const word& groupName = selector.name();
778  const dictionary& dict = selector.dict();
779  const dictionary* patchesDictPtr = dict.findDict("patches");
780 
781  const label zoneID = mesh.faceZones().findZoneID(groupName);
782  const faceZone& fZone = mesh.faceZones()[zoneID];
783 
784  DynamicList<label> newMasterPatches;
785  DynamicList<label> newSlavePatches;
786 
787  if (patchesDictPtr)
788  {
789  bool master = true;
790 
791  for (const entry& dEntry : *patchesDictPtr)
792  {
793  const word patchName(dEntry.dict().get<word>("name"));
794 
795  const label patchi = pbm.findPatchID(patchName);
796 
797  if (master)
798  {
799  newMasterPatches.push_back(patchi);
800  }
801  else
802  {
803  newSlavePatches.push_back(patchi);
804  }
805  master = !master;
806  }
807  }
808  else
809  {
810  const word masterName(groupName + "_master");
811  newMasterPatches.push_back(pbm.findPatchID(masterName));
812 
813  const word slaveName(groupName + "_slave");
814  newSlavePatches.push_back(pbm.findPatchID(slaveName));
815  }
816 
817 
818  createFaces
819  (
820  internalFacesOnly,
821  mesh,
822  fZone,
823  newMasterPatches,
824  newSlavePatches,
825  meshMod,
826  modifiedFace,
827  nModified
828  );
829  }
830 
831 
832  Info<< "Converted " << returnReduce(nModified, sumOp<label>())
833  << " faces into boundary faces in patches "
834  << bafflePatches.sortedToc() << nl << endl;
835 
836  if (!overwrite)
837  {
838  ++runTime;
839  }
840 
841  // Change the mesh. Change points directly (no inflation).
842  autoPtr<mapPolyMesh> map = meshMod.changeMesh(mesh, false);
843 
844  // Update fields
845  mesh.updateMesh(map());
846 
847 
848  // Correct boundary faces mapped-out-of-nothing.
849  // This is just a hack to correct the value field.
850  {
851  fvMeshMapper mapper(mesh, map());
852  bool hasWarned = false;
853 
854  for (const word& patchName : bafflePatches)
855  {
856  label patchi = mesh.boundaryMesh().findPatchID(patchName);
857 
858  const fvPatchMapper& pm = mapper.boundaryMap()[patchi];
859 
860  if (pm.sizeBeforeMapping() == 0)
861  {
862  if (!hasWarned)
863  {
864  hasWarned = true;
866  << "Setting field on boundary faces to zero." << endl
867  << "You might have to edit these fields." << endl;
868  }
869 
871  }
872  }
873  }
874 
875 
876  // Pass 2: change patchFields
877  // ~~~~~~~~~~~~~~~~~~~~~~~~~~
878 
879  {
881  for (const faceSelection& selector : selectors)
882  {
883  const dictionary& dict = selector.dict();
884  const word& groupName = selector.name();
885  const dictionary* patchesDictPtr = dict.findDict("patches");
886 
887  if (patchesDictPtr)
888  {
889  for (const entry& dEntry : *patchesDictPtr)
890  {
891  const dictionary& dict = dEntry.dict();
892  const word patchName(dict.get<word>("name"));
893  const label patchi = pbm.findPatchID(patchName);
894 
895  const dictionary* patchFieldsDictPtr
896  = dict.findDict("patchFields");
897 
898  if (patchFieldsDictPtr)
899  {
900  fvMeshTools::setPatchFields
901  (
902  mesh,
903  patchi,
904  *patchFieldsDictPtr
905  );
906  }
907  }
908  }
909  else
910  {
911  const dictionary& patchSource = dict.subDict("patchPairs");
912 
913  const bool sameGroup =
914  patchSource.getOrDefault("sameGroup", true);
915 
916  const dictionary* patchFieldsDictPtr
917  = patchSource.findDict("patchFields");
918 
919  if (patchFieldsDictPtr)
920  {
921  if (sameGroup)
922  {
923  // Make a copy
924  dictionary patchFieldsDict(*patchFieldsDictPtr);
925 
926  // Add coupleGroup to all entries
927  for (entry& dEntry : patchFieldsDict)
928  {
929  if (dEntry.isDict())
930  {
931  dictionary& dict = dEntry.dict();
932  dict.set("coupleGroup", groupName);
933  }
934  }
935 
936  for
937  (
938  const label patchi
939  : pbm.groupPatchIDs()[groupName]
940  )
941  {
942  fvMeshTools::setPatchFields
943  (
944  mesh,
945  patchi,
946  patchFieldsDict
947  );
948  }
949  }
950  else
951  {
952  const word masterPatchName(groupName + "_master");
953  const word slavePatchName(groupName + "_slave");
954 
955  label patchiMaster = pbm.findPatchID(masterPatchName);
956  label patchiSlave = pbm.findPatchID(slavePatchName);
957 
958  fvMeshTools::setPatchFields
959  (
960  mesh,
961  patchiMaster,
962  *patchFieldsDictPtr
963  );
964 
965  fvMeshTools::setPatchFields
966  (
967  mesh,
968  patchiSlave,
969  *patchFieldsDictPtr
970  );
971  }
972  }
973  }
974  }
975  }
976 
977 
978  // Move mesh (since morphing might not do this)
979  if (map().hasMotionPoints())
980  {
981  mesh.movePoints(map().preMotionPoints());
982  }
983 
984 
985  // Remove any now zero-sized patches
986  filterPatches(mesh, bafflePatches);
987 
988 
989  if (overwrite)
990  {
991  mesh.setInstance(oldInstance);
992  }
993 
994  Info<< "Writing mesh to " << runTime.timeName() << endl;
995 
996  mesh.write();
999 
1000  Info<< "End\n" << endl;
1001 
1002  return 0;
1003 }
1004 
1005 
1006 // ************************************************************************* //
label findPatchID(const word &patchName, const bool allowNotFound=true) const
Find patch index given a name, return -1 if not found.
Foam::surfaceFields.
bool checkParallelSync(const bool report=false) const
Check whether all procs have all zones and in same order.
Definition: ZoneMesh.C:743
virtual label sizeBeforeMapping() const
Return size of field before mapping.
const polyBoundaryMesh & pbm
dictionary dict
static void noFunctionObjects(bool addWithOption=false)
Remove &#39;-noFunctionObjects&#39; option and ignore any occurrences.
Definition: argList.C:547
static void reorderPatches(fvMesh &, const labelList &oldToNew, const label nPatches, const bool validBoundary)
Reorder and remove trailing patches.
Definition: fvMeshTools.C:317
static void addNote(const string &note)
Add extra notes for the usage information.
Definition: argList.C:462
static void removeFiles(const polyMesh &mesh)
Helper: remove all procAddressing files from mesh instance.
void clearAddressing()
Clear addressing.
Definition: ZoneMesh.C:700
void set(const bitSet &bitset)
Set specified bits from another bitset.
Definition: bitSetI.H:571
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
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1117
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:162
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:49
bool coupled(solutionDict.getOrDefault("coupledEnergyField", false))
const word dictName("faMeshDefinition")
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:487
virtual void movePoints(const pointField &)
Move points, returns volumes swept by faces in motion.
Definition: fvMesh.C:885
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
T & emplace_back(Args &&... args)
Construct and append an element to the end of the list, return reference to the new list element...
Definition: PtrListI.H:97
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
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:414
static void zeroPatchFields(fvMesh &mesh, const label patchi)
Change patchField to zero on registered fields.
Definition: fvMeshTools.C:223
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:848
virtual void updateMesh(const mapPolyMesh &mpm)
Update mesh corresponding to the given map.
Definition: fvMesh.C:961
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:52
bool contains(const Key &key) const
True if hashed key is contained (found) in table.
Definition: HashTableI.H:72
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
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition: polyMesh.H:584
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:38
label size() const noexcept
The number of entries in the list.
Definition: UPtrListI.H:113
Face selection method for createBaffles.
Definition: faceSelection.H:54
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1111
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:385
label nInternalFaces() const noexcept
Number of internal faces.
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1098
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:159
virtual bool write(const bool writeOnProc=true) const
Write mesh using IO settings from time.
Definition: fvMesh.C:1069
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO...
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:770
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 push_back(const T &val)
Copy append an element to the end of this list.
Definition: DynamicListI.H:518
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.
const word & name() const noexcept
The zone name.
static const word & calculatedType() noexcept
The type name for calculated patch fields.
Definition: fvPatchField.H:201
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...
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
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
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.
const dictionary * findDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary pointer if present (and a sub-dictionary) otherwise return nullptr...
Definition: dictionaryI.H:120
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:133