All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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-2021 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.
29 
30  Note: if any coupled patch face is selected for baffling the opposite
31  member has to be selected for baffling as well.
32 
33  - if the patch already exists will not override it nor its fields
34  - if the patch does not exist it will be created together with 'calculated'
35  patchfields unless the field is mentioned in the patchFields section.
36  - any 0-sized patches (since faces have been moved out) will get removed
37 
38 \*---------------------------------------------------------------------------*/
39 
40 #include "argList.H"
41 #include "Time.H"
42 #include "polyTopoChange.H"
43 #include "polyModifyFace.H"
44 #include "polyAddFace.H"
45 #include "ReadFields.H"
46 #include "volFields.H"
47 #include "surfaceFields.H"
48 #include "pointFields.H"
49 #include "fvMeshMapper.H"
50 #include "faceSelection.H"
51 #include "searchableSurface.H"
52 #include "fvMeshTools.H"
53 #include "systemDict.H"
54 
55 using namespace Foam;
56 
57 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
58 
59 label addPatch
60 (
61  fvMesh& mesh,
62  const word& patchName,
63  const word& groupName,
64  const dictionary& patchDict
65 )
66 {
67  const polyBoundaryMesh& pbm = mesh.boundaryMesh();
68 
69  if (pbm.findPatchID(patchName) == -1)
70  {
71  autoPtr<polyPatch> ppPtr
72  (
74  (
75  patchName,
76  patchDict,
77  0,
78  pbm
79  )
80  );
81  polyPatch& pp = ppPtr();
82 
83  if (!groupName.empty() && !pp.inGroup(groupName))
84  {
85  pp.inGroups().append(groupName);
86  }
87 
88 
89  // Add patch, create calculated everywhere
91  (
92  mesh,
93  pp,
94  dictionary(), // do not set specialised patchFields
96  true // parallel sync'ed addition
97  );
98  }
99  else
100  {
101  Info<< "Patch '" << patchName
102  << "' already exists. Only "
103  << "moving patch faces - type will remain the same"
104  << endl;
105  }
106 
107  return pbm.findPatchID(patchName);
108 }
109 
110 
111 // Filter out the empty patches.
112 void filterPatches(fvMesh& mesh, const HashSet<word>& addedPatchNames)
113 {
114  // Remove any zero-sized ones. Assumes
115  // - processor patches are already only there if needed
116  // - all other patches are available on all processors
117  // - but coupled ones might still be needed, even if zero-size
118  // (e.g. processorCyclic)
119  // See also logic in createPatch.
120  const polyBoundaryMesh& pbm = mesh.boundaryMesh();
121 
122  labelList oldToNew(pbm.size(), -1);
123  label newPatchi = 0;
124  forAll(pbm, patchi)
125  {
126  const polyPatch& pp = pbm[patchi];
127 
128  if (!isA<processorPolyPatch>(pp))
129  {
130  if
131  (
132  isA<coupledPolyPatch>(pp)
133  || returnReduce(pp.size(), sumOp<label>())
134  || addedPatchNames.found(pp.name())
135  )
136  {
137  // Coupled (and unknown size) or uncoupled and used
138  oldToNew[patchi] = newPatchi++;
139  }
140  }
141  }
142 
143  forAll(pbm, patchi)
144  {
145  const polyPatch& pp = pbm[patchi];
146 
147  if (isA<processorPolyPatch>(pp))
148  {
149  oldToNew[patchi] = newPatchi++;
150  }
151  }
152 
153 
154  const label nKeepPatches = newPatchi;
155 
156  // Shuffle unused ones to end
157  if (nKeepPatches != pbm.size())
158  {
159  Info<< endl
160  << "Removing zero-sized patches:" << endl << incrIndent;
161 
162  forAll(oldToNew, patchi)
163  {
164  if (oldToNew[patchi] == -1)
165  {
166  Info<< indent << pbm[patchi].name()
167  << " type " << pbm[patchi].type()
168  << " at position " << patchi << endl;
169  oldToNew[patchi] = newPatchi++;
170  }
171  }
172  Info<< decrIndent;
173 
174  fvMeshTools::reorderPatches(mesh, oldToNew, nKeepPatches, true);
175  Info<< endl;
176  }
177 }
178 
179 
180 void modifyOrAddFace
181 (
182  polyTopoChange& meshMod,
183  const face& f,
184  const label facei,
185  const label own,
186  const bool flipFaceFlux,
187  const label newPatchi,
188  const label zoneID,
189  const bool zoneFlip,
190 
191  PackedBoolList& modifiedFace
192 )
193 {
194  if (!modifiedFace[facei])
195  {
196  // First usage of face. Modify.
197  meshMod.setAction
198  (
200  (
201  f, // modified face
202  facei, // label of face
203  own, // owner
204  -1, // neighbour
205  flipFaceFlux, // face flip
206  newPatchi, // patch for face
207  false, // remove from zone
208  zoneID, // zone for face
209  zoneFlip // face flip in zone
210  )
211  );
212  modifiedFace[facei] = 1;
213  }
214  else
215  {
216  // Second or more usage of face. Add.
217  meshMod.setAction
218  (
220  (
221  f, // modified face
222  own, // owner
223  -1, // neighbour
224  -1, // master point
225  -1, // master edge
226  facei, // master face
227  flipFaceFlux, // face flip
228  newPatchi, // patch for face
229  zoneID, // zone for face
230  zoneFlip // face flip in zone
231  )
232  );
233  }
234 }
235 
236 
237 // Create faces for fZone faces. Usually newMasterPatches, newSlavePatches
238 // only size one but can be more for duplicate baffle sets
239 void createFaces
240 (
241  const bool internalFacesOnly,
242  const fvMesh& mesh,
243  const faceZone& fZone,
244  const labelList& newMasterPatches,
245  const labelList& newSlavePatches,
246  polyTopoChange& meshMod,
247  PackedBoolList& modifiedFace,
248  label& nModified
249 )
250 {
251  const polyBoundaryMesh& pbm = mesh.boundaryMesh();
252 
253  forAll(newMasterPatches, i)
254  {
255  // Pass 1. Do selected side of zone
256  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
257 
258  for (label facei = 0; facei < mesh.nInternalFaces(); facei++)
259  {
260  label zoneFacei = fZone.whichFace(facei);
261 
262  if (zoneFacei != -1)
263  {
264  if (!fZone.flipMap()[zoneFacei])
265  {
266  // Use owner side of face
267  modifyOrAddFace
268  (
269  meshMod,
270  mesh.faces()[facei], // modified face
271  facei, // label of face
272  mesh.faceOwner()[facei],// owner
273  false, // face flip
274  newMasterPatches[i], // patch for face
275  fZone.index(), // zone for face
276  false, // face flip in zone
277  modifiedFace // modify or add status
278  );
279  }
280  else
281  {
282  // Use neighbour side of face.
283  // To keep faceZone pointing out of original neighbour
284  // we don't need to set faceFlip since that cell
285  // now becomes the owner
286  modifyOrAddFace
287  (
288  meshMod,
289  mesh.faces()[facei].reverseFace(), // modified face
290  facei, // label of face
291  mesh.faceNeighbour()[facei],// owner
292  true, // face flip
293  newMasterPatches[i], // patch for face
294  fZone.index(), // zone for face
295  false, // face flip in zone
296  modifiedFace // modify or add status
297  );
298  }
299 
300  nModified++;
301  }
302  }
303 
304 
305  // Pass 2. Do other side of zone
306  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
307 
308  for (label facei = 0; facei < mesh.nInternalFaces(); facei++)
309  {
310  label zoneFacei = fZone.whichFace(facei);
311 
312  if (zoneFacei != -1)
313  {
314  if (!fZone.flipMap()[zoneFacei])
315  {
316  // Use neighbour side of face
317  modifyOrAddFace
318  (
319  meshMod,
320  mesh.faces()[facei].reverseFace(), // modified face
321  facei, // label of face
322  mesh.faceNeighbour()[facei], // owner
323  true, // face flip
324  newSlavePatches[i], // patch for face
325  fZone.index(), // zone for face
326  true, // face flip in zone
327  modifiedFace // modify or add
328  );
329  }
330  else
331  {
332  // Use owner side of face
333  modifyOrAddFace
334  (
335  meshMod,
336  mesh.faces()[facei], // modified face
337  facei, // label of face
338  mesh.faceOwner()[facei],// owner
339  false, // face flip
340  newSlavePatches[i], // patch for face
341  fZone.index(), // zone for face
342  true, // face flip in zone
343  modifiedFace // modify or add status
344  );
345  }
346  }
347  }
348 
349 
350  // Modify any boundary faces
351  // ~~~~~~~~~~~~~~~~~~~~~~~~~
352 
353  // Normal boundary:
354  // - move to new patch. Might already be back-to-back baffle
355  // you want to add cyclic to. Do warn though.
356  //
357  // Processor boundary:
358  // - do not move to cyclic
359  // - add normal patches though.
360 
361  // For warning once per patch.
362  labelHashSet patchWarned;
363 
364  forAll(pbm, patchi)
365  {
366  const polyPatch& pp = pbm[patchi];
367 
368  const label newMasterPatchi = newMasterPatches[i];
369  const label newSlavePatchi = newSlavePatches[i];
370 
371  if
372  (
373  pp.coupled()
374  && (
375  pbm[newMasterPatchi].coupled()
376  || pbm[newSlavePatchi].coupled()
377  )
378  )
379  {
380  // Do not allow coupled faces to be moved to different
381  // coupled patches.
382  }
383  else if (pp.coupled() || !internalFacesOnly)
384  {
385  forAll(pp, i)
386  {
387  label facei = pp.start()+i;
388 
389  label zoneFacei = fZone.whichFace(facei);
390 
391  if (zoneFacei != -1)
392  {
393  if (patchWarned.insert(patchi))
394  {
396  << "Found boundary face (in patch "
397  << pp.name()
398  << ") in faceZone " << fZone.name()
399  << " to convert to baffle patches "
400  << pbm[newMasterPatchi].name() << "/"
401  << pbm[newSlavePatchi].name()
402  << endl
403  << " Set internalFacesOnly to true in the"
404  << " createBaffles control dictionary if you"
405  << " don't wish to convert boundary faces."
406  << endl;
407  }
408 
409  modifyOrAddFace
410  (
411  meshMod,
412  mesh.faces()[facei], // modified face
413  facei, // label of face
414  mesh.faceOwner()[facei], // owner
415  false, // face flip
416  fZone.flipMap()[zoneFacei]
417  ? newSlavePatchi
418  : newMasterPatchi, // patch for face
419  fZone.index(), // zone for face
420  fZone.flipMap()[zoneFacei], // face flip in zone
421  modifiedFace // modify or add
422  );
423 
424  nModified++;
425  }
426  }
427  }
428  }
429  }
430 }
431 
432 
433 int main(int argc, char *argv[])
434 {
436  (
437  "Makes internal faces into boundary faces.\n"
438  "Does not duplicate points."
439  );
440  #include "addDictOption.H"
441  #include "addOverwriteOption.H"
442  #include "addDictOption.H"
443  #include "addRegionOption.H"
444  #include "setRootCase.H"
445  #include "createTime.H"
447  #include "createNamedMesh.H"
448 
449  const bool overwrite = args.optionFound("overwrite");
450 
451  const word oldInstance = mesh.pointsInstance();
452 
453  Switch internalFacesOnly(false);
454 
455  Switch fields(false);
456 
457  PtrList<faceSelection> selectors;
458  {
459  const dictionary dict(systemDict("createBafflesDict", args, mesh));
460 
461  dict.lookup("internalFacesOnly") >> internalFacesOnly;
462  fields = dict.lookupOrDefault("fields", false);
463 
464  const dictionary& selectionsDict = dict.subDict("baffles");
465 
466  label n = 0;
467  forAllConstIter(dictionary, selectionsDict, iter)
468  {
469  if (iter().isDict())
470  {
471  n++;
472  }
473  }
474  selectors.setSize(n);
475  n = 0;
476  forAllConstIter(dictionary, selectionsDict, iter)
477  {
478  if (iter().isDict())
479  {
480  selectors.set
481  (
482  n++,
483  faceSelection::New(iter().keyword(), mesh, iter().dict())
484  );
485  }
486  }
487  }
488 
489 
490  if (internalFacesOnly)
491  {
492  Info<< "Not converting faces on non-coupled patches." << nl << endl;
493  }
494 
495 
496  // Read objects in time directory
498 
499  if (fields) Info<< "Reading geometric fields" << nl << endl;
500 
501  #include "readVolFields.H"
502  #include "readSurfaceFields.H"
503  #include "readPointFields.H"
504 
505 
506  // Creating (if necessary) faceZones
507  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
508 
509  forAll(selectors, selectorI)
510  {
511  const word& name = selectors[selectorI].name();
512 
513  if (mesh.faceZones().findZoneID(name) == -1)
514  {
515  mesh.faceZones().clearAddressing();
516  label sz = mesh.faceZones().size();
517 
518  labelList addr(0);
519  boolList flip(0);
520  mesh.faceZones().setSize(sz+1);
521  mesh.faceZones().set
522  (
523  sz,
524  new faceZone(name, addr, flip, sz, mesh.faceZones())
525  );
526  }
527  }
528 
529 
530  // Select faces
531  // ~~~~~~~~~~~~
532 
533  //- Per face zoneID it is in and flip status.
534  labelList faceToZoneID(mesh.nFaces(), -1);
535  boolList faceToFlip(mesh.nFaces(), false);
536  forAll(selectors, selectorI)
537  {
538  const word& name = selectors[selectorI].name();
539  label zoneID = mesh.faceZones().findZoneID(name);
540 
541  selectors[selectorI].select(zoneID, faceToZoneID, faceToFlip);
542  }
543 
544  // Add faces to faceZones
545  labelList nFaces(mesh.faceZones().size(), 0);
546  forAll(faceToZoneID, facei)
547  {
548  label zoneID = faceToZoneID[facei];
549  if (zoneID != -1)
550  {
551  nFaces[zoneID]++;
552  }
553  }
554 
555  forAll(selectors, selectorI)
556  {
557  const word& name = selectors[selectorI].name();
558  label zoneID = mesh.faceZones().findZoneID(name);
559 
560  label& n = nFaces[zoneID];
561  labelList addr(n);
562  boolList flip(n);
563  n = 0;
564  forAll(faceToZoneID, facei)
565  {
566  label zone = faceToZoneID[facei];
567  if (zone == zoneID)
568  {
569  addr[n] = facei;
570  flip[n] = faceToFlip[facei];
571  n++;
572  }
573  }
574 
575  Info<< "Created zone " << name
576  << " at index " << zoneID
577  << " with " << n << " faces" << endl;
578 
579  mesh.faceZones().set
580  (
581  zoneID,
582  new faceZone(name, addr, flip, zoneID, mesh.faceZones())
583  );
584  }
585 
586 
587 
588  // Count patches to add
589  // ~~~~~~~~~~~~~~~~~~~~
590  HashSet<word> bafflePatches;
591  {
592  forAll(selectors, selectorI)
593  {
594  const dictionary& dict = selectors[selectorI].dict();
595 
596  if (dict.found("patches"))
597  {
598  const dictionary& patchSources = dict.subDict("patches");
599  forAllConstIter(dictionary, patchSources, iter)
600  {
601  const word patchName(iter().dict()["name"]);
602  bafflePatches.insert(patchName);
603  }
604  }
605  else
606  {
607  const word masterName = selectors[selectorI].name() + "_master";
608  bafflePatches.insert(masterName);
609  const word slaveName = selectors[selectorI].name() + "_slave";
610  bafflePatches.insert(slaveName);
611  }
612  }
613  }
614 
615 
616  // Create baffles
617  // ~~~~~~~~~~~~~~
618  // Is done in multiple steps
619  // - create patches with 'calculated' patchFields
620  // - move faces into these patches
621  // - change the patchFields to the wanted type
622  // This order is done so e.g. fixedJump works:
623  // - you cannot create patchfields at the same time as patches since
624  // they do an evaluate upon construction
625  // - you want to create the patchField only after you have faces
626  // so you don't get the 'create-from-nothing' mapping problem.
627 
628 
629  // Pass 1: add patches
630  // ~~~~~~~~~~~~~~~~~~~
631 
632  // HashSet<word> addedPatches;
633  {
634  const polyBoundaryMesh& pbm = mesh.boundaryMesh();
635  forAll(selectors, selectorI)
636  {
637  const dictionary& dict = selectors[selectorI].dict();
638  const word& groupName = selectors[selectorI].name();
639 
640  if (dict.found("patches"))
641  {
642  const dictionary& patchSources = dict.subDict("patches");
643  forAllConstIter(dictionary, patchSources, iter)
644  {
645  const word patchName(iter().dict()["name"]);
646 
647  if (pbm.findPatchID(patchName) == -1)
648  {
649  dictionary patchDict = iter().dict();
650  patchDict.set("nFaces", 0);
651  patchDict.set("startFace", 0);
652 
653  // Note: do not set coupleGroup if constructed from
654  // baffles so you have freedom specifying it
655  // yourself.
656  // patchDict.set("coupleGroup", groupName);
657 
658  addPatch(mesh, patchName, groupName, patchDict);
659  }
660  else
661  {
662  Info<< "Patch '" << patchName
663  << "' already exists. Only "
664  << "moving patch faces - type will remain the same"
665  << endl;
666  }
667  }
668  }
669  else
670  {
671  const dictionary& patchSource = dict.subDict("patchPairs");
672  const word masterName = groupName + "_master";
673  const word slaveName = groupName + "_slave";
674 
675  word groupNameMaster = groupName;
676  word groupNameSlave = groupName;
677 
678 
679  dictionary patchDictMaster(patchSource);
680  patchDictMaster.set("nFaces", 0);
681  patchDictMaster.set("startFace", 0);
682  patchDictMaster.set("coupleGroup", groupName);
683 
684  dictionary patchDictSlave(patchDictMaster);
685 
686  // Note: This is added for the particular case where we want
687  // master and slave in different groupNames
688  // (ie 3D thermal baffles)
689 
690  Switch sameGroup
691  (
692  patchSource.lookupOrDefault("sameGroup", true)
693  );
694  if (!sameGroup)
695  {
696  groupNameMaster = groupName + "Group_master";
697  groupNameSlave = groupName + "Group_slave";
698  patchDictMaster.set("coupleGroup", groupNameMaster);
699  patchDictSlave.set("coupleGroup", groupNameSlave);
700  }
701 
702  addPatch(mesh, masterName, groupNameMaster, patchDictMaster);
703  addPatch(mesh, slaveName, groupNameSlave, patchDictSlave);
704  }
705  }
706  }
707 
708 
709  // Make sure patches and zoneFaces are synchronised across couples
710  mesh.boundaryMesh().checkParallelSync(true);
711  mesh.faceZones().checkParallelSync(true);
712 
713 
714 
715  // Mesh change container
716  polyTopoChange meshMod(mesh);
717 
718  const polyBoundaryMesh& pbm = mesh.boundaryMesh();
719 
720 
721  // Do the actual changes. Note:
722  // - loop in incrementing face order (not necessary if faceZone ordered).
723  // Preserves any existing ordering on patch faces.
724  // - two passes, do non-flip faces first and flip faces second. This
725  // guarantees that when e.g. creating a cyclic all faces from one
726  // side come first and faces from the other side next.
727 
728  // Whether first use of face (modify) or consecutive (add)
729  PackedBoolList modifiedFace(mesh.nFaces());
730  label nModified = 0;
731 
732  forAll(selectors, selectorI)
733  {
734  const word& name = selectors[selectorI].name();
735  label zoneID = mesh.faceZones().findZoneID(name);
736  const faceZone& fZone = mesh.faceZones()[zoneID];
737 
738  const dictionary& dict = selectors[selectorI].dict();
739 
740  DynamicList<label> newMasterPatches;
741  DynamicList<label> newSlavePatches;
742 
743  if (dict.found("patches"))
744  {
745  const dictionary& patchSources = dict.subDict("patches");
746 
747  bool master = true;
748  forAllConstIter(dictionary, patchSources, iter)
749  {
750  const word patchName(iter().dict()["name"]);
751  label patchi = pbm.findPatchID(patchName);
752  if (master)
753  {
754  newMasterPatches.append(patchi);
755  }
756  else
757  {
758  newSlavePatches.append(patchi);
759  }
760  master = !master;
761  }
762  }
763  else
764  {
765  const word masterName = selectors[selectorI].name() + "_master";
766  newMasterPatches.append(pbm.findPatchID(masterName));
767 
768  const word slaveName = selectors[selectorI].name() + "_slave";
769  newSlavePatches.append(pbm.findPatchID(slaveName));
770  }
771 
772 
773 
774  createFaces
775  (
776  internalFacesOnly,
777  mesh,
778  fZone,
779  newMasterPatches,
780  newSlavePatches,
781  meshMod,
782  modifiedFace,
783  nModified
784  );
785  }
786 
787 
788  Info<< "Converted " << returnReduce(nModified, sumOp<label>())
789  << " faces into boundary faces in patches "
790  << bafflePatches.sortedToc() << nl << endl;
791 
792  if (!overwrite)
793  {
794  runTime++;
795  }
796 
797  // Change the mesh. Change points directly (no inflation).
798  autoPtr<mapPolyMesh> map = meshMod.changeMesh(mesh, false);
799 
800  // Update fields
801  mesh.updateMesh(map);
802 
803 
804 
805  // Correct boundary faces mapped-out-of-nothing.
806  // This is just a hack to correct the value field.
807  {
808  fvMeshMapper mapper(mesh, map);
809  bool hasWarned = false;
810 
811  forAllConstIter(HashSet<word>, bafflePatches, iter)
812  {
813  label patchi = mesh.boundaryMesh().findPatchID(iter.key());
814 
815  const fvPatchMapper& pm = mapper.boundaryMap()[patchi];
816 
817  if (pm.sizeBeforeMapping() == 0)
818  {
819  if (!hasWarned)
820  {
821  hasWarned = true;
823  << "Setting field on boundary faces to zero." << endl
824  << "You might have to edit these fields." << endl;
825  }
826 
827  fvMeshTools::zeroPatchFields(mesh, patchi);
828  }
829  }
830  }
831 
832 
833  // Pass 2: change patchFields
834  // ~~~~~~~~~~~~~~~~~~~~~~~~~~
835 
836  {
837  const polyBoundaryMesh& pbm = mesh.boundaryMesh();
838  forAll(selectors, selectorI)
839  {
840  const dictionary& dict = selectors[selectorI].dict();
841  if (dict.found("patches"))
842  {
843  const dictionary& patchSources = dict.subDict("patches");
844 
845  forAllConstIter(dictionary, patchSources, iter)
846  {
847  const word patchName(iter().dict()["name"]);
848  label patchi = pbm.findPatchID(patchName);
849 
850  if (iter().dict().found("patchFields"))
851  {
852  const dictionary& patchFieldsDict =
853  iter().dict().subDict
854  (
855  "patchFields"
856  );
857 
858  fvMeshTools::setPatchFields
859  (
860  mesh,
861  patchi,
862  patchFieldsDict
863  );
864  }
865  }
866  }
867  else
868  {
869  const dictionary& patchSource = dict.subDict("patchPairs");
870 
871  Switch sameGroup
872  (
873  patchSource.lookupOrDefault("sameGroup", true)
874  );
875 
876  const word& groupName = selectors[selectorI].name();
877 
878  if (patchSource.found("patchFields"))
879  {
880  dictionary patchFieldsDict = patchSource.subDict
881  (
882  "patchFields"
883  );
884 
885  if (sameGroup)
886  {
887  // Add coupleGroup to all entries
888  forAllIter(dictionary, patchFieldsDict, iter)
889  {
890  if (iter().isDict())
891  {
892  dictionary& dict = iter().dict();
893  dict.set("coupleGroup", groupName);
894  }
895  }
896 
897  const labelList& patchIDs =
898  pbm.groupPatchIDs()[groupName];
899 
900  forAll(patchIDs, i)
901  {
902  fvMeshTools::setPatchFields
903  (
904  mesh,
905  patchIDs[i],
906  patchFieldsDict
907  );
908  }
909  }
910  else
911  {
912  const word masterPatchName(groupName + "_master");
913  const word slavePatchName(groupName + "_slave");
914 
915  label patchiMaster = pbm.findPatchID(masterPatchName);
916  label patchiSlave = pbm.findPatchID(slavePatchName);
917 
918  fvMeshTools::setPatchFields
919  (
920  mesh,
921  patchiMaster,
922  patchFieldsDict
923  );
924 
925  fvMeshTools::setPatchFields
926  (
927  mesh,
928  patchiSlave,
929  patchFieldsDict
930  );
931  }
932  }
933  }
934  }
935  }
936 
937 
938  // Move mesh (since morphing might not do this)
939  if (map().hasMotionPoints())
940  {
941  mesh.movePoints(map().preMotionPoints());
942  }
943 
944 
945  // Remove any now zero-sized patches
946  filterPatches(mesh, bafflePatches);
947 
948 
949  if (overwrite)
950  {
951  mesh.setInstance(oldInstance);
952  }
953 
954  Info<< "Writing mesh to " << runTime.timeName() << endl;
955 
956  mesh.write();
957 
958  Info<< "End\n" << endl;
959 
960  return 0;
961 }
962 
963 
964 // ************************************************************************* //
Foam::surfaceFields.
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:434
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:140
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:643
label findZoneID(const word &zoneName) const
Find zone index given a name.
Definition: MeshZones.C:341
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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:303
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:221
bool set(const label) const
Is element set.
Definition: PtrListI.H:65
List of IOobjects with searching and retrieving facilities.
Definition: IOobjectList.H:50
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:156
label nInternalFaces() const
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1175
label nFaces() const
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:459
bool inGroup(const word &) const
Test if in group.
const boolList & flipMap() const
Return face flip map.
Definition: faceZone.H:248
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
bool checkParallelSync(const bool report=false) const
Check whether all procs have all zones and in same order. Return.
Definition: MeshZones.C:428
bool optionFound(const word &opt) const
Return true if the named option is found.
Definition: argListI.H:114
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/any.
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:636
const wordList & inGroups() const
Return the optional groups patch belongs to.
Info<< "Calculating turbulent flame speed field St\"<< endl;volScalarField St(IOobject("St", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), flameWrinkling->Xi() *Su);multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:228
virtual bool write(const bool write=true) const
Write mesh using IO settings from time.
Definition: fvMesh.C:1093
bool insert(const Key &key)
Insert a new entry.
Definition: HashSet.H:111
IOdictionary systemDict(const word &dictName, const argList &args, const objectRegistry &ob, const word &regionName=polyMesh::defaultRegion)
Definition: systemDict.C:92
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
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:982
static void zeroPatchFields(fvMesh &mesh, const label patchi)
Change patchField to zero on registered fields.
Definition: fvMeshTools.C:113
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:802
virtual void updateMesh(const mapPolyMesh &mpm)
Update mesh corresponding to the given map.
Definition: fvMesh.C:795
virtual bool coupled() const
Return true if this patch is geometrically coupled (i.e. faces and.
Definition: polyPatch.H:319
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:33
void append(const T &)
Append an element at the end of the list.
Definition: ListI.H:178
This boundary condition is not designed to be evaluated; it is assumed that the value is assigned via...
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1169
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
Definition: DynamicListI.H:296
const word & name() const
Return name.
Definition: zone.H:147
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1156
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:260
objects
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition: Ostream.H:235
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:716
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:70
const meshFaceZones & faceZones() const
Return face zones.
Definition: polyMesh.H:476
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:309
void set(entry *)
Assign a new entry, overwrite any existing entry.
Definition: dictionary.C:1271
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:159
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:228
void clearAddressing()
Clear addressing.
Definition: MeshZones.C:387
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:844
const HashTable< labelList, word > & groupPatchIDs() const
Per patch group the patch indices.
label setAction(const topoAction &action)
For compatibility with polyTopoChange: set topological action.