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-2021 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 "processorMeshes.H"
61 
62 using namespace Foam;
63 
64 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
65 
66 label addPatch
67 (
68  fvMesh& mesh,
69  const word& patchName,
70  const word& groupName,
71  const dictionary& patchDict
72 )
73 {
74  const polyBoundaryMesh& pbm = mesh.boundaryMesh();
75 
76  if (pbm.findPatchID(patchName) == -1)
77  {
78  autoPtr<polyPatch> ppPtr
79  (
81  (
82  patchName,
83  patchDict,
84  0,
85  pbm
86  )
87  );
88  auto& pp = *ppPtr;
89 
90  if (!groupName.empty())
91  {
92  pp.inGroups().appendUniq(groupName);
93  }
94 
95 
96  // Add patch, create calculated everywhere
98  (
99  mesh,
100  pp,
101  dictionary(), // do not set specialised patchFields
103  true // parallel sync'ed addition
104  );
105  }
106  else
107  {
108  Info<< "Patch '" << patchName
109  << "' already exists. Only "
110  << "moving patch faces - type will remain the same"
111  << endl;
112  }
113 
114  return pbm.findPatchID(patchName);
115 }
116 
117 
118 // Filter out the empty patches.
119 void filterPatches(fvMesh& mesh, const wordHashSet& addedPatchNames)
120 {
121  // Remove any zero-sized ones. Assumes
122  // - processor patches are already only there if needed
123  // - all other patches are available on all processors
124  // - but coupled ones might still be needed, even if zero-size
125  // (e.g. processorCyclic)
126  // See also logic in createPatch.
127  const polyBoundaryMesh& pbm = mesh.boundaryMesh();
128 
129  labelList oldToNew(pbm.size(), -1);
130  label newPatchi = 0;
131  forAll(pbm, patchi)
132  {
133  const polyPatch& pp = pbm[patchi];
134 
135  if (!isA<processorPolyPatch>(pp))
136  {
137  if
138  (
139  isA<coupledPolyPatch>(pp)
140  || returnReduce(pp.size(), sumOp<label>())
141  || addedPatchNames.found(pp.name())
142  )
143  {
144  // Coupled (and unknown size) or uncoupled and used
145  oldToNew[patchi] = newPatchi++;
146  }
147  }
148  }
149 
150  forAll(pbm, patchi)
151  {
152  const polyPatch& pp = pbm[patchi];
153 
154  if (isA<processorPolyPatch>(pp))
155  {
156  oldToNew[patchi] = newPatchi++;
157  }
158  }
159 
160 
161  const label nKeepPatches = newPatchi;
162 
163  // Shuffle unused ones to end
164  if (nKeepPatches != pbm.size())
165  {
166  Info<< endl
167  << "Removing zero-sized patches:" << endl << incrIndent;
168 
169  forAll(oldToNew, patchi)
170  {
171  if (oldToNew[patchi] == -1)
172  {
173  Info<< indent << pbm[patchi].name()
174  << " type " << pbm[patchi].type()
175  << " at position " << patchi << endl;
176  oldToNew[patchi] = newPatchi++;
177  }
178  }
179  Info<< decrIndent;
180 
181  fvMeshTools::reorderPatches(mesh, oldToNew, nKeepPatches, true);
182  Info<< endl;
183  }
184 }
185 
186 
187 void modifyOrAddFace
188 (
189  polyTopoChange& meshMod,
190  const face& f,
191  const label facei,
192  const label own,
193  const bool flipFaceFlux,
194  const label newPatchi,
195  const label zoneID,
196  const bool zoneFlip,
197 
198  bitSet& modifiedFace
199 )
200 {
201  if (modifiedFace.set(facei))
202  {
203  // First usage of face. Modify.
204  meshMod.setAction
205  (
207  (
208  f, // modified face
209  facei, // label of face
210  own, // owner
211  -1, // neighbour
212  flipFaceFlux, // face flip
213  newPatchi, // patch for face
214  false, // remove from zone
215  zoneID, // zone for face
216  zoneFlip // face flip in zone
217  )
218  );
219  }
220  else
221  {
222  // Second or more usage of face. Add.
223  meshMod.setAction
224  (
226  (
227  f, // modified face
228  own, // owner
229  -1, // neighbour
230  -1, // master point
231  -1, // master edge
232  facei, // master face
233  flipFaceFlux, // face flip
234  newPatchi, // patch for face
235  zoneID, // zone for face
236  zoneFlip // face flip in zone
237  )
238  );
239  }
240 }
241 
242 
243 // Create faces for fZone faces. Usually newMasterPatches, newSlavePatches
244 // only size one but can be more for duplicate baffle sets
245 void createFaces
246 (
247  const bool internalFacesOnly,
248  const fvMesh& mesh,
249  const faceZone& fZone,
250  const labelList& newMasterPatches,
251  const labelList& newSlavePatches,
252  polyTopoChange& meshMod,
253  bitSet& modifiedFace,
254  label& nModified
255 )
256 {
257  const polyBoundaryMesh& pbm = mesh.boundaryMesh();
258 
259  forAll(newMasterPatches, i)
260  {
261  // Pass 1. Do selected side of zone
262  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
263 
264  for (label facei = 0; facei < mesh.nInternalFaces(); facei++)
265  {
266  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  newMasterPatches[i], // 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  newMasterPatches[i], // 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  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  newSlavePatches[i], // 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  newSlavePatches[i], // 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 
374  const label newMasterPatchi = newMasterPatches[i];
375  const label newSlavePatchi = newSlavePatches[i];
376 
377  if
378  (
379  pp.coupled()
380  && (
381  pbm[newMasterPatchi].coupled()
382  || pbm[newSlavePatchi].coupled()
383  )
384  )
385  {
386  // Do not allow coupled faces to be moved to different
387  // coupled patches.
388  }
389  else if (pp.coupled() || !internalFacesOnly)
390  {
391  forAll(pp, i)
392  {
393  label facei = pp.start()+i;
394 
395  label zoneFacei = fZone.whichFace(facei);
396 
397  if (zoneFacei != -1)
398  {
399  if (patchWarned.insert(patchi))
400  {
402  << "Found boundary face (in patch "
403  << pp.name()
404  << ") in faceZone " << fZone.name()
405  << " to convert to baffle patches "
406  << pbm[newMasterPatchi].name() << "/"
407  << pbm[newSlavePatchi].name()
408  << endl
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  forAll(selectors, selectorI)
550  {
551  const word& name = selectors[selectorI].name();
552 
553  if (mesh.faceZones().findZoneID(name) == -1)
554  {
556  const label zoneID = mesh.faceZones().size();
557 
558  mesh.faceZones().setSize(zoneID+1);
559  mesh.faceZones().set
560  (
561  zoneID,
562  new faceZone(name, labelList(), false, zoneID, mesh.faceZones())
563  );
564  }
565  }
566 
567 
568  // Select faces
569  // ~~~~~~~~~~~~
570 
571  //- Per face zoneID it is in and flip status.
572  labelList faceToZoneID(mesh.nFaces(), -1);
573  boolList faceToFlip(mesh.nFaces(), false);
574  forAll(selectors, selectorI)
575  {
576  const word& name = selectors[selectorI].name();
577  label zoneID = mesh.faceZones().findZoneID(name);
578 
579  selectors[selectorI].select(zoneID, faceToZoneID, faceToFlip);
580  }
581 
582  // Add faces to faceZones
583  labelList nFaces(mesh.faceZones().size(), Zero);
584  forAll(faceToZoneID, facei)
585  {
586  label zoneID = faceToZoneID[facei];
587  if (zoneID != -1)
588  {
589  nFaces[zoneID]++;
590  }
591  }
592 
593  forAll(selectors, selectorI)
594  {
595  const word& name = selectors[selectorI].name();
596  label zoneID = mesh.faceZones().findZoneID(name);
597 
598  label& n = nFaces[zoneID];
599  labelList addr(n);
600  boolList flip(n);
601  n = 0;
602  forAll(faceToZoneID, facei)
603  {
604  label zone = faceToZoneID[facei];
605  if (zone == zoneID)
606  {
607  addr[n] = facei;
608  flip[n] = faceToFlip[facei];
609  ++n;
610  }
611  }
612 
613  Info<< "Created zone " << name
614  << " at index " << zoneID
615  << " with " << returnReduce(n, sumOp<label>()) << " faces" << endl;
616 
617  mesh.faceZones().set
618  (
619  zoneID,
620  new faceZone
621  (
622  name,
623  std::move(addr),
624  std::move(flip),
625  zoneID,
626  mesh.faceZones()
627  )
628  );
629  }
630 
631 
632 
633  // Count patches to add
634  // ~~~~~~~~~~~~~~~~~~~~
635  wordHashSet bafflePatches;
636  {
637  forAll(selectors, selectorI)
638  {
639  const dictionary& dict = selectors[selectorI].dict();
640 
641  if (dict.found("patches"))
642  {
643  for (const entry& dEntry : dict.subDict("patches"))
644  {
645  const word patchName(dEntry.dict().get<word>("name"));
646 
647  bafflePatches.insert(patchName);
648  }
649  }
650  else
651  {
652  const word masterName = selectors[selectorI].name() + "_master";
653  bafflePatches.insert(masterName);
654 
655  const word slaveName = selectors[selectorI].name() + "_slave";
656  bafflePatches.insert(slaveName);
657  }
658  }
659  }
660 
661 
662  // Create baffles
663  // ~~~~~~~~~~~~~~
664  // Is done in multiple steps
665  // - create patches with 'calculated' patchFields
666  // - move faces into these patches
667  // - change the patchFields to the wanted type
668  // This order is done so e.g. fixedJump works:
669  // - you cannot create patchfields at the same time as patches since
670  // they do an evaluate upon construction
671  // - you want to create the patchField only after you have faces
672  // so you don't get the 'create-from-nothing' mapping problem.
673 
674 
675  // Pass 1: add patches
676  // ~~~~~~~~~~~~~~~~~~~
677 
678  // wordHashSet addedPatches;
679  {
680  const polyBoundaryMesh& pbm = mesh.boundaryMesh();
681  forAll(selectors, selectorI)
682  {
683  const dictionary& dict = selectors[selectorI].dict();
684  const word& groupName = selectors[selectorI].name();
685 
686  if (dict.found("patches"))
687  {
688  for (const entry& dEntry : dict.subDict("patches"))
689  {
690  const dictionary& dict = dEntry.dict();
691 
692  const word patchName(dict.get<word>("name"));
693 
694  if (pbm.findPatchID(patchName) == -1)
695  {
696  dictionary patchDict = dict;
697  patchDict.set("nFaces", 0);
698  patchDict.set("startFace", 0);
699 
700  // Note: do not set coupleGroup if constructed from
701  // baffles so you have freedom specifying it
702  // yourself.
703  //patchDict.set("coupleGroup", groupName);
704 
705  addPatch(mesh, patchName, groupName, patchDict);
706  }
707  else
708  {
709  Info<< "Patch '" << patchName
710  << "' already exists. Only "
711  << "moving patch faces - type will remain the same"
712  << endl;
713  }
714  }
715  }
716  else
717  {
718  const dictionary& patchSource = dict.subDict("patchPairs");
719  const word masterName = groupName + "_master";
720  const word slaveName = groupName + "_slave";
721 
722  word groupNameMaster = groupName;
723  word groupNameSlave = groupName;
724 
725 
726  dictionary patchDictMaster(patchSource);
727  patchDictMaster.set("nFaces", 0);
728  patchDictMaster.set("startFace", 0);
729  patchDictMaster.set("coupleGroup", groupName);
730 
731  dictionary patchDictSlave(patchDictMaster);
732 
733  // Note: This is added for the particular case where we want
734  // master and slave in different groupNames
735  // (ie 3D thermal baffles)
736 
737  const bool sameGroup =
738  patchSource.getOrDefault("sameGroup", true);
739 
740  if (!sameGroup)
741  {
742  groupNameMaster = groupName + "Group_master";
743  groupNameSlave = groupName + "Group_slave";
744  patchDictMaster.set("coupleGroup", groupNameMaster);
745  patchDictSlave.set("coupleGroup", groupNameSlave);
746  }
747 
748  addPatch(mesh, masterName, groupNameMaster, patchDictMaster);
749  addPatch(mesh, slaveName, groupNameSlave, patchDictSlave);
750  }
751  }
752  }
753 
754 
755  // Make sure patches and zoneFaces are synchronised across couples
758 
759 
760 
761  // Mesh change container
762  polyTopoChange meshMod(mesh);
763 
764  const polyBoundaryMesh& pbm = mesh.boundaryMesh();
765 
766 
767  // Do the actual changes. Note:
768  // - loop in incrementing face order (not necessary if faceZone ordered).
769  // Preserves any existing ordering on patch faces.
770  // - two passes, do non-flip faces first and flip faces second. This
771  // guarantees that when e.g. creating a cyclic all faces from one
772  // side come first and faces from the other side next.
773 
774  // Whether first use of face (modify) or consecutive (add)
775  bitSet modifiedFace(mesh.nFaces());
776  label nModified = 0;
777 
778  forAll(selectors, selectorI)
779  {
780  const word& name = selectors[selectorI].name();
781  label zoneID = mesh.faceZones().findZoneID(name);
782  const faceZone& fZone = mesh.faceZones()[zoneID];
783 
784  const dictionary& dict = selectors[selectorI].dict();
785 
786  DynamicList<label> newMasterPatches;
787  DynamicList<label> newSlavePatches;
788 
789  if (dict.found("patches"))
790  {
791  bool master = true;
792 
793  for (const entry& dEntry : dict.subDict("patches"))
794  {
795  const word patchName(dEntry.dict().get<word>("name"));
796 
797  const label patchi = pbm.findPatchID(patchName);
798 
799  if (master)
800  {
801  newMasterPatches.append(patchi);
802  }
803  else
804  {
805  newSlavePatches.append(patchi);
806  }
807  master = !master;
808  }
809  }
810  else
811  {
812  const word masterName = selectors[selectorI].name() + "_master";
813  newMasterPatches.append(pbm.findPatchID(masterName));
814 
815  const word slaveName = selectors[selectorI].name() + "_slave";
816  newSlavePatches.append(pbm.findPatchID(slaveName));
817  }
818 
819 
820  createFaces
821  (
822  internalFacesOnly,
823  mesh,
824  fZone,
825  newMasterPatches,
826  newSlavePatches,
827  meshMod,
828  modifiedFace,
829  nModified
830  );
831  }
832 
833 
834  Info<< "Converted " << returnReduce(nModified, sumOp<label>())
835  << " faces into boundary faces in patches "
836  << bafflePatches.sortedToc() << nl << endl;
837 
838  if (!overwrite)
839  {
840  ++runTime;
841  }
842 
843  // Change the mesh. Change points directly (no inflation).
844  autoPtr<mapPolyMesh> map = meshMod.changeMesh(mesh, false);
845 
846  // Update fields
847  mesh.updateMesh(map());
848 
849 
850  // Correct boundary faces mapped-out-of-nothing.
851  // This is just a hack to correct the value field.
852  {
853  fvMeshMapper mapper(mesh, map());
854  bool hasWarned = false;
855 
856  for (const word& patchName : bafflePatches)
857  {
858  label patchi = mesh.boundaryMesh().findPatchID(patchName);
859 
860  const fvPatchMapper& pm = mapper.boundaryMap()[patchi];
861 
862  if (pm.sizeBeforeMapping() == 0)
863  {
864  if (!hasWarned)
865  {
866  hasWarned = true;
868  << "Setting field on boundary faces to zero." << endl
869  << "You might have to edit these fields." << endl;
870  }
871 
873  }
874  }
875  }
876 
877 
878  // Pass 2: change patchFields
879  // ~~~~~~~~~~~~~~~~~~~~~~~~~~
880 
881  {
882  const polyBoundaryMesh& pbm = mesh.boundaryMesh();
883  forAll(selectors, selectorI)
884  {
885  const dictionary& dict = selectors[selectorI].dict();
886  if (dict.found("patches"))
887  {
888  for (const entry& dEntry : dict.subDict("patches"))
889  {
890  const dictionary& dict = dEntry.dict();
891 
892  const word patchName(dict.get<word>("name"));
893 
894  label patchi = pbm.findPatchID(patchName);
895 
896  if (dEntry.dict().found("patchFields"))
897  {
898  const dictionary& patchFieldsDict =
899  dEntry.dict().subDict
900  (
901  "patchFields"
902  );
903 
904  fvMeshTools::setPatchFields
905  (
906  mesh,
907  patchi,
908  patchFieldsDict
909  );
910  }
911  }
912  }
913  else
914  {
915  const dictionary& patchSource = dict.subDict("patchPairs");
916 
917  const bool sameGroup =
918  patchSource.getOrDefault("sameGroup", true);
919 
920  const word& groupName = selectors[selectorI].name();
921 
922  if (patchSource.found("patchFields"))
923  {
924  dictionary patchFieldsDict = patchSource.subDict
925  (
926  "patchFields"
927  );
928 
929  if (sameGroup)
930  {
931  // Add coupleGroup to all entries
932  for (entry& dEntry : patchFieldsDict)
933  {
934  if (dEntry.isDict())
935  {
936  dictionary& dict = dEntry.dict();
937  dict.set("coupleGroup", groupName);
938  }
939  }
940 
941  const labelList& patchIDs =
942  pbm.groupPatchIDs()[groupName];
943 
944  forAll(patchIDs, i)
945  {
946  fvMeshTools::setPatchFields
947  (
948  mesh,
949  patchIDs[i],
950  patchFieldsDict
951  );
952  }
953  }
954  else
955  {
956  const word masterPatchName(groupName + "_master");
957  const word slavePatchName(groupName + "_slave");
958 
959  label patchiMaster = pbm.findPatchID(masterPatchName);
960  label patchiSlave = pbm.findPatchID(slavePatchName);
961 
962  fvMeshTools::setPatchFields
963  (
964  mesh,
965  patchiMaster,
966  patchFieldsDict
967  );
968 
969  fvMeshTools::setPatchFields
970  (
971  mesh,
972  patchiSlave,
973  patchFieldsDict
974  );
975  }
976  }
977  }
978  }
979  }
980 
981 
982  // Move mesh (since morphing might not do this)
983  if (map().hasMotionPoints())
984  {
985  mesh.movePoints(map().preMotionPoints());
986  }
987 
988 
989  // Remove any now zero-sized patches
990  filterPatches(mesh, bafflePatches);
991 
992 
993  if (overwrite)
994  {
995  mesh.setInstance(oldInstance);
996  }
997 
998  Info<< "Writing mesh to " << runTime.timeName() << endl;
999 
1000  mesh.write();
1003 
1004  Info<< "End\n" << endl;
1005 
1006  return 0;
1007 }
1008 
1009 
1010 // ************************************************************************* //
faceSelection.H
Foam::entry
A keyword and a list of tokens is an 'entry'.
Definition: entry.H:63
Foam::fvPatchMapper
Mapping class for a fvPatchField.
Definition: fvPatchMapper.H:53
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:63
volFields.H
Foam::IOdictionary
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:50
runTime
engineTime & runTime
Definition: createEngineTime.H:13
Foam::faceSelection::New
static autoPtr< faceSelection > New(const word &name, const fvMesh &mesh, const dictionary &dict)
Foam::fvMeshTools::zeroPatchFields
static void zeroPatchFields(fvMesh &mesh, const label patchi)
Definition: fvMeshTools.C:236
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:63
Foam::faceZone::flipMap
const boolList & flipMap() const noexcept
Definition: faceZone.H:268
Foam::fvMesh::write
virtual bool write(const bool valid=true) const
Definition: fvMesh.C:1034
Foam::returnReduce
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Definition: PstreamReduceOps.H:88
Foam::bitSet
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition: bitSet.H:59
Foam::polyBoundaryMesh
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
Definition: polyBoundaryMesh.H:59
Foam::ReadFields
wordList ReadFields(const typename GeoMesh::Mesh &mesh, const IOobjectList &objects, PtrList< GeometricField< Type, PatchField, GeoMesh >> &fields, const bool syncPar=true, const bool readOldTime=false)
topoSet.H
Foam::Zero
static constexpr const zero Zero
Definition: zero.H:131
Foam::DynamicList< label >
Foam::dictionary::found
bool found(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryI.H:80
polyAddFace.H
addOverwriteOption.H
dictName
const word dictName("faMeshDefinition")
Foam::polyPatch::coupled
virtual bool coupled() const
Definition: polyPatch.H:373
Foam::polyBoundaryMesh::groupPatchIDs
const HashTable< labelList > & groupPatchIDs() const
Definition: polyBoundaryMesh.C:474
polyTopoChange.H
Foam::Time::timeName
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:773
Foam::polyTopoChange
Direct mesh changes based on v1.3 polyTopoChange syntax.
Definition: polyTopoChange.H:95
Foam::zone
Base class for mesh zones.
Definition: zone.H:59
Foam::argList::addNote
static void addNote(const string &note)
Definition: argList.C:405
Foam::bitSet::set
void set(const bitSet &bitset)
Definition: bitSetI.H:567
Foam::polyModifyFace
Class describing modification of a face.
Definition: polyModifyFace.H:44
Foam::fvMesh::movePoints
virtual tmp< scalarField > movePoints(const pointField &)
Definition: fvMesh.C:858
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Definition: polyMesh.H:440
Foam::endl
Ostream & endl(Ostream &os)
Definition: Ostream.H:381
surfaceFields.H
Foam::surfaceFields.
Foam::PtrList::set
const T * set(const label i) const
Definition: PtrList.H:134
Foam::dictionary::get
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:100
Foam::dictionary::set
entry * set(entry *entryPtr)
Definition: dictionary.C:773
Foam::fvPatchMapper::sizeBeforeMapping
virtual label sizeBeforeMapping() const
Definition: fvPatchMapper.H:126
Foam::HashSet
A HashTable with keys but without contents that is similar to std::unordered_set.
Definition: HashSet.H:73
processorMeshes.H
Foam::incrIndent
Ostream & incrIndent(Ostream &os)
Definition: Ostream.H:352
setSystemMeshDictionaryIO.H
Foam::sumOp
Definition: ops.H:207
forAll
#define forAll(list, i)
Definition: stdFoam.H:349
Foam::polyMesh::pointsInstance
const fileName & pointsInstance() const
Definition: polyMesh.C:839
n
label n
Definition: TABSMDCalcMethod2.H:31
Foam::polyBoundaryMesh::checkParallelSync
bool checkParallelSync(const bool report=false) const
Definition: polyBoundaryMesh.C:965
Foam::argList::noFunctionObjects
static void noFunctionObjects(bool addWithOption=false)
Definition: argList.C:466
fvMeshMapper.H
Foam::patchIdentifier::inGroups
const wordList & inGroups() const noexcept
Definition: patchIdentifier.H:167
Foam::faceZone
A subset of mesh faces organised as a primitive patch.
Definition: faceZone.H:60
Foam::Info
messageStream Info
Foam::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:64
Foam::DynamicList::append
DynamicList< T, SizeMin > & append(const T &val)
Definition: DynamicListI.H:504
argList.H
Foam::polyMesh::faceOwner
virtual const labelList & faceOwner() const
Definition: polyMesh.C:1100
Foam::dictionary::subDict
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionary.C:453
addRegionOption.H
Foam::fvMesh::updateMesh
virtual void updateMesh(const mapPolyMesh &mpm)
Definition: fvMesh.C:946
Foam::polyMesh::faceZones
const faceZoneMesh & faceZones() const noexcept
Definition: polyMesh.H:482
Foam::PtrList
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: List.H:55
Foam::fvMeshTools::reorderPatches
static void reorderPatches(fvMesh &, const labelList &oldToNew, const label nPatches, const bool validBoundary)
Definition: fvMeshTools.C:322
Foam::polyBoundaryMesh::findPatchID
label findPatchID(const word &patchName, const bool allowNotFound=true) const
Definition: polyBoundaryMesh.C:758
Foam::List::appendUniq
label appendUniq(const T &val)
Definition: ListI.H:225
Foam::polyTopoChange::changeMesh
autoPtr< mapPolyMesh > changeMesh(polyMesh &mesh, const labelUList &patchMap, const bool inflate, const bool syncParallel=true, const bool orderCells=false, const bool orderPoints=false)
Definition: polyTopoChange.C:2973
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:119
createNamedMesh.H
Required Variables.
Foam::polyPatch::New
static autoPtr< polyPatch > New(const word &patchType, const word &name, const label size, const label start, const label index, const polyBoundaryMesh &bm)
Definition: polyPatchNew.C:28
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
zoneID
const labelIOList & zoneID
Definition: interpolatedFaces.H:22
dictIO
IOobject dictIO
Definition: setConstantMeshDictionaryIO.H:1
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:81
Foam::calculatedFvPatchField
This boundary condition is not designed to be evaluated; it is assmued that the value is assigned via...
Definition: calculatedFvPatchField.H:62
Foam
Definition: atmBoundaryLayer.C:26
Foam::decrIndent
Ostream & decrIndent(Ostream &os)
Definition: Ostream.H:361
Foam::faceZone::whichFace
label whichFace(const label globalCellID) const
Definition: faceZone.C:365
Foam::ZoneMesh::findZoneID
label findZoneID(const word &zoneName) const
Definition: ZoneMesh.C:512
Foam::indent
Ostream & indent(Ostream &os)
Definition: Ostream.H:343
Foam::polyPatch::start
label start() const
Definition: polyPatch.H:357
Foam::IOobjectList
List of IOobjects with searching and retrieving facilities.
Definition: IOobjectList.H:51
polyModifyFace.H
Foam::PtrList::resize
void resize(const label newLen)
Definition: PtrList.C:96
Foam::zoneIdentifier::index
label index() const noexcept
Definition: zoneIdentifier.H:131
Foam::IOobject::name
const word & name() const noexcept
Definition: IOobjectI.H:58
Time.H
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:49
setRootCase.H
Foam::polyMesh::faces
virtual const faceList & faces() const
Definition: polyMesh.C:1087
Foam::polyTopoChange::setAction
label setAction(const topoAction &action)
Definition: polyTopoChange.C:2474
ReadFields.H
Field reading functions for post-processing utilities.
Foam::nl
constexpr char nl
Definition: Ostream.H:424
Foam::zoneIdentifier::name
const word & name() const noexcept
Definition: zoneIdentifier.H:119
fvMeshTools.H
f
labelList f(nPoints)
Foam::fvMeshMapper
Class holds all the necessary information for mapping fields associated with fvMesh.
Definition: fvMeshMapper.H:53
Foam::primitiveMesh::nInternalFaces
label nInternalFaces() const noexcept
Definition: primitiveMeshI.H:71
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: BitOps.H:58
Foam::fvMeshTools::addPatch
static label addPatch(fvMesh &mesh, const polyPatch &patch, const dictionary &patchFieldDict, const word &defaultPatchFieldType, const bool validBoundary)
Definition: fvMeshTools.C:30
Foam::processorMeshes::removeFiles
static void removeFiles(const polyMesh &mesh)
Definition: processorMeshes.C:267
Foam::HashSet::insert
bool insert(const Key &key)
Definition: HashSet.H:191
createTime.H
Foam::name
word name(const expressions::valueTypeCode typeCode)
Definition: exprTraits.C:52
Foam::polyAddFace
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
Foam::primitiveMesh::nFaces
label nFaces() const noexcept
Definition: primitiveMeshI.H:83
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:68
Foam::patchIdentifier::name
const word & name() const noexcept
Definition: patchIdentifier.H:131
Foam::ZoneMesh::clearAddressing
void clearAddressing()
Definition: ZoneMesh.C:702
Foam::dictionary::getOrDefault
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:141
Foam::argList::addOption
static void addOption(const word &optName, const string &param="", const string &usage="", bool advanced=false)
Definition: argList.C:328
args
Foam::argList args(argc, argv)
Foam::topoSet::removeFiles
static void removeFiles(const polyMesh &)
Definition: topoSet.C:628
Foam::ZoneMesh::checkParallelSync
bool checkParallelSync(const bool report=false) const
Definition: ZoneMesh.C:745
WarningInFunction
#define WarningInFunction
Definition: messageStream.H:353
Foam::polyMesh::setInstance
void setInstance(const fileName &instance, const IOobject::writeOption wOpt=IOobject::AUTO_WRITE)
Definition: polyMeshIO.C:29
Foam::polyMesh::faceNeighbour
virtual const labelList & faceNeighbour() const
Definition: polyMesh.C:1106
Foam::argList::found
bool found(const word &optName) const
Definition: argListI.H:171