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 | Copyright (C) 2011-2015 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 
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 "fvMeshMapper.H"
49 #include "faceSelection.H"
50 
51 #include "fvMeshTools.H"
52 
53 using namespace Foam;
54 
55 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
56 
57 label addPatch
58 (
59  fvMesh& mesh,
60  const word& patchName,
61  const word& groupName,
62  const dictionary& patchDict
63 )
64 {
65  const polyBoundaryMesh& pbm = mesh.boundaryMesh();
66 
67  if (pbm.findPatchID(patchName) == -1)
68  {
69  autoPtr<polyPatch> ppPtr
70  (
72  (
73  patchName,
74  patchDict,
75  0,
76  pbm
77  )
78  );
79  polyPatch& pp = ppPtr();
80 
81  if (!groupName.empty() && !pp.inGroup(groupName))
82  {
83  pp.inGroups().append(groupName);
84  }
85 
86 
87  // Add patch, create calculated everywhere
89  (
90  mesh,
91  pp,
92  dictionary(), // do not set specialised patchFields
94  true // parallel sync'ed addition
95  );
96  }
97  else
98  {
99  Info<< "Patch '" << patchName
100  << "' already exists. Only "
101  << "moving patch faces - type will remain the same"
102  << endl;
103  }
104 
105  return pbm.findPatchID(patchName);
106 }
107 
108 
109 void modifyOrAddFace
110 (
111  polyTopoChange& meshMod,
112  const face& f,
113  const label faceI,
114  const label own,
115  const bool flipFaceFlux,
116  const label newPatchI,
117  const label zoneID,
118  const bool zoneFlip,
119 
120  PackedBoolList& modifiedFace
121 )
122 {
123  if (!modifiedFace[faceI])
124  {
125  // First usage of face. Modify.
126  meshMod.setAction
127  (
129  (
130  f, // modified face
131  faceI, // label of face
132  own, // owner
133  -1, // neighbour
134  flipFaceFlux, // face flip
135  newPatchI, // patch for face
136  false, // remove from zone
137  zoneID, // zone for face
138  zoneFlip // face flip in zone
139  )
140  );
141  modifiedFace[faceI] = 1;
142  }
143  else
144  {
145  // Second or more usage of face. Add.
146  meshMod.setAction
147  (
149  (
150  f, // modified face
151  own, // owner
152  -1, // neighbour
153  -1, // master point
154  -1, // master edge
155  faceI, // master face
156  flipFaceFlux, // face flip
157  newPatchI, // patch for face
158  zoneID, // zone for face
159  zoneFlip // face flip in zone
160  )
161  );
162  }
163 }
164 
165 
166 // Create faces for fZone faces. Usually newMasterPatches, newSlavePatches
167 // only size one but can be more for duplicate baffle sets
168 void createFaces
169 (
170  const bool internalFacesOnly,
171  const fvMesh& mesh,
172  const faceZone& fZone,
173  const labelList& newMasterPatches,
174  const labelList& newSlavePatches,
175  polyTopoChange& meshMod,
176  PackedBoolList& modifiedFace,
177  label& nModified
178 )
179 {
180  const polyBoundaryMesh& pbm = mesh.boundaryMesh();
181 
182  forAll(newMasterPatches, i)
183  {
184  // Pass 1. Do selected side of zone
185  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
186 
187  for (label faceI = 0; faceI < mesh.nInternalFaces(); faceI++)
188  {
189  label zoneFaceI = fZone.whichFace(faceI);
190 
191  if (zoneFaceI != -1)
192  {
193  if (!fZone.flipMap()[zoneFaceI])
194  {
195  // Use owner side of face
196  modifyOrAddFace
197  (
198  meshMod,
199  mesh.faces()[faceI], // modified face
200  faceI, // label of face
201  mesh.faceOwner()[faceI],// owner
202  false, // face flip
203  newMasterPatches[i], // patch for face
204  fZone.index(), // zone for face
205  false, // face flip in zone
206  modifiedFace // modify or add status
207  );
208  }
209  else
210  {
211  // Use neighbour side of face.
212  // To keep faceZone pointing out of original neighbour
213  // we don't need to set faceFlip since that cell
214  // now becomes the owner
215  modifyOrAddFace
216  (
217  meshMod,
218  mesh.faces()[faceI].reverseFace(), // modified face
219  faceI, // label of face
220  mesh.faceNeighbour()[faceI],// owner
221  true, // face flip
222  newMasterPatches[i], // patch for face
223  fZone.index(), // zone for face
224  false, // face flip in zone
225  modifiedFace // modify or add status
226  );
227  }
228 
229  nModified++;
230  }
231  }
232 
233 
234  // Pass 2. Do other side of zone
235  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
236 
237  for (label faceI = 0; faceI < mesh.nInternalFaces(); faceI++)
238  {
239  label zoneFaceI = fZone.whichFace(faceI);
240 
241  if (zoneFaceI != -1)
242  {
243  if (!fZone.flipMap()[zoneFaceI])
244  {
245  // Use neighbour side of face
246  modifyOrAddFace
247  (
248  meshMod,
249  mesh.faces()[faceI].reverseFace(), // modified face
250  faceI, // label of face
251  mesh.faceNeighbour()[faceI], // owner
252  true, // face flip
253  newSlavePatches[i], // patch for face
254  fZone.index(), // zone for face
255  true, // face flip in zone
256  modifiedFace // modify or add
257  );
258  }
259  else
260  {
261  // Use owner side of face
262  modifyOrAddFace
263  (
264  meshMod,
265  mesh.faces()[faceI], // modified face
266  faceI, // label of face
267  mesh.faceOwner()[faceI],// owner
268  false, // face flip
269  newSlavePatches[i], // patch for face
270  fZone.index(), // zone for face
271  true, // face flip in zone
272  modifiedFace // modify or add status
273  );
274  }
275  }
276  }
277 
278 
279  // Modify any boundary faces
280  // ~~~~~~~~~~~~~~~~~~~~~~~~~
281 
282  // Normal boundary:
283  // - move to new patch. Might already be back-to-back baffle
284  // you want to add cyclic to. Do warn though.
285  //
286  // Processor boundary:
287  // - do not move to cyclic
288  // - add normal patches though.
289 
290  // For warning once per patch.
291  labelHashSet patchWarned;
292 
293  forAll(pbm, patchI)
294  {
295  const polyPatch& pp = pbm[patchI];
296 
297  label newPatchI = newMasterPatches[i];
298 
299  if (pp.coupled() && pbm[newPatchI].coupled())
300  {
301  // Do not allow coupled faces to be moved to different
302  // coupled patches.
303  }
304  else if (pp.coupled() || !internalFacesOnly)
305  {
306  forAll(pp, i)
307  {
308  label faceI = pp.start()+i;
309 
310  label zoneFaceI = fZone.whichFace(faceI);
311 
312  if (zoneFaceI != -1)
313  {
314  if (patchWarned.insert(patchI))
315  {
317  << "Found boundary face (in patch "
318  << pp.name()
319  << ") in faceZone " << fZone.name()
320  << " to convert to baffle patch "
321  << pbm[newPatchI].name()
322  << endl
323  << " Run with -internalFacesOnly option"
324  << " if you don't wish to convert"
325  << " boundary faces." << endl;
326  }
327 
328  modifyOrAddFace
329  (
330  meshMod,
331  mesh.faces()[faceI], // modified face
332  faceI, // label of face
333  mesh.faceOwner()[faceI], // owner
334  false, // face flip
335  newPatchI, // patch for face
336  fZone.index(), // zone for face
337  fZone.flipMap()[zoneFaceI], // face flip in zone
338  modifiedFace // modify or add
339  );
340  nModified++;
341  }
342  }
343  }
344  }
345  }
346 }
347 
348 
349 int main(int argc, char *argv[])
350 {
352  (
353  "Makes internal faces into boundary faces.\n"
354  "Does not duplicate points."
355  );
356  #include "addDictOption.H"
357  #include "addOverwriteOption.H"
358  #include "addDictOption.H"
359  #include "addRegionOption.H"
360  #include "setRootCase.H"
361  #include "createTime.H"
362  runTime.functionObjects().off();
363  #include "createNamedMesh.H"
364 
365 
366  const bool overwrite = args.optionFound("overwrite");
367 
368  const word oldInstance = mesh.pointsInstance();
369 
370  const word dictName("createBafflesDict");
371  #include "setSystemMeshDictionaryIO.H"
372 
373  Switch internalFacesOnly(false);
374 
375  Switch noFields(false);
376 
377  PtrList<faceSelection> selectors;
378  {
379  Info<< "Reading baffle criteria from " << dictName << nl << endl;
381 
382  dict.lookup("internalFacesOnly") >> internalFacesOnly;
383  noFields = dict.lookupOrDefault("noFields", false);
384 
385  const dictionary& selectionsDict = dict.subDict("baffles");
386 
387  label n = 0;
388  forAllConstIter(dictionary, selectionsDict, iter)
389  {
390  if (iter().isDict())
391  {
392  n++;
393  }
394  }
395  selectors.setSize(n);
396  n = 0;
397  forAllConstIter(dictionary, selectionsDict, iter)
398  {
399  if (iter().isDict())
400  {
401  selectors.set
402  (
403  n++,
404  faceSelection::New(iter().keyword(), mesh, iter().dict())
405  );
406  }
407  }
408  }
409 
410 
411  if (internalFacesOnly)
412  {
413  Info<< "Not converting faces on non-coupled patches." << nl << endl;
414  }
415 
416 
417  // Read objects in time directory
418  IOobjectList objects(mesh, runTime.timeName());
419 
420  // Read vol fields.
421  Info<< "Reading geometric fields" << nl << endl;
423  if (!noFields) ReadFields(mesh, objects, vsFlds);
424 
426  if (!noFields) ReadFields(mesh, objects, vvFlds);
427 
429  if (!noFields) ReadFields(mesh, objects, vstFlds);
430 
431  PtrList<volSymmTensorField> vsymtFlds;
432  if (!noFields) ReadFields(mesh, objects, vsymtFlds);
433 
435  if (!noFields) ReadFields(mesh, objects, vtFlds);
436 
437  // Read surface fields.
438 
440  if (!noFields) ReadFields(mesh, objects, ssFlds);
441 
443  if (!noFields) ReadFields(mesh, objects, svFlds);
444 
446  if (!noFields) ReadFields(mesh, objects, sstFlds);
447 
449  if (!noFields) ReadFields(mesh, objects, ssymtFlds);
450 
452  if (!noFields) ReadFields(mesh, objects, stFlds);
453 
454 
455 
456 
457  // Creating (if necessary) faceZones
458  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
459 
460  forAll(selectors, selectorI)
461  {
462  const word& name = selectors[selectorI].name();
463 
464  if (mesh.faceZones().findZoneID(name) == -1)
465  {
467  label sz = mesh.faceZones().size();
468 
469  labelList addr(0);
470  boolList flip(0);
471  mesh.faceZones().setSize(sz+1);
472  mesh.faceZones().set
473  (
474  sz,
475  new faceZone(name, addr, flip, sz, mesh.faceZones())
476  );
477  }
478  }
479 
480 
481  // Select faces
482  // ~~~~~~~~~~~~
483 
484  //- Per face zoneID it is in and flip status.
485  labelList faceToZoneID(mesh.nFaces(), -1);
486  boolList faceToFlip(mesh.nFaces(), false);
487  forAll(selectors, selectorI)
488  {
489  const word& name = selectors[selectorI].name();
490  label zoneID = mesh.faceZones().findZoneID(name);
491 
492  selectors[selectorI].select(zoneID, faceToZoneID, faceToFlip);
493  }
494 
495  // Add faces to faceZones
496  labelList nFaces(mesh.faceZones().size(), 0);
497  forAll(faceToZoneID, faceI)
498  {
499  label zoneID = faceToZoneID[faceI];
500  if (zoneID != -1)
501  {
502  nFaces[zoneID]++;
503  }
504  }
505 
506  forAll(selectors, selectorI)
507  {
508  const word& name = selectors[selectorI].name();
509  label zoneID = mesh.faceZones().findZoneID(name);
510 
511  label& n = nFaces[zoneID];
512  labelList addr(n);
513  boolList flip(n);
514  n = 0;
515  forAll(faceToZoneID, faceI)
516  {
517  label zone = faceToZoneID[faceI];
518  if (zone == zoneID)
519  {
520  addr[n] = faceI;
521  flip[n] = faceToFlip[faceI];
522  n++;
523  }
524  }
525 
526  Info<< "Created zone " << name
527  << " at index " << zoneID
528  << " with " << returnReduce(n, sumOp<label>()) << " faces" << endl;
529 
530  mesh.faceZones().set
531  (
532  zoneID,
533  new faceZone(name, addr, flip, zoneID, mesh.faceZones())
534  );
535  }
536 
537 
538 
539  // Count patches to add
540  // ~~~~~~~~~~~~~~~~~~~~
541  HashSet<word> bafflePatches;
542  {
543  forAll(selectors, selectorI)
544  {
545  const dictionary& dict = selectors[selectorI].dict();
546 
547  if (dict.found("patches"))
548  {
549  const dictionary& patchSources = dict.subDict("patches");
550  forAllConstIter(dictionary, patchSources, iter)
551  {
552  const word patchName(iter().dict()["name"]);
553  bafflePatches.insert(patchName);
554  }
555  }
556  else
557  {
558  const word masterName = selectors[selectorI].name() + "_master";
559  bafflePatches.insert(masterName);
560  const word slaveName = selectors[selectorI].name() + "_slave";
561  bafflePatches.insert(slaveName);
562  }
563  }
564  }
565 
566 
567  // Create baffles
568  // ~~~~~~~~~~~~~~
569  // Is done in multiple steps
570  // - create patches with 'calculated' patchFields
571  // - move faces into these patches
572  // - change the patchFields to the wanted type
573  // This order is done so e.g. fixedJump works:
574  // - you cannot create patchfields at the same time as patches since
575  // they do an evaluate upon construction
576  // - you want to create the patchField only after you have faces
577  // so you don't get the 'create-from-nothing' mapping problem.
578 
579 
580  // Pass 1: add patches
581  // ~~~~~~~~~~~~~~~~~~~
582 
583  //HashSet<word> addedPatches;
584  {
585  const polyBoundaryMesh& pbm = mesh.boundaryMesh();
586  forAll(selectors, selectorI)
587  {
588  const dictionary& dict = selectors[selectorI].dict();
589  const word& groupName = selectors[selectorI].name();
590 
591  if (dict.found("patches"))
592  {
593  const dictionary& patchSources = dict.subDict("patches");
594  forAllConstIter(dictionary, patchSources, iter)
595  {
596  const word patchName(iter().dict()["name"]);
597 
598  if (pbm.findPatchID(patchName) == -1)
599  {
600  dictionary patchDict = iter().dict();
601  patchDict.set("nFaces", 0);
602  patchDict.set("startFace", 0);
603 
604  // Note: do not set coupleGroup if constructed from
605  // baffles so you have freedom specifying it
606  // yourself.
607  //patchDict.set("coupleGroup", groupName);
608 
609  addPatch(mesh, patchName, groupName, patchDict);
610  }
611  else
612  {
613  Info<< "Patch '" << patchName
614  << "' already exists. Only "
615  << "moving patch faces - type will remain the same"
616  << endl;
617  }
618  }
619  }
620  else
621  {
622  const dictionary& patchSource = dict.subDict("patchPairs");
623  const word masterName = groupName + "_master";
624  const word slaveName = groupName + "_slave";
625 
626  word groupNameMaster = groupName;
627  word groupNameSlave = groupName;
628 
629 
630  dictionary patchDictMaster(patchSource);
631  patchDictMaster.set("nFaces", 0);
632  patchDictMaster.set("startFace", 0);
633  patchDictMaster.set("coupleGroup", groupName);
634 
635  dictionary patchDictSlave(patchDictMaster);
636 
637  // Note: This is added for the particular case where we want
638  // master and slave in different groupNames
639  // (ie 3D thermal baffles)
640 
641  Switch sameGroup
642  (
643  patchSource.lookupOrDefault("sameGroup", true)
644  );
645  if (!sameGroup)
646  {
647  groupNameMaster = groupName + "Group_master";
648  groupNameSlave = groupName + "Group_slave";
649  patchDictMaster.set("coupleGroup", groupNameMaster);
650  patchDictSlave.set("coupleGroup", groupNameSlave);
651  }
652 
653  addPatch(mesh, masterName, groupNameMaster, patchDictMaster);
654  addPatch(mesh, slaveName, groupNameSlave, patchDictSlave);
655  }
656  }
657  }
658 
659 
660  // Make sure patches and zoneFaces are synchronised across couples
663 
664 
665 
666  // Mesh change container
667  polyTopoChange meshMod(mesh);
668 
669  const polyBoundaryMesh& pbm = mesh.boundaryMesh();
670 
671 
672  // Do the actual changes. Note:
673  // - loop in incrementing face order (not necessary if faceZone ordered).
674  // Preserves any existing ordering on patch faces.
675  // - two passes, do non-flip faces first and flip faces second. This
676  // guarantees that when e.g. creating a cyclic all faces from one
677  // side come first and faces from the other side next.
678 
679  // Whether first use of face (modify) or consecutive (add)
680  PackedBoolList modifiedFace(mesh.nFaces());
681  label nModified = 0;
682 
683  forAll(selectors, selectorI)
684  {
685  const word& name = selectors[selectorI].name();
686  label zoneID = mesh.faceZones().findZoneID(name);
687  const faceZone& fZone = mesh.faceZones()[zoneID];
688 
689  const dictionary& dict = selectors[selectorI].dict();
690 
691  DynamicList<label> newMasterPatches;
692  DynamicList<label> newSlavePatches;
693 
694  if (dict.found("patches"))
695  {
696  const dictionary& patchSources = dict.subDict("patches");
697 
698  bool master = true;
699  forAllConstIter(dictionary, patchSources, iter)
700  {
701  const word patchName(iter().dict()["name"]);
702  label patchI = pbm.findPatchID(patchName);
703  if (master)
704  {
705  newMasterPatches.append(patchI);
706  }
707  else
708  {
709  newSlavePatches.append(patchI);
710  }
711  master = !master;
712  }
713  }
714  else
715  {
716  const word masterName = selectors[selectorI].name() + "_master";
717  newMasterPatches.append(pbm.findPatchID(masterName));
718 
719  const word slaveName = selectors[selectorI].name() + "_slave";
720  newSlavePatches.append(pbm.findPatchID(slaveName));
721  }
722 
723 
724 
725  createFaces
726  (
727  internalFacesOnly,
728  mesh,
729  fZone,
730  newMasterPatches,
731  newSlavePatches,
732  meshMod,
733  modifiedFace,
734  nModified
735  );
736  }
737 
738 
739  Info<< "Converted " << returnReduce(nModified, sumOp<label>())
740  << " faces into boundary faces in patches "
741  << bafflePatches.sortedToc() << nl << endl;
742 
743  if (!overwrite)
744  {
745  runTime++;
746  }
747 
748  // Change the mesh. Change points directly (no inflation).
749  autoPtr<mapPolyMesh> map = meshMod.changeMesh(mesh, false);
750 
751  // Update fields
752  mesh.updateMesh(map);
753 
754 
755 
756  // Correct boundary faces mapped-out-of-nothing.
757  // This is just a hack to correct the value field.
758  {
759  fvMeshMapper mapper(mesh, map);
760  bool hasWarned = false;
761 
762  forAllConstIter(HashSet<word>, bafflePatches, iter)
763  {
764  label patchI = mesh.boundaryMesh().findPatchID(iter.key());
765 
766  const fvPatchMapper& pm = mapper.boundaryMap()[patchI];
767 
768  if (pm.sizeBeforeMapping() == 0)
769  {
770  if (!hasWarned)
771  {
772  hasWarned = true;
774  << "Setting field on boundary faces to zero." << endl
775  << "You might have to edit these fields." << endl;
776  }
777 
779  }
780  }
781  }
782 
783 
784  // Pass 2: change patchFields
785  // ~~~~~~~~~~~~~~~~~~~~~~~~~~
786 
787  {
788  const polyBoundaryMesh& pbm = mesh.boundaryMesh();
789  forAll(selectors, selectorI)
790  {
791  const dictionary& dict = selectors[selectorI].dict();
792  if (dict.found("patches"))
793  {
794  const dictionary& patchSources = dict.subDict("patches");
795 
796  forAllConstIter(dictionary, patchSources, iter)
797  {
798  const word patchName(iter().dict()["name"]);
799  label patchI = pbm.findPatchID(patchName);
800 
801  if (iter().dict().found("patchFields"))
802  {
803  const dictionary& patchFieldsDict =
804  iter().dict().subDict
805  (
806  "patchFields"
807  );
808 
810  (
811  mesh,
812  patchI,
813  patchFieldsDict
814  );
815  }
816  }
817  }
818  else
819  {
820  const dictionary& patchSource = dict.subDict("patchPairs");
821 
822  Switch sameGroup
823  (
824  patchSource.lookupOrDefault("sameGroup", true)
825  );
826 
827  const word& groupName = selectors[selectorI].name();
828 
829  if (patchSource.found("patchFields"))
830  {
831  dictionary patchFieldsDict = patchSource.subDict
832  (
833  "patchFields"
834  );
835 
836  if (sameGroup)
837  {
838  // Add coupleGroup to all entries
839  forAllIter(dictionary, patchFieldsDict, iter)
840  {
841  if (iter().isDict())
842  {
843  dictionary& dict = iter().dict();
844  dict.set("coupleGroup", groupName);
845  }
846  }
847 
848  const labelList& patchIDs =
849  pbm.groupPatchIDs()[groupName];
850 
851  forAll(patchIDs, i)
852  {
854  (
855  mesh,
856  patchIDs[i],
857  patchFieldsDict
858  );
859  }
860  }
861  else
862  {
863  const word masterPatchName(groupName + "_master");
864  const word slavePatchName(groupName + "_slave");
865 
866  label patchIMaster = pbm.findPatchID(masterPatchName);
867  label patchISlave = pbm.findPatchID(slavePatchName);
868 
870  (
871  mesh,
872  patchIMaster,
873  patchFieldsDict
874  );
875 
877  (
878  mesh,
879  patchISlave,
880  patchFieldsDict
881  );
882  }
883  }
884  }
885  }
886  }
887 
888 
889  // Move mesh (since morphing might not do this)
890  if (map().hasMotionPoints())
891  {
892  mesh.movePoints(map().preMotionPoints());
893  }
894 
895  if (overwrite)
896  {
897  mesh.setInstance(oldInstance);
898  }
899 
900  Info<< "Writing mesh to " << runTime.timeName() << endl;
901 
902  mesh.write();
903 
904  Info<< "End\n" << endl;
905 
906  return 0;
907 }
908 
909 
910 // ************************************************************************* //
faceSelection.H
Foam::fvPatchMapper
Mapping class for a fvPatchField.
Definition: fvPatchMapper.H:55
volFields.H
Foam::IOdictionary
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:53
Foam::polyBoundaryMesh::groupPatchIDs
const HashTable< labelList, word > & groupPatchIDs() const
Per patch group the patch indices.
Definition: polyBoundaryMesh.C:429
Foam::PackedBoolList
A bit-packed bool list.
Definition: PackedBoolList.H:63
Foam::faceSelection::New
static autoPtr< faceSelection > New(const word &name, const fvMesh &mesh, const dictionary &dict)
Return a reference to the selected faceSelection.
Foam::Switch
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition: Switch.H:60
Foam::polyBoundaryMesh::findPatchID
label findPatchID(const word &patchName) const
Find patch index given a name.
Definition: polyBoundaryMesh.C:678
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
forAllIter
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:431
Foam::returnReduce
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Definition: PstreamReduceOps.H:86
Foam::polyBoundaryMesh
Foam::polyBoundaryMesh.
Definition: polyBoundaryMesh.H:60
Foam::argList::addNote
static void addNote(const string &)
Add extra notes for the usage information.
Definition: argList.C:139
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::polyTopoChange::changeMesh
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.
Definition: polyTopoChange.C:3039
Foam::DynamicList< label >
polyAddFace.H
addOverwriteOption.H
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)
Helper routine to read Geometric fields.
Foam::polyPatch::coupled
virtual bool coupled() const
Return true if this patch is geometrically coupled (i.e. faces and.
Definition: polyPatch.H:322
polyTopoChange.H
Foam::dictionary::lookup
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:449
Foam::polyTopoChange
Direct mesh changes based on v1.3 polyTopoChange syntax.
Definition: polyTopoChange.H:97
Foam::zone
Base class for zones.
Definition: zone.H:57
Foam::dictionary::lookupOrDefault
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
Definition: dictionaryTemplates.C:33
Foam::fvMeshTools::zeroPatchFields
static void zeroPatchFields(fvMesh &mesh, const label patchI)
Change patchField to zero on registered fields.
Definition: fvMeshTools.C:241
Foam::polyModifyFace
Class describing modification of a face.
Definition: polyModifyFace.H:48
Foam::fvMesh::movePoints
virtual tmp< scalarField > movePoints(const pointField &)
Move points, returns volumes swept by faces in motion.
Definition: fvMesh.C:726
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:421
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
surfaceFields.H
Foam::surfaceFields.
Foam::fvPatchMapper::sizeBeforeMapping
virtual label sizeBeforeMapping() const
Return size of field before mapping.
Definition: fvPatchMapper.H:131
Foam::HashSet< label, Hash< label > >
setSystemMeshDictionaryIO.H
forAllConstIter
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:39
Foam::polyMesh::faceZones
const faceZoneMesh & faceZones() const
Return face zone mesh.
Definition: polyMesh.H:463
addDictOption.H
Foam::fvMesh::write
virtual bool write() const
Write mesh using IO settings from time.
Definition: fvMesh.C:873
dictName
const word dictName("particleTrackDict")
Foam::polyMesh::pointsInstance
const fileName & pointsInstance() const
Return the current instance directory for points.
Definition: polyMesh.C:772
n
label n
Definition: TABSMDCalcMethod2.H:31
Foam::polyBoundaryMesh::checkParallelSync
bool checkParallelSync(const bool report=false) const
Check whether all procs have all patches and in same order. Return.
Definition: polyBoundaryMesh.C:859
Foam::PtrList::set
bool set(const label) const
Is element set.
fvMeshMapper.H
Foam::label
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
Foam::patchIdentifier::inGroup
bool inGroup(const word &) const
Test if in group.
Definition: patchIdentifier.C:83
Foam::dictionary::found
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:304
Foam::patchIdentifier::inGroups
const wordList & inGroups() const
Return the optional groups patch belongs to.
Definition: patchIdentifier.H:145
Foam::faceZone
A subset of mesh faces organised as a primitive patch.
Definition: faceZone.H:64
Foam::List::append
void append(const T &)
Append an element at the end of the list.
Foam::nl
static const char nl
Definition: Ostream.H:260
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:66
argList.H
Foam::polyMesh::faceOwner
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1017
addRegionOption.H
Foam::fvMesh::updateMesh
virtual void updateMesh(const mapPolyMesh &mpm)
Update mesh corresponding to the given map.
Definition: fvMesh.C:802
Foam::IOobject::name
const word & name() const
Return name.
Definition: IOobject.H:273
main
int main(int argc, char *argv[])
Definition: postCalc.C:54
Foam::PtrList
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:61
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
Foam::primitiveMesh::nInternalFaces
label nInternalFaces() const
Definition: primitiveMeshI.H:52
createNamedMesh.H
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)
Return a pointer to a new patch created on freestore from.
Definition: polyPatchNew.C:32
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Foam::calculatedFvPatchField
This boundary condition is not designed to be evaluated; it is assmued that the value is assigned via...
Definition: calculatedFvPatchField.H:65
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::faceZone::whichFace
label whichFace(const label globalCellID) const
Helper function to re-direct to zone::localID(...)
Definition: faceZone.C:314
Foam::zone::name
const word & name() const
Return name.
Definition: zone.H:150
Foam::ZoneMesh::findZoneID
label findZoneID(const word &zoneName) const
Find zone index given a name.
Definition: ZoneMesh.C:348
Foam::HashTable< nil, word, string::hash >::sortedToc
List< Key > sortedToc() const
Return the table of contents as a sorted list.
Definition: HashTable.C:216
Foam::polyPatch::start
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:312
Foam::IOobjectList
List of IOobjects with searching and retrieving facilities.
Definition: IOobjectList.H:50
polyModifyFace.H
Foam::DynamicList::append
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
Foam::polyMesh::setInstance
void setInstance(const fileName &)
Set the instance for mesh files.
Definition: polyMeshIO.C:32
found
bool found
Definition: TABSMDCalcMethod2.H:32
Foam::autoPtr
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:117
setRootCase.H
Foam::polyMesh::faces
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1004
Foam::polyTopoChange::setAction
label setAction(const topoAction &action)
For compatibility with polyTopoChange: set topological action.
Definition: polyTopoChange.C:2530
Foam::primitiveMesh::nFaces
label nFaces() const
Definition: primitiveMeshI.H:58
ReadFields.H
Helper routine to read fields.
Foam::sumOp
Definition: ops.H:162
fvMeshTools.H
f
labelList f(nPoints)
Foam::fvMeshMapper
Class holds all the necessary information for mapping fields associated with fvMesh.
Definition: fvMeshMapper.H:55
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
Foam::fvMeshTools::addPatch
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:35
Foam::fvMeshTools::setPatchFields
static void setPatchFields(fvMesh &mesh, const label patchI, const dictionary &patchFieldDict)
Set patchFields according to dictionary.
Definition: fvMeshToolsTemplates.C:89
Foam::HashSet::insert
bool insert(const Key &key)
Insert a new entry.
Definition: HashSet.H:116
Foam::PtrList::setSize
void setSize(const label)
Reset size of PtrList. If extending the PtrList, new entries are.
Definition: PtrList.C:142
createTime.H
Foam::argList::optionFound
bool optionFound(const word &opt) const
Return true if the named option is found.
Definition: argListI.H:108
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:51
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
Foam::ZoneMesh::clearAddressing
void clearAddressing()
Clear addressing.
Definition: ZoneMesh.C:394
Foam::dictionary::subDict
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:631
args
Foam::argList args(argc, argv)
Foam::faceZone::flipMap
const boolList & flipMap() const
Return face flip map.
Definition: faceZone.H:249
Foam::zone::index
label index() const
Return the index of this zone in zone list.
Definition: zone.H:161
Foam::ZoneMesh::checkParallelSync
bool checkParallelSync(const bool report=false) const
Check whether all procs have all zones and in same order. Return.
Definition: ZoneMesh.C:436
Foam::patchIdentifier::name
const word & name() const
Return name.
Definition: patchIdentifier.H:109
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:259
dictIO
IOobject dictIO(dictName, runTime.constant(), mesh, IOobject::MUST_READ_IF_MODIFIED, IOobject::NO_WRITE)
Foam::name
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
Foam::polyMesh::faceNeighbour
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1023
Foam::dictionary::set
void set(entry *)
Assign a new entry, overwrite any existing entry.
Definition: dictionary.C:856