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