createPatch.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  createPatch
26 
27 Description
28  Utility to create patches out of selected boundary faces. Faces come either
29  from existing patches or from a faceSet.
30 
31  More specifically it:
32  - creates new patches (from selected boundary faces). Synchronise faces
33  on coupled patches.
34  - synchronises points on coupled boundaries
35  - remove patches with 0 faces in them
36 
37 \*---------------------------------------------------------------------------*/
38 
39 #include "cyclicPolyPatch.H"
40 #include "syncTools.H"
41 #include "argList.H"
42 #include "polyMesh.H"
43 #include "Time.H"
44 #include "SortableList.H"
45 #include "OFstream.H"
46 #include "meshTools.H"
47 #include "faceSet.H"
48 #include "IOPtrList.H"
49 #include "polyTopoChange.H"
50 #include "polyModifyFace.H"
51 #include "wordReList.H"
52 
53 using namespace Foam;
54 
55 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
56 
57 namespace Foam
58 {
60 }
61 
62 void changePatchID
63 (
64  const polyMesh& mesh,
65  const label faceID,
66  const label patchID,
67  polyTopoChange& meshMod
68 )
69 {
70  const label zoneID = mesh.faceZones().whichZone(faceID);
71 
72  bool zoneFlip = false;
73 
74  if (zoneID >= 0)
75  {
76  const faceZone& fZone = mesh.faceZones()[zoneID];
77 
78  zoneFlip = fZone.flipMap()[fZone.whichFace(faceID)];
79  }
80 
81  meshMod.setAction
82  (
84  (
85  mesh.faces()[faceID], // face
86  faceID, // face ID
87  mesh.faceOwner()[faceID], // owner
88  -1, // neighbour
89  false, // flip flux
90  patchID, // patch ID
91  false, // remove from zone
92  zoneID, // zone ID
93  zoneFlip // zone flip
94  )
95  );
96 }
97 
98 
99 // Filter out the empty patches.
100 void filterPatches(polyMesh& mesh, const HashSet<word>& addedPatchNames)
101 {
103 
104  // Patches to keep
105  DynamicList<polyPatch*> allPatches(patches.size());
106 
107  label nOldPatches = returnReduce(patches.size(), sumOp<label>());
108 
109  // Copy old patches.
110  forAll(patches, patchI)
111  {
112  const polyPatch& pp = patches[patchI];
113 
114  // Note: reduce possible since non-proc patches guaranteed in same order
115  if (!isA<processorPolyPatch>(pp))
116  {
117 
118  // Add if
119  // - non zero size
120  // - or added from the createPatchDict
121  // - or cyclic (since referred to by other cyclic half or
122  // proccyclic)
123 
124  if
125  (
126  addedPatchNames.found(pp.name())
127  || returnReduce(pp.size(), sumOp<label>()) > 0
128  || isA<coupledPolyPatch>(pp)
129  )
130  {
131  allPatches.append
132  (
133  pp.clone
134  (
135  patches,
136  allPatches.size(),
137  pp.size(),
138  pp.start()
139  ).ptr()
140  );
141  }
142  else
143  {
144  Info<< "Removing zero-sized patch " << pp.name()
145  << " type " << pp.type()
146  << " at position " << patchI << endl;
147  }
148  }
149  }
150  // Copy non-empty processor patches
151  forAll(patches, patchI)
152  {
153  const polyPatch& pp = patches[patchI];
154 
155  if (isA<processorPolyPatch>(pp))
156  {
157  if (pp.size())
158  {
159  allPatches.append
160  (
161  pp.clone
162  (
163  patches,
164  allPatches.size(),
165  pp.size(),
166  pp.start()
167  ).ptr()
168  );
169  }
170  else
171  {
172  Info<< "Removing empty processor patch " << pp.name()
173  << " at position " << patchI << endl;
174  }
175  }
176  }
177 
178  label nAllPatches = returnReduce(allPatches.size(), sumOp<label>());
179  if (nAllPatches != nOldPatches)
180  {
181  Info<< "Removing patches." << endl;
182  allPatches.shrink();
184  mesh.addPatches(allPatches);
185  }
186  else
187  {
188  Info<< "No patches removed." << endl;
189  forAll(allPatches, i)
190  {
191  delete allPatches[i];
192  }
193  }
194 }
195 
196 
197 // Dump for all patches the current match
198 void dumpCyclicMatch(const fileName& prefix, const polyMesh& mesh)
199 {
201 
202  forAll(patches, patchI)
203  {
204  if
205  (
206  isA<cyclicPolyPatch>(patches[patchI])
207  && refCast<const cyclicPolyPatch>(patches[patchI]).owner()
208  )
209  {
210  const cyclicPolyPatch& cycPatch =
211  refCast<const cyclicPolyPatch>(patches[patchI]);
212 
213  // Dump patches
214  {
215  OFstream str(prefix+cycPatch.name()+".obj");
216  Pout<< "Dumping " << cycPatch.name()
217  << " faces to " << str.name() << endl;
219  (
220  str,
221  cycPatch,
222  cycPatch.points()
223  );
224  }
225 
226  const cyclicPolyPatch& nbrPatch = cycPatch.neighbPatch();
227  {
228  OFstream str(prefix+nbrPatch.name()+".obj");
229  Pout<< "Dumping " << nbrPatch.name()
230  << " faces to " << str.name() << endl;
232  (
233  str,
234  nbrPatch,
235  nbrPatch.points()
236  );
237  }
238 
239 
240  // Lines between corresponding face centres
241  OFstream str(prefix+cycPatch.name()+nbrPatch.name()+"_match.obj");
242  label vertI = 0;
243 
244  Pout<< "Dumping cyclic match as lines between face centres to "
245  << str.name() << endl;
246 
247  forAll(cycPatch, faceI)
248  {
249  const point& fc0 = mesh.faceCentres()[cycPatch.start()+faceI];
250  meshTools::writeOBJ(str, fc0);
251  vertI++;
252  const point& fc1 = mesh.faceCentres()[nbrPatch.start()+faceI];
253  meshTools::writeOBJ(str, fc1);
254  vertI++;
255 
256  str<< "l " << vertI-1 << ' ' << vertI << nl;
257  }
258  }
259  }
260 }
261 
262 
263 void separateList
264 (
265  const vectorField& separation,
266  UList<vector>& field
267 )
268 {
269  if (separation.size() == 1)
270  {
271  // Single value for all.
272 
273  forAll(field, i)
274  {
275  field[i] += separation[0];
276  }
277  }
278  else if (separation.size() == field.size())
279  {
280  forAll(field, i)
281  {
282  field[i] += separation[i];
283  }
284  }
285  else
286  {
288  << "Sizes of field and transformation not equal. field:"
289  << field.size() << " transformation:" << separation.size()
290  << abort(FatalError);
291  }
292 }
293 
294 
295 // Synchronise points on both sides of coupled boundaries.
296 template<class CombineOp>
297 void syncPoints
298 (
299  const polyMesh& mesh,
301  const CombineOp& cop,
302  const point& nullValue
303 )
304 {
305  if (points.size() != mesh.nPoints())
306  {
308  << "Number of values " << points.size()
309  << " is not equal to the number of points in the mesh "
310  << mesh.nPoints() << abort(FatalError);
311  }
312 
314 
315  // Is there any coupled patch with transformation?
316  bool hasTransformation = false;
317 
318  if (Pstream::parRun())
319  {
320  // Send
321 
322  forAll(patches, patchI)
323  {
324  const polyPatch& pp = patches[patchI];
325 
326  if
327  (
328  isA<processorPolyPatch>(pp)
329  && pp.nPoints() > 0
330  && refCast<const processorPolyPatch>(pp).owner()
331  )
332  {
333  const processorPolyPatch& procPatch =
334  refCast<const processorPolyPatch>(pp);
335 
336  // Get data per patchPoint in neighbouring point numbers.
337  pointField patchInfo(procPatch.nPoints(), nullValue);
338 
339  const labelList& meshPts = procPatch.meshPoints();
340  const labelList& nbrPts = procPatch.neighbPoints();
341 
342  forAll(nbrPts, pointI)
343  {
344  label nbrPointI = nbrPts[pointI];
345  if (nbrPointI >= 0 && nbrPointI < patchInfo.size())
346  {
347  patchInfo[nbrPointI] = points[meshPts[pointI]];
348  }
349  }
350 
351  OPstream toNbr(Pstream::blocking, procPatch.neighbProcNo());
352  toNbr << patchInfo;
353  }
354  }
355 
356 
357  // Receive and set.
358 
359  forAll(patches, patchI)
360  {
361  const polyPatch& pp = patches[patchI];
362 
363  if
364  (
365  isA<processorPolyPatch>(pp)
366  && pp.nPoints() > 0
367  && !refCast<const processorPolyPatch>(pp).owner()
368  )
369  {
370  const processorPolyPatch& procPatch =
371  refCast<const processorPolyPatch>(pp);
372 
373  pointField nbrPatchInfo(procPatch.nPoints());
374  {
375  // We do not know the number of points on the other side
376  // so cannot use Pstream::read.
377  IPstream fromNbr
378  (
380  procPatch.neighbProcNo()
381  );
382  fromNbr >> nbrPatchInfo;
383  }
384  // Null any value which is not on neighbouring processor
385  nbrPatchInfo.setSize(procPatch.nPoints(), nullValue);
386 
387  if (!procPatch.parallel())
388  {
389  hasTransformation = true;
390  transformList(procPatch.forwardT(), nbrPatchInfo);
391  }
392  else if (procPatch.separated())
393  {
394  hasTransformation = true;
395  separateList(-procPatch.separation(), nbrPatchInfo);
396  }
397 
398  const labelList& meshPts = procPatch.meshPoints();
399 
400  forAll(meshPts, pointI)
401  {
402  label meshPointI = meshPts[pointI];
403  points[meshPointI] = nbrPatchInfo[pointI];
404  }
405  }
406  }
407  }
408 
409  // Do the cyclics.
410  forAll(patches, patchI)
411  {
412  const polyPatch& pp = patches[patchI];
413 
414  if
415  (
416  isA<cyclicPolyPatch>(pp)
417  && refCast<const cyclicPolyPatch>(pp).owner()
418  )
419  {
420  const cyclicPolyPatch& cycPatch =
421  refCast<const cyclicPolyPatch>(pp);
422 
423  const edgeList& coupledPoints = cycPatch.coupledPoints();
424  const labelList& meshPts = cycPatch.meshPoints();
425  const cyclicPolyPatch& nbrPatch = cycPatch.neighbPatch();
426  const labelList& nbrMeshPts = nbrPatch.meshPoints();
427 
428  pointField half0Values(coupledPoints.size());
429 
430  forAll(coupledPoints, i)
431  {
432  const edge& e = coupledPoints[i];
433  label point0 = meshPts[e[0]];
434  half0Values[i] = points[point0];
435  }
436 
437  if (!cycPatch.parallel())
438  {
439  hasTransformation = true;
440  transformList(cycPatch.reverseT(), half0Values);
441  }
442  else if (cycPatch.separated())
443  {
444  hasTransformation = true;
445  separateList(cycPatch.separation(), half0Values);
446  }
447 
448  forAll(coupledPoints, i)
449  {
450  const edge& e = coupledPoints[i];
451  label point1 = nbrMeshPts[e[1]];
452  points[point1] = half0Values[i];
453  }
454  }
455  }
456 
457  //- Note: hasTransformation is only used for warning messages so
458  // reduction not strictly nessecary.
459  //reduce(hasTransformation, orOp<bool>());
460 
461  // Synchronize multiple shared points.
462  const globalMeshData& pd = mesh.globalData();
463 
464  if (pd.nGlobalPoints() > 0)
465  {
466  if (hasTransformation)
467  {
469  << "There are decomposed cyclics in this mesh with"
470  << " transformations." << endl
471  << "This is not supported. The result will be incorrect"
472  << endl;
473  }
474 
475 
476  // Values on shared points.
477  pointField sharedPts(pd.nGlobalPoints(), nullValue);
478 
479  forAll(pd.sharedPointLabels(), i)
480  {
481  label meshPointI = pd.sharedPointLabels()[i];
482  // Fill my entries in the shared points
483  sharedPts[pd.sharedPointAddr()[i]] = points[meshPointI];
484  }
485 
486  // Combine on master.
487  Pstream::listCombineGather(sharedPts, cop);
488  Pstream::listCombineScatter(sharedPts);
489 
490  // Now we will all have the same information. Merge it back with
491  // my local information.
492  forAll(pd.sharedPointLabels(), i)
493  {
494  label meshPointI = pd.sharedPointLabels()[i];
495  points[meshPointI] = sharedPts[pd.sharedPointAddr()[i]];
496  }
497  }
498 }
499 
500 
501 
502 int main(int argc, char *argv[])
503 {
504  #include "addOverwriteOption.H"
505  #include "addRegionOption.H"
506  #include "addDictOption.H"
508  (
509  "writeObj",
510  "write obj files showing the cyclic matching process"
511  );
512  #include "setRootCase.H"
513  #include "createTime.H"
514  runTime.functionObjects().off();
515 
516  Foam::word meshRegionName = polyMesh::defaultRegion;
517  args.optionReadIfPresent("region", meshRegionName);
518 
519  const bool overwrite = args.optionFound("overwrite");
520 
521  #include "createNamedPolyMesh.H"
522 
523  const bool writeObj = args.optionFound("writeObj");
524 
525  const word oldInstance = mesh.pointsInstance();
526 
527  const word dictName("createPatchDict");
528  #include "setSystemMeshDictionaryIO.H"
529  Info<< "Reading " << dictIO.instance()/dictIO.name() << nl << endl;
530 
532 
533  // Whether to synchronise points
534  const Switch pointSync(dict.lookup("pointSync"));
535 
536 
538 
539  // If running parallel check same patches everywhere
540  patches.checkParallelSync(true);
541 
542 
543  if (writeObj)
544  {
545  dumpCyclicMatch("initial_", mesh);
546  }
547 
548  // Read patch construct info from dictionary
549  PtrList<dictionary> patchSources(dict.lookup("patches"));
550 
551  HashSet<word> addedPatchNames;
552  forAll(patchSources, addedI)
553  {
554  const dictionary& dict = patchSources[addedI];
555  addedPatchNames.insert(dict.lookup("name"));
556  }
557 
558 
559  // 1. Add all new patches
560  // ~~~~~~~~~~~~~~~~~~~~~~
561 
562  if (patchSources.size())
563  {
564  // Old and new patches.
565  DynamicList<polyPatch*> allPatches(patches.size()+patchSources.size());
566 
567  label startFaceI = mesh.nInternalFaces();
568 
569  // Copy old patches.
570  forAll(patches, patchI)
571  {
572  const polyPatch& pp = patches[patchI];
573 
574  if (!isA<processorPolyPatch>(pp))
575  {
576  allPatches.append
577  (
578  pp.clone
579  (
580  patches,
581  patchI,
582  pp.size(),
583  startFaceI
584  ).ptr()
585  );
586  startFaceI += pp.size();
587  }
588  }
589 
590  forAll(patchSources, addedI)
591  {
592  const dictionary& dict = patchSources[addedI];
593 
594  word patchName(dict.lookup("name"));
595 
596  label destPatchI = patches.findPatchID(patchName);
597 
598  if (destPatchI == -1)
599  {
600  dictionary patchDict(dict.subDict("patchInfo"));
601 
602  destPatchI = allPatches.size();
603 
604  Info<< "Adding new patch " << patchName
605  << " as patch " << destPatchI
606  << " from " << patchDict << endl;
607 
608  patchDict.set("nFaces", 0);
609  patchDict.set("startFace", startFaceI);
610 
611  // Add an empty patch.
612  allPatches.append
613  (
615  (
616  patchName,
617  patchDict,
618  destPatchI,
619  patches
620  ).ptr()
621  );
622  }
623  else
624  {
625  Info<< "Patch '" << patchName << "' already exists. Only "
626  << "moving patch faces - type will remain the same" << endl;
627  }
628  }
629 
630  // Copy old patches.
631  forAll(patches, patchI)
632  {
633  const polyPatch& pp = patches[patchI];
634 
635  if (isA<processorPolyPatch>(pp))
636  {
637  allPatches.append
638  (
639  pp.clone
640  (
641  patches,
642  patchI,
643  pp.size(),
644  startFaceI
645  ).ptr()
646  );
647  startFaceI += pp.size();
648  }
649  }
650 
651  allPatches.shrink();
653  mesh.addPatches(allPatches);
654 
655  Info<< endl;
656  }
657 
658 
659 
660  // 2. Repatch faces
661  // ~~~~~~~~~~~~~~~~
662 
663  polyTopoChange meshMod(mesh);
664 
665 
666  forAll(patchSources, addedI)
667  {
668  const dictionary& dict = patchSources[addedI];
669 
670  const word patchName(dict.lookup("name"));
671  label destPatchI = patches.findPatchID(patchName);
672 
673  if (destPatchI == -1)
674  {
676  << "patch " << patchName << " not added. Problem."
677  << abort(FatalError);
678  }
679 
680  const word sourceType(dict.lookup("constructFrom"));
681 
682  if (sourceType == "patches")
683  {
684  labelHashSet patchSources
685  (
686  patches.patchSet
687  (
688  wordReList(dict.lookup("patches"))
689  )
690  );
691 
692  // Repatch faces of the patches.
693  forAllConstIter(labelHashSet, patchSources, iter)
694  {
695  const polyPatch& pp = patches[iter.key()];
696 
697  Info<< "Moving faces from patch " << pp.name()
698  << " to patch " << destPatchI << endl;
699 
700  forAll(pp, i)
701  {
703  (
704  mesh,
705  pp.start() + i,
706  destPatchI,
707  meshMod
708  );
709  }
710  }
711  }
712  else if (sourceType == "set")
713  {
714  const word setName(dict.lookup("set"));
715 
716  faceSet faces(mesh, setName);
717 
718  Info<< "Read " << returnReduce(faces.size(), sumOp<label>())
719  << " faces from faceSet " << faces.name() << endl;
720 
721  // Sort (since faceSet contains faces in arbitrary order)
722  labelList faceLabels(faces.toc());
723 
724  SortableList<label> patchFaces(faceLabels);
725 
726  forAll(patchFaces, i)
727  {
728  label faceI = patchFaces[i];
729 
730  if (mesh.isInternalFace(faceI))
731  {
733  << "Face " << faceI << " specified in set "
734  << faces.name()
735  << " is not an external face of the mesh." << endl
736  << "This application can only repatch existing boundary"
737  << " faces." << exit(FatalError);
738  }
739 
741  (
742  mesh,
743  faceI,
744  destPatchI,
745  meshMod
746  );
747  }
748  }
749  else
750  {
752  << "Invalid source type " << sourceType << endl
753  << "Valid source types are 'patches' 'set'" << exit(FatalError);
754  }
755  }
756  Info<< endl;
757 
758 
759  // Change mesh, use inflation to reforce calculation of transformation
760  // tensors.
761  Info<< "Doing topology modification to order faces." << nl << endl;
762  autoPtr<mapPolyMesh> map = meshMod.changeMesh(mesh, true);
763  mesh.movePoints(map().preMotionPoints());
764 
765  if (writeObj)
766  {
767  dumpCyclicMatch("coupled_", mesh);
768  }
769 
770  // Synchronise points.
771  if (!pointSync)
772  {
773  Info<< "Not synchronising points." << nl << endl;
774  }
775  else
776  {
777  Info<< "Synchronising points." << nl << endl;
778 
779  // This is a bit tricky. Both normal and position might be out and
780  // current separation also includes the normal
781  // ( separation_ = (nf&(Cr - Cf))*nf ).
782 
783  // For cyclic patches:
784  // - for separated ones use user specified offset vector
785 
786  forAll(mesh.boundaryMesh(), patchI)
787  {
788  const polyPatch& pp = mesh.boundaryMesh()[patchI];
789 
790  if (pp.size() && isA<coupledPolyPatch>(pp))
791  {
792  const coupledPolyPatch& cpp =
793  refCast<const coupledPolyPatch>(pp);
794 
795  if (cpp.separated())
796  {
797  Info<< "On coupled patch " << pp.name()
798  << " separation[0] was "
799  << cpp.separation()[0] << endl;
800 
801  if (isA<cyclicPolyPatch>(pp) && pp.size())
802  {
803  const cyclicPolyPatch& cycpp =
804  refCast<const cyclicPolyPatch>(pp);
805 
807  {
808  // Force to wanted separation
809  Info<< "On cyclic translation patch " << pp.name()
810  << " forcing uniform separation of "
811  << cycpp.separationVector() << endl;
812  const_cast<vectorField&>(cpp.separation()) =
813  pointField(1, cycpp.separationVector());
814  }
815  else
816  {
817  const cyclicPolyPatch& nbr = cycpp.neighbPatch();
818  const_cast<vectorField&>(cpp.separation()) =
819  pointField
820  (
821  1,
822  nbr[0].centre(mesh.points())
823  - cycpp[0].centre(mesh.points())
824  );
825  }
826  }
827  Info<< "On coupled patch " << pp.name()
828  << " forcing uniform separation of "
829  << cpp.separation() << endl;
830  }
831  else if (!cpp.parallel())
832  {
833  Info<< "On coupled patch " << pp.name()
834  << " forcing uniform rotation of "
835  << cpp.forwardT()[0] << endl;
836 
837  const_cast<tensorField&>
838  (
839  cpp.forwardT()
840  ).setSize(1);
841  const_cast<tensorField&>
842  (
843  cpp.reverseT()
844  ).setSize(1);
845 
846  Info<< "On coupled patch " << pp.name()
847  << " forcing uniform rotation of "
848  << cpp.forwardT() << endl;
849  }
850  }
851  }
852 
853  Info<< "Synchronising points." << endl;
854 
855  pointField newPoints(mesh.points());
856 
857  syncPoints
858  (
859  mesh,
860  newPoints,
862  point(GREAT, GREAT, GREAT)
863  );
864 
865  scalarField diff(mag(newPoints-mesh.points()));
866  Info<< "Points changed by average:" << gAverage(diff)
867  << " max:" << gMax(diff) << nl << endl;
868 
869  mesh.movePoints(newPoints);
870  }
871 
872  // 3. Remove zeros-sized patches
873  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
874 
875  Info<< "Removing patches with no faces in them." << nl<< endl;
876  filterPatches(mesh, addedPatchNames);
877 
878 
879  if (writeObj)
880  {
881  dumpCyclicMatch("final_", mesh);
882  }
883 
884 
885  // Set the precision of the points data to 10
887 
888  if (!overwrite)
889  {
890  runTime++;
891  }
892  else
893  {
894  mesh.setInstance(oldInstance);
895  }
896 
897  // Write resulting mesh
898  Info<< "Writing repatched mesh to " << runTime.timeName() << nl << endl;
899  mesh.write();
900 
901  Info<< "End\n" << endl;
902 
903  return 0;
904 }
905 
906 
907 // ************************************************************************* //
Foam::pointField
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
Foam::IOdictionary
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:53
Foam::polyMesh::points
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:979
syncPoints
void syncPoints(const polyMesh &mesh, pointField &points, const CombineOp &cop, const point &nullValue)
Definition: createPatch.C:298
dumpCyclicMatch
void dumpCyclicMatch(const fileName &prefix, const polyMesh &mesh)
Definition: createPatch.C:198
setSize
points setSize(newPointi)
meshTools.H
Foam::coupledPolyPatch::separation
virtual const vectorField & separation() const
If the planes are separated the separation vector.
Definition: coupledPolyPatch.H:282
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::PrimitivePatch::points
const Field< PointType > & points() const
Return reference to global points.
Definition: PrimitivePatchTemplate.H:282
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::fileName
A class for handling file names.
Definition: fileName.H:69
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::IOPtrList
A PtrList of objects of type <T> with automated input and output.
Definition: IOPtrList.H:50
Foam::polyBoundaryMesh
Foam::polyBoundaryMesh.
Definition: polyBoundaryMesh.H:60
Foam::polyMesh::defaultRegion
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:306
cyclicPolyPatch.H
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::cyclicPolyPatch
Cyclic plane patch.
Definition: cyclicPolyPatch.H:63
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
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:56
Foam::gAverage
Type gAverage(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:571
Foam::meshTools::writeOBJ
void writeOBJ(Ostream &os, const point &pt)
Write obj representation of point.
Definition: meshTools.C:203
Foam::OPstream
Output inter-processor communications stream.
Definition: OPstream.H:50
Foam::cyclicPolyPatch::neighbPatch
const cyclicPolyPatch & neighbPatch() const
Definition: cyclicPolyPatch.H:328
addOverwriteOption.H
Foam::coupledPolyPatch::forwardT
virtual const tensorField & forwardT() const
Return face transformation tensor.
Definition: coupledPolyPatch.H:294
main
int main(int argc, char *argv[])
Definition: createPatch.C:561
polyTopoChange.H
Foam::UPstream::parRun
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:377
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::edge
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:58
Foam::coupledPolyPatch::reverseT
virtual const tensorField & reverseT() const
Return neighbour-cell transformation tensor.
Definition: coupledPolyPatch.H:300
Foam::argList::addBoolOption
static void addBoolOption(const word &opt, const string &usage="")
Add to a bool option to validOptions with usage information.
Definition: argList.C:98
Foam::processorPolyPatch::neighbProcNo
int neighbProcNo() const
Return neigbour processor number.
Definition: processorPolyPatch.H:255
Foam::faceSet
A list of face labels.
Definition: faceSet.H:48
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
Foam::coupledPolyPatch
The coupledPolyPatch is an abstract base class for patches that couple regions of the computational d...
Definition: coupledPolyPatch.H:51
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
separateList
void separateList(const vectorField &separation, UList< vector > &field)
Definition: createPatch.C:264
Foam::IOobject::instance
const fileName & instance() const
Definition: IOobject.H:350
polyMesh.H
patchFaces
labelList patchFaces(const polyBoundaryMesh &patches, const wordList &names)
Definition: extrudeMesh.C:148
Foam::HashSet
A HashTable with keys but without contents.
Definition: HashSet.H:59
syncTools.H
setSystemMeshDictionaryIO.H
createNamedPolyMesh.H
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
forAllConstIter
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:39
OFstream.H
Foam::polyMesh::faceZones
const faceZoneMesh & faceZones() const
Return face zone mesh.
Definition: polyMesh.H:463
Foam::cyclicPolyPatch::coupledPoints
const edgeList & coupledPoints() const
Return connected points (from patch local to neighbour patch local)
Definition: cyclicPolyPatch.C:1012
addDictOption.H
Foam::coupledPolyPatch::parallel
virtual bool parallel() const
Are the cyclic planes parallel.
Definition: coupledPolyPatch.H:288
Foam::fvMesh::write
virtual bool write() const
Write mesh using IO settings from time.
Definition: fvMesh.C:873
Foam::diff
scalar diff(const triad &A, const triad &B)
Return a quantity of the difference between two triads.
Definition: triad.C:407
dictName
const word dictName("particleTrackDict")
Foam::polyMesh::pointsInstance
const fileName & pointsInstance() const
Return the current instance directory for points.
Definition: polyMesh.C:772
Foam::UPstream::blocking
@ blocking
Definition: UPstream.H:66
SortableList.H
Foam::coupledPolyPatch::separated
virtual bool separated() const
Are the planes separated.
Definition: coupledPolyPatch.H:276
Foam::primitiveMesh::nPoints
label nPoints() const
Definition: primitiveMeshI.H:35
Foam::globalMeshData::nGlobalPoints
label nGlobalPoints() const
Return number of globally shared points.
Definition: globalMeshData.C:1986
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::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::faceZone
A subset of mesh faces organised as a primitive patch.
Definition: faceZone.H:64
IOPtrList.H
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
Foam::OSstream::name
virtual const fileName & name() const
Return the name of the stream.
Definition: OSstream.H:88
argList.H
Foam::polyMesh::faceOwner
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1017
faceSet.H
addRegionOption.H
Foam::IOobject::name
const word & name() const
Return name.
Definition: IOobject.H:273
Foam::processorPolyPatch
Neighbour processor patch.
Definition: processorPolyPatch.H:55
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
Foam::ZoneMesh::whichZone
label whichZone(const label objectIndex) const
Given a global object index, return the zone it is in.
Definition: ZoneMesh.C:226
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::PrimitivePatch::nPoints
label nPoints() const
Return number of points supporting patch faces.
Definition: PrimitivePatchTemplate.H:293
Foam::FatalError
error FatalError
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
Foam::globalMeshData
Various mesh related information for a parallel run. Upon construction, constructs all info using par...
Definition: globalMeshData.H:106
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
Foam::SortableList
A list that is sorted upon construction or when explicitly requested with the sort() method.
Definition: List.H:65
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
Foam::HashTable< nil, word, string::hash >::found
bool found(const Key &) const
Return true if hashedEntry is found in table.
Definition: HashTable.C:109
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
writeObj
void writeObj(Ostream &os, const pointField &points)
Definition: Test-PrimitivePatch.C:44
Foam::faceZone::whichFace
label whichFace(const label globalCellID) const
Helper function to re-direct to zone::localID(...)
Definition: faceZone.C:314
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:131
Foam::e
const double e
Elementary charge.
Definition: doubleFloat.H:94
Foam::polyPatch::start
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:312
Foam::coupledPolyPatch::transform
virtual transformType transform() const
Type of transform.
Definition: coupledPolyPatch.H:256
polyModifyFace.H
Foam::polyPatch::clone
virtual autoPtr< polyPatch > clone(const polyBoundaryMesh &bm) const
Construct and return a clone, resetting the boundary mesh.
Definition: polyPatch.H:231
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Foam::polyMesh::removeBoundary
void removeBoundary()
Remove boundary patches.
Definition: polyMeshClear.C:36
Foam::OFstream
Output to file stream.
Definition: OFstream.H:81
Foam::polyMesh::setInstance
void setInstance(const fileName &)
Set the instance for mesh files.
Definition: polyMeshIO.C:32
Foam::Pstream::listCombineGather
static void listCombineGather(const List< commsStruct > &comms, List< T > &Value, const CombineOp &cop, const int tag, const label comm)
Definition: combineGatherScatter.C:282
Foam::max
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
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
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
Foam::processorPolyPatch::neighbPoints
const labelList & neighbPoints() const
Return neighbour point labels. WIP.
Definition: processorPolyPatch.C:462
Foam::polyTopoChange::setAction
label setAction(const topoAction &action)
For compatibility with polyTopoChange: set topological action.
Definition: polyTopoChange.C:2530
Foam::sumOp
Definition: ops.H:162
Foam::IOstream::defaultPrecision
static unsigned int defaultPrecision()
Return the default precision.
Definition: IOstream.H:461
Foam::Pout
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
Foam::Vector< scalar >
wordReList.H
Foam::polyMesh::addPatches
void addPatches(const List< polyPatch * > &, const bool validBoundary=true)
Add boundary patches.
Definition: polyMesh.C:878
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::primitiveMesh::isInternalFace
bool isInternalFace(const label faceIndex) const
Return true if given face label is internal to the mesh.
Definition: primitiveMeshI.H:70
Foam::primitiveMesh::faceCentres
const vectorField & faceCentres() const
Definition: primitiveMeshFaceCentresAndAreas.C:130
Foam::UList
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: HashTable.H:60
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::HashSet::insert
bool insert(const Key &key)
Insert a new entry.
Definition: HashSet.H:116
createTime.H
Foam::argList::optionFound
bool optionFound(const word &opt) const
Return true if the named option is found.
Definition: argListI.H:108
Foam::globalMeshData::sharedPointLabels
const labelList & sharedPointLabels() const
Return indices of local points that are globally shared.
Definition: globalMeshData.C:1996
Foam::globalMeshData::sharedPointAddr
const labelList & sharedPointAddr() const
Return addressing into the complete globally shared points.
Definition: globalMeshData.C:2006
patches
patches[0]
Definition: createSingleCellMesh.H:36
Foam::IPstream
Input inter-processor communications stream.
Definition: IPstream.H:50
Foam::minMagSqrEqOp
Definition: ops.H:79
changePatchID
void changePatchID(const polyMesh &mesh, const label faceID, const label patchID, polyTopoChange &meshMod)
Definition: createPatch.C:63
Foam::List::size
void size(const label)
Override size to be inconsistent with allocated storage.
Foam::transformList
void transformList(const tensor &, UList< T > &)
Apply transformation to list. Either single transformation tensor.
Definition: transformList.C:49
Foam::polyMesh::globalData
const globalMeshData & globalData() const
Return parallel info.
Definition: polyMesh.C:1143
Foam::Pstream::listCombineScatter
static void listCombineScatter(const List< commsStruct > &comms, List< T > &Value, const int tag, const label comm)
Scatter data. Reverse of combineGather.
Definition: combineGatherScatter.C:417
Foam::point
vector point
Point is a vector.
Definition: point.H:41
Foam::UList::size
label size() const
Return the number of elements in the UList.
Definition: UListI.H:299
Foam::PrimitivePatch::meshPoints
const labelList & meshPoints() const
Return labelList of mesh points in patch.
Definition: PrimitivePatchTemplate.C:392
Foam::defineTemplateTypeNameAndDebug
defineTemplateTypeNameAndDebug(IOPtrList< ensightPart >, 0)
Foam::dictionary::subDict
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:631
Foam::coupledPolyPatch::TRANSLATIONAL
@ TRANSLATIONAL
Definition: coupledPolyPatch.H:61
args
Foam::argList args(argc, argv)
filterPatches
void filterPatches(polyMesh &mesh, const HashSet< word > &addedPatchNames)
Definition: createPatch.C:100
Foam::faceZone::flipMap
const boolList & flipMap() const
Return face flip map.
Definition: faceZone.H:249
Foam::argList::optionReadIfPresent
bool optionReadIfPresent(const word &opt, T &) const
Read a value from the named option if present.
Definition: argListI.H:198
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
Foam::cyclicPolyPatch::separationVector
const vector & separationVector() const
Translation vector for translational cyclics.
Definition: cyclicPolyPatch.H:384
Foam::gMax
Type gMax(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:562
dictIO
IOobject dictIO(dictName, runTime.constant(), mesh, IOobject::MUST_READ_IF_MODIFIED, IOobject::NO_WRITE)