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  pp.addGroup(groupName);
91 
92 
93  // Add patch, create calculated everywhere
95  (
96  mesh,
97  pp,
98  dictionary(), // do not set specialised patchFields
100  true // parallel sync'ed addition
101  );
102  }
103  else
104  {
105  Info<< "Patch '" << patchName
106  << "' already exists. Only "
107  << "moving patch faces - type will remain the same"
108  << endl;
109  }
110 
111  return pbm.findPatchID(patchName);
112 }
113 
114 
115 // Filter out the empty patches.
116 void filterPatches(fvMesh& mesh, const wordHashSet& addedPatchNames)
117 {
118  // Remove any zero-sized ones. Assumes
119  // - processor patches are already only there if needed
120  // - all other patches are available on all processors
121  // - but coupled ones might still be needed, even if zero-size
122  // (e.g. processorCyclic)
123  // See also logic in createPatch.
125 
126  labelList oldToNew(pbm.size(), -1);
127  label newPatchi = 0;
128  forAll(pbm, patchi)
129  {
130  const polyPatch& pp = pbm[patchi];
131 
132  if (!isA<processorPolyPatch>(pp))
133  {
134  if
135  (
136  isA<coupledPolyPatch>(pp)
137  || returnReduceOr(pp.size())
138  || addedPatchNames.contains(pp.name())
139  )
140  {
141  // Coupled (and unknown size) or uncoupled and used
142  oldToNew[patchi] = newPatchi++;
143  }
144  }
145  }
146 
147  forAll(pbm, patchi)
148  {
149  const polyPatch& pp = pbm[patchi];
150 
151  if (isA<processorPolyPatch>(pp))
152  {
153  oldToNew[patchi] = newPatchi++;
154  }
155  }
156 
157 
158  const label nKeepPatches = newPatchi;
159 
160  // Shuffle unused ones to end
161  if (nKeepPatches != pbm.size())
162  {
163  Info<< endl
164  << "Removing zero-sized patches:" << endl << incrIndent;
165 
166  forAll(oldToNew, patchi)
167  {
168  if (oldToNew[patchi] == -1)
169  {
170  Info<< indent << pbm[patchi].name()
171  << " type " << pbm[patchi].type()
172  << " at position " << patchi << endl;
173  oldToNew[patchi] = newPatchi++;
174  }
175  }
176  Info<< decrIndent;
177 
178  fvMeshTools::reorderPatches(mesh, oldToNew, nKeepPatches, true);
179  Info<< endl;
180  }
181 }
182 
183 
184 void modifyOrAddFace
185 (
186  polyTopoChange& meshMod,
187  const face& f,
188  const label facei,
189  const label own,
190  const bool flipFaceFlux,
191  const label newPatchi,
192  const label zoneID,
193  const bool zoneFlip,
194 
195  bitSet& modifiedFace
196 )
197 {
198  if (modifiedFace.set(facei))
199  {
200  // First usage of face. Modify.
201  meshMod.setAction
202  (
204  (
205  f, // modified face
206  facei, // label of face
207  own, // owner
208  -1, // neighbour
209  flipFaceFlux, // face flip
210  newPatchi, // patch for face
211  false, // remove from zone
212  zoneID, // zone for face
213  zoneFlip // face flip in zone
214  )
215  );
216  }
217  else
218  {
219  // Second or more usage of face. Add.
220  meshMod.setAction
221  (
223  (
224  f, // modified face
225  own, // owner
226  -1, // neighbour
227  -1, // master point
228  -1, // master edge
229  facei, // master face
230  flipFaceFlux, // face flip
231  newPatchi, // patch for face
232  zoneID, // zone for face
233  zoneFlip // face flip in zone
234  )
235  );
236  }
237 }
238 
239 
240 // Create faces for fZone faces. Usually newMasterPatches, newSlavePatches
241 // only size one but can be more for duplicate baffle sets
242 void createFaces
243 (
244  const bool internalFacesOnly,
245  const fvMesh& mesh,
246  const faceZone& fZone,
247  const labelUList& newMasterPatches,
248  const labelUList& newSlavePatches,
249  polyTopoChange& meshMod,
250  bitSet& modifiedFace,
251  label& nModified
252 )
253 {
255 
256  forAll(newMasterPatches, i)
257  {
258  const label newMasterPatchi = newMasterPatches[i];
259  const label newSlavePatchi = newSlavePatches[i];
260 
261  // Pass 1. Do selected side of zone
262  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
263 
264  for (label facei = 0; facei < mesh.nInternalFaces(); ++facei)
265  {
266  const label zoneFacei = fZone.whichFace(facei);
267 
268  if (zoneFacei != -1)
269  {
270  if (!fZone.flipMap()[zoneFacei])
271  {
272  // Use owner side of face
273  modifyOrAddFace
274  (
275  meshMod,
276  mesh.faces()[facei], // modified face
277  facei, // label of face
278  mesh.faceOwner()[facei],// owner
279  false, // face flip
280  newMasterPatchi, // patch for face
281  fZone.index(), // zone for face
282  false, // face flip in zone
283  modifiedFace // modify or add status
284  );
285  }
286  else
287  {
288  // Use neighbour side of face.
289  // To keep faceZone pointing out of original neighbour
290  // we don't need to set faceFlip since that cell
291  // now becomes the owner
292  modifyOrAddFace
293  (
294  meshMod,
295  mesh.faces()[facei].reverseFace(), // modified face
296  facei, // label of face
297  mesh.faceNeighbour()[facei],// owner
298  true, // face flip
299  newMasterPatchi, // patch for face
300  fZone.index(), // zone for face
301  false, // face flip in zone
302  modifiedFace // modify or add status
303  );
304  }
305 
306  ++nModified;
307  }
308  }
309 
310 
311  // Pass 2. Do other side of zone
312  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
313 
314  for (label facei = 0; facei < mesh.nInternalFaces(); ++facei)
315  {
316  const label zoneFacei = fZone.whichFace(facei);
317 
318  if (zoneFacei != -1)
319  {
320  if (!fZone.flipMap()[zoneFacei])
321  {
322  // Use neighbour side of face
323  modifyOrAddFace
324  (
325  meshMod,
326  mesh.faces()[facei].reverseFace(), // modified face
327  facei, // label of face
328  mesh.faceNeighbour()[facei], // owner
329  true, // face flip
330  newSlavePatchi, // patch for face
331  fZone.index(), // zone for face
332  true, // face flip in zone
333  modifiedFace // modify or add
334  );
335  }
336  else
337  {
338  // Use owner side of face
339  modifyOrAddFace
340  (
341  meshMod,
342  mesh.faces()[facei], // modified face
343  facei, // label of face
344  mesh.faceOwner()[facei],// owner
345  false, // face flip
346  newSlavePatchi, // patch for face
347  fZone.index(), // zone for face
348  true, // face flip in zone
349  modifiedFace // modify or add status
350  );
351  }
352  }
353  }
354 
355 
356  // Modify any boundary faces
357  // ~~~~~~~~~~~~~~~~~~~~~~~~~
358 
359  // Normal boundary:
360  // - move to new patch. Might already be back-to-back baffle
361  // you want to add cyclic to. Do warn though.
362  //
363  // Processor boundary:
364  // - do not move to cyclic
365  // - add normal patches though.
366 
367  // For warning once per patch.
368  labelHashSet patchWarned;
369 
370  forAll(pbm, patchi)
371  {
372  const polyPatch& pp = pbm[patchi];
373  const bool isCoupled = pp.coupled();
374 
375  if
376  (
377  isCoupled
378  && (
379  pbm[newMasterPatchi].coupled()
380  || pbm[newSlavePatchi].coupled()
381  )
382  )
383  {
384  // Do not allow coupled faces to be moved to different
385  // coupled patches.
386  }
387  else if (isCoupled || !internalFacesOnly)
388  {
389  forAll(pp, i)
390  {
391  label facei = pp.start()+i;
392 
393  label zoneFacei = fZone.whichFace(facei);
394 
395  if (zoneFacei != -1)
396  {
397  if (!isCoupled && patchWarned.insert(patchi))
398  {
400  << "Found boundary face (patch " << pp.name()
401  << ") in faceZone " << fZone.name()
402  << " to convert to baffle patches "
403  << pbm[newMasterPatchi].name() << "/"
404  << pbm[newSlavePatchi].name() << nl
405  << " Set internalFacesOnly to true in the"
406  << " createBaffles control dictionary if you"
407  << " don't wish to convert boundary faces."
408  << endl;
409  }
410 
411  modifyOrAddFace
412  (
413  meshMod,
414  mesh.faces()[facei], // modified face
415  facei, // label of face
416  mesh.faceOwner()[facei], // owner
417  false, // face flip
418  fZone.flipMap()[zoneFacei]
419  ? newSlavePatchi
420  : newMasterPatchi, // patch for face
421  fZone.index(), // zone for face
422  fZone.flipMap()[zoneFacei], // face flip in zone
423  modifiedFace // modify or add
424  );
425 
426  ++nModified;
427  }
428  }
429  }
430  }
431  }
432 }
433 
434 
435 int main(int argc, char *argv[])
436 {
438  (
439  "Makes internal faces into boundary faces.\n"
440  "Does not duplicate points."
441  );
442 
443  argList::addOption("dict", "file", "Alternative createBafflesDict");
444  #include "addOverwriteOption.H"
445  #include "addRegionOption.H"
446 
447  argList::noFunctionObjects(); // Never use function objects
448 
449  #include "setRootCase.H"
450  #include "createTime.H"
451  #include "createNamedMesh.H"
452 
453 
454  const bool overwrite = args.found("overwrite");
455 
456  const word oldInstance = mesh.pointsInstance();
457 
458  const word dictName("createBafflesDict");
459  #include "setSystemMeshDictionaryIO.H"
460 
461  bool internalFacesOnly(false);
462 
463  bool noFields(false);
464 
465  PtrList<faceSelection> selectors;
466  {
467  Info<< "Reading baffle criteria from " << dictIO.name() << nl << endl;
469 
470  internalFacesOnly = dict.get<bool>("internalFacesOnly");
471  noFields = dict.getOrDefault("noFields", false);
472 
473  const dictionary& selectionsDict = dict.subDict("baffles");
474 
475  selectors.resize(selectionsDict.size());
476 
477  label nselect = 0;
478  for (const entry& dEntry : selectionsDict)
479  {
480  if (dEntry.isDict())
481  {
482  selectors.set
483  (
484  nselect,
485  faceSelection::New(dEntry.keyword(), mesh, dEntry.dict())
486  );
487 
488  ++nselect;
489  }
490  }
491 
492  selectors.resize(nselect);
493  }
494 
495 
496  if (internalFacesOnly)
497  {
498  Info<< "Not converting faces on non-coupled patches." << nl << endl;
499  }
500 
501 
502  // Read objects in time directory
503  IOobjectList objects(mesh, runTime.timeName());
504 
505  // Read vol fields.
506  Info<< "Reading geometric fields" << nl << endl;
508  if (!noFields) ReadFields(mesh, objects, vsFlds);
509 
511  if (!noFields) ReadFields(mesh, objects, vvFlds);
512 
514  if (!noFields) ReadFields(mesh, objects, vstFlds);
515 
516  PtrList<volSymmTensorField> vsymtFlds;
517  if (!noFields) ReadFields(mesh, objects, vsymtFlds);
518 
520  if (!noFields) ReadFields(mesh, objects, vtFlds);
521 
522  // Read surface fields.
523 
525  if (!noFields) ReadFields(mesh, objects, ssFlds);
526 
528  if (!noFields) ReadFields(mesh, objects, svFlds);
529 
531  if (!noFields) ReadFields(mesh, objects, sstFlds);
532 
534  if (!noFields) ReadFields(mesh, objects, ssymtFlds);
535 
537  if (!noFields) ReadFields(mesh, objects, stFlds);
538 
539 
540 
541 
542  // Creating (if necessary) faceZones
543  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
544 
545  for (const faceSelection& selector : selectors)
546  {
547  const word& name = selector.name();
548 
549  if (mesh.faceZones().findZoneID(name) == -1)
550  {
552  const label zoneID = mesh.faceZones().size();
553 
555  }
556  }
557 
558 
559  // Select faces
560  // ~~~~~~~~~~~~
561 
562  //- Per face zoneID it is in and flip status.
563  labelList faceToZoneID(mesh.nFaces(), -1);
564  boolList faceToFlip(mesh.nFaces(), false);
565  for (faceSelection& selector : selectors)
566  {
567  const word& name = selector.name();
568  const label zoneID = mesh.faceZones().findZoneID(name);
569 
570  selector.select(zoneID, faceToZoneID, faceToFlip);
571  }
572 
573  // Add faces to faceZones
574  labelList nFaces(mesh.faceZones().size(), Zero);
575  for (const label zoneID : faceToZoneID)
576  {
577  if (zoneID != -1)
578  {
579  ++nFaces[zoneID];
580  }
581  }
582 
583  for (const faceSelection& selector : selectors)
584  {
585  const word& name = selector.name();
586  const label zoneID = mesh.faceZones().findZoneID(name);
587 
588  label& n = nFaces[zoneID];
589  labelList addr(n);
590  boolList flip(n);
591  n = 0;
592  forAll(faceToZoneID, facei)
593  {
594  label zone = faceToZoneID[facei];
595  if (zone == zoneID)
596  {
597  addr[n] = facei;
598  flip[n] = faceToFlip[facei];
599  ++n;
600  }
601  }
602 
603  Info<< "Created zone " << name
604  << " at index " << zoneID
605  << " with " << returnReduce(n, sumOp<label>()) << " faces" << endl;
606 
607  mesh.faceZones().set
608  (
609  zoneID,
610  new faceZone
611  (
612  name,
613  std::move(addr),
614  std::move(flip),
615  zoneID,
616  mesh.faceZones()
617  )
618  );
619  }
620 
621 
622 
623  // Count patches to add
624  // ~~~~~~~~~~~~~~~~~~~~
625  wordHashSet bafflePatches;
626  {
627  for (const faceSelection& selector : selectors)
628  {
629  const dictionary& dict = selector.dict();
630  const word& groupName = selector.name();
631  const dictionary* patchesDictPtr = dict.findDict("patches");
632 
633  if (patchesDictPtr)
634  {
635  for (const entry& dEntry : *patchesDictPtr)
636  {
637  const word patchName(dEntry.dict().get<word>("name"));
638 
639  bafflePatches.insert(patchName);
640  }
641  }
642  else
643  {
644  const word masterName(groupName + "_master");
645  bafflePatches.insert(masterName);
646 
647  const word slaveName(groupName + "_slave");
648  bafflePatches.insert(slaveName);
649  }
650  }
651  }
652 
653 
654  // Create baffles
655  // ~~~~~~~~~~~~~~
656  // Is done in multiple steps
657  // - create patches with 'calculated' patchFields
658  // - move faces into these patches
659  // - change the patchFields to the wanted type
660  // This order is done so e.g. fixedJump works:
661  // - you cannot create patchfields at the same time as patches since
662  // they do an evaluate upon construction
663  // - you want to create the patchField only after you have faces
664  // so you don't get the 'create-from-nothing' mapping problem.
665 
666 
667  // Pass 1: add patches
668  // ~~~~~~~~~~~~~~~~~~~
669 
670  // wordHashSet addedPatches;
671  {
673  for (const faceSelection& selector : selectors)
674  {
675  const dictionary& dict = selector.dict();
676  const word& groupName = selector.name();
677  const dictionary* patchesDictPtr = dict.findDict("patches");
678 
679  if (patchesDictPtr)
680  {
681  for (const entry& dEntry : *patchesDictPtr)
682  {
683  const dictionary& dict = dEntry.dict();
684 
685  const word patchName(dict.get<word>("name"));
686 
687  if (pbm.findPatchID(patchName) == -1)
688  {
689  dictionary patchDict = dict;
690  patchDict.set("nFaces", 0);
691  patchDict.set("startFace", 0);
692 
693  // Note: do not set coupleGroup if constructed from
694  // baffles so you have freedom specifying it
695  // yourself.
696  //patchDict.set("coupleGroup", groupName);
697 
698  addPatch(mesh, patchName, groupName, patchDict);
699  }
700  else
701  {
702  Info<< "Patch '" << patchName
703  << "' already exists. Only "
704  << "moving patch faces - type will remain the same"
705  << endl;
706  }
707  }
708  }
709  else
710  {
711  const dictionary& patchSource = dict.subDict("patchPairs");
712  const word masterName(groupName + "_master");
713  const word slaveName(groupName + "_slave");
714 
715  word groupNameMaster = groupName;
716  word groupNameSlave = groupName;
717 
718 
719  dictionary patchDictMaster(patchSource);
720  patchDictMaster.set("nFaces", 0);
721  patchDictMaster.set("startFace", 0);
722  patchDictMaster.set("coupleGroup", groupName);
723 
724  dictionary patchDictSlave(patchDictMaster);
725 
726  // Note: This is added for the particular case where we want
727  // master and slave in different groupNames
728  // (ie 3D thermal baffles)
729 
730  const bool sameGroup =
731  patchSource.getOrDefault("sameGroup", true);
732 
733  if (!sameGroup)
734  {
735  groupNameMaster = groupName + "Group_master";
736  groupNameSlave = groupName + "Group_slave";
737  patchDictMaster.set("coupleGroup", groupNameMaster);
738  patchDictSlave.set("coupleGroup", groupNameSlave);
739  }
740 
741  addPatch(mesh, masterName, groupNameMaster, patchDictMaster);
742  addPatch(mesh, slaveName, groupNameSlave, patchDictSlave);
743  }
744  }
745  }
746 
747 
748  // Make sure patches and zoneFaces are synchronised across couples
751 
752 
753 
754  // Mesh change container
755  polyTopoChange meshMod(mesh);
756 
758 
759 
760  // Do the actual changes. Note:
761  // - loop in incrementing face order (not necessary if faceZone ordered).
762  // Preserves any existing ordering on patch faces.
763  // - two passes, do non-flip faces first and flip faces second. This
764  // guarantees that when e.g. creating a cyclic all faces from one
765  // side come first and faces from the other side next.
766 
767  // Whether first use of face (modify) or consecutive (add)
768  bitSet modifiedFace(mesh.nFaces());
769  label nModified = 0;
770 
771  for (const faceSelection& selector : selectors)
772  {
773  const word& groupName = selector.name();
774  const dictionary& dict = selector.dict();
775  const dictionary* patchesDictPtr = dict.findDict("patches");
776 
777  const label zoneID = mesh.faceZones().findZoneID(groupName);
778  const faceZone& fZone = mesh.faceZones()[zoneID];
779 
780  DynamicList<label> newMasterPatches;
781  DynamicList<label> newSlavePatches;
782 
783  if (patchesDictPtr)
784  {
785  bool master = true;
786 
787  for (const entry& dEntry : *patchesDictPtr)
788  {
789  const word patchName(dEntry.dict().get<word>("name"));
790 
791  const label patchi = pbm.findPatchID(patchName);
792 
793  if (master)
794  {
795  newMasterPatches.push_back(patchi);
796  }
797  else
798  {
799  newSlavePatches.push_back(patchi);
800  }
801  master = !master;
802  }
803  }
804  else
805  {
806  const word masterName(groupName + "_master");
807  newMasterPatches.push_back(pbm.findPatchID(masterName));
808 
809  const word slaveName(groupName + "_slave");
810  newSlavePatches.push_back(pbm.findPatchID(slaveName));
811  }
812 
813 
814  createFaces
815  (
816  internalFacesOnly,
817  mesh,
818  fZone,
819  newMasterPatches,
820  newSlavePatches,
821  meshMod,
822  modifiedFace,
823  nModified
824  );
825  }
826 
827 
828  Info<< "Converted " << returnReduce(nModified, sumOp<label>())
829  << " faces into boundary faces in patches "
830  << bafflePatches.sortedToc() << nl << endl;
831 
832  if (!overwrite)
833  {
834  ++runTime;
835  }
836 
837  // Change the mesh. Change points directly (no inflation).
838  autoPtr<mapPolyMesh> map = meshMod.changeMesh(mesh, false);
839 
840  // Update fields
841  mesh.updateMesh(map());
842 
843 
844  // Correct boundary faces mapped-out-of-nothing.
845  // This is just a hack to correct the value field.
846  {
847  fvMeshMapper mapper(mesh, map());
848  bool hasWarned = false;
849 
850  for (const word& patchName : bafflePatches)
851  {
852  label patchi = mesh.boundaryMesh().findPatchID(patchName);
853 
854  const fvPatchMapper& pm = mapper.boundaryMap()[patchi];
855 
856  if (pm.sizeBeforeMapping() == 0)
857  {
858  if (!hasWarned)
859  {
860  hasWarned = true;
862  << "Setting field on boundary faces to zero." << endl
863  << "You might have to edit these fields." << endl;
864  }
865 
867  }
868  }
869  }
870 
871 
872  // Pass 2: change patchFields
873  // ~~~~~~~~~~~~~~~~~~~~~~~~~~
874 
875  {
877  for (const faceSelection& selector : selectors)
878  {
879  const dictionary& dict = selector.dict();
880  const word& groupName = selector.name();
881  const dictionary* patchesDictPtr = dict.findDict("patches");
882 
883  if (patchesDictPtr)
884  {
885  for (const entry& dEntry : *patchesDictPtr)
886  {
887  const dictionary& dict = dEntry.dict();
888  const word patchName(dict.get<word>("name"));
889  const label patchi = pbm.findPatchID(patchName);
890 
891  const dictionary* patchFieldsDictPtr
892  = dict.findDict("patchFields");
893 
894  if (patchFieldsDictPtr)
895  {
896  fvMeshTools::setPatchFields
897  (
898  mesh,
899  patchi,
900  *patchFieldsDictPtr
901  );
902  }
903  }
904  }
905  else
906  {
907  const dictionary& patchSource = dict.subDict("patchPairs");
908 
909  const bool sameGroup =
910  patchSource.getOrDefault("sameGroup", true);
911 
912  const dictionary* patchFieldsDictPtr
913  = patchSource.findDict("patchFields");
914 
915  if (patchFieldsDictPtr)
916  {
917  if (sameGroup)
918  {
919  // Make a copy
920  dictionary patchFieldsDict(*patchFieldsDictPtr);
921 
922  // Add coupleGroup to all entries
923  for (entry& dEntry : patchFieldsDict)
924  {
925  if (dEntry.isDict())
926  {
927  dictionary& dict = dEntry.dict();
928  dict.set("coupleGroup", groupName);
929  }
930  }
931 
932  for
933  (
934  const label patchi
935  : pbm.groupPatchIDs()[groupName]
936  )
937  {
938  fvMeshTools::setPatchFields
939  (
940  mesh,
941  patchi,
942  patchFieldsDict
943  );
944  }
945  }
946  else
947  {
948  const word masterPatchName(groupName + "_master");
949  const word slavePatchName(groupName + "_slave");
950 
951  label patchiMaster = pbm.findPatchID(masterPatchName);
952  label patchiSlave = pbm.findPatchID(slavePatchName);
953 
954  fvMeshTools::setPatchFields
955  (
956  mesh,
957  patchiMaster,
958  *patchFieldsDictPtr
959  );
960 
961  fvMeshTools::setPatchFields
962  (
963  mesh,
964  patchiSlave,
965  *patchFieldsDictPtr
966  );
967  }
968  }
969  }
970  }
971  }
972 
973 
974  // Move mesh (since morphing might not do this)
975  if (map().hasMotionPoints())
976  {
977  mesh.movePoints(map().preMotionPoints());
978  }
979 
980 
981  // Remove any now zero-sized patches
982  filterPatches(mesh, bafflePatches);
983 
984 
985  if (overwrite)
986  {
987  mesh.setInstance(oldInstance);
988  }
989 
990  Info<< "Writing mesh to " << runTime.timeName() << endl;
991 
992  mesh.write();
995 
996  Info<< "End\n" << endl;
997 
998  return 0;
999 }
1000 
1001 
1002 // ************************************************************************* //
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:959
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:905
void set(const bitSet &bitset)
Set specified bits from another bitset.
Definition: bitSetI.H:498
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:493
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:129
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1122
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:195
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
const word dictName("faMeshDefinition")
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
virtual void movePoints(const pointField &)
Move points, returns volumes swept by faces in motion.
Definition: fvMesh.C:929
Required Classes.
label whichFace(const label globalCellID) const
Helper function to re-direct to zone::localID(...)
Definition: faceZone.C:397
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:232
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:441
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:105
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:421
static void zeroPatchFields(fvMesh &mesh, const label patchi)
Change patchField to zero on registered fields.
Definition: fvMeshTools.C:223
Required Classes.
static void removeFiles(const polyMesh &)
Helper: remove all sets files from mesh instance.
Definition: topoSet.C:693
const fileName & pointsInstance() const
Return the current instance directory for points.
Definition: polyMesh.C:853
virtual void updateMesh(const mapPolyMesh &mpm)
Update mesh corresponding to the given map.
Definition: fvMesh.C:1005
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
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:718
Base class for mesh zones.
Definition: zone.H:59
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition: polyMesh.H:609
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:106
Face selection method for createBaffles.
Definition: faceSelection.H:54
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1116
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:1103
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:1113
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 a time name for the given scalar time value formatted with the given precision.
Definition: Time.C:714
const faceZoneMesh & faceZones() const noexcept
Return face zone mesh.
Definition: polyMesh.H:671
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition: Ostream.H:511
labelList f(nPoints)
void push_back(const T &val)
Copy append an element to the end of this list.
Definition: DynamicListI.H:555
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:78
List< Key > sortedToc() const
The table of contents (the keys) in sorted order.
Definition: HashTable.C:163
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:204
const boolList & flipMap() const noexcept
Return face flip map.
Definition: faceZone.H:367
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:765
A subset of mesh faces organised as a primitive patch.
Definition: faceZone.H:60
Mapping class for a fvPatchField.
Definition: fvPatchMapper.H:54
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
bool coupled
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:502
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 it is a dictionary) otherwise return nullptr...
Definition: dictionaryI.H:124
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127