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