All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
createPatch.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2022 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
33  - Synchronises faces on coupled patches
34  - Synchronises points on coupled patches (optional)
35  - Removes non-constraint patches with no faces
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 "systemDict.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  // Add if
119  // - listed in createPatchDict
120  // - or constraint type
121  // - or non zero size
122  if
123  (
124  addedPatchNames.found(pp.name())
125  || polyPatch::constraintType(pp.type())
126  || returnReduce(pp.size(), sumOp<label>()) > 0
127  )
128  {
129  allPatches.append
130  (
131  pp.clone
132  (
133  patches,
134  allPatches.size(),
135  pp.size(),
136  pp.start()
137  ).ptr()
138  );
139  }
140  else
141  {
142  Info<< "Removing zero-sized patch " << pp.name()
143  << " type " << pp.type()
144  << " at position " << patchi << endl;
145  }
146  }
147  }
148 
149  // Copy non-empty processor patches
150  forAll(patches, patchi)
151  {
152  const polyPatch& pp = patches[patchi];
153 
154  if (isA<processorPolyPatch>(pp))
155  {
156  if (pp.size())
157  {
158  allPatches.append
159  (
160  pp.clone
161  (
162  patches,
163  allPatches.size(),
164  pp.size(),
165  pp.start()
166  ).ptr()
167  );
168  }
169  else
170  {
171  Info<< "Removing empty processor patch " << pp.name()
172  << " at position " << patchi << endl;
173  }
174  }
175  }
176 
177  label nAllPatches = returnReduce(allPatches.size(), sumOp<label>());
178  if (nAllPatches != nOldPatches)
179  {
180  Info<< "Removing patches." << endl;
181  allPatches.shrink();
182  mesh.removeBoundary();
183  mesh.addPatches(allPatches);
184  }
185  else
186  {
187  Info<< "No patches removed." << endl;
188  forAll(allPatches, i)
189  {
190  delete allPatches[i];
191  }
192  }
193 }
194 
195 
196 // Write current match for all patches the as OBJ files
197 void writeCyclicMatchObjs(const fileName& prefix, const polyMesh& mesh)
198 {
199  const polyBoundaryMesh& patches = mesh.boundaryMesh();
200 
201  forAll(patches, patchi)
202  {
203  if
204  (
205  isA<cyclicPolyPatch>(patches[patchi])
206  && refCast<const cyclicPolyPatch>(patches[patchi]).owner()
207  )
208  {
209  const cyclicPolyPatch& cycPatch =
210  refCast<const cyclicPolyPatch>(patches[patchi]);
211 
212  // Write patches
213  {
214  OFstream str(prefix+cycPatch.name() + ".obj");
215  Pout<< "Writing " << cycPatch.name()
216  << " faces to " << str.name() << endl;
218  (
219  str,
220  cycPatch,
221  cycPatch.points()
222  );
223  }
224 
225  const cyclicPolyPatch& nbrPatch = cycPatch.nbrPatch();
226  {
227  OFstream str(prefix+nbrPatch.name()+".obj");
228  Pout<< "Writing " << nbrPatch.name()
229  << " faces to " << str.name() << endl;
231  (
232  str,
233  nbrPatch,
234  nbrPatch.points()
235  );
236  }
237 
238 
239  // Lines between corresponding face centres
240  OFstream str(prefix+cycPatch.name()+nbrPatch.name()+"_match.obj");
241  label vertI = 0;
242 
243  Pout<< "Writing cyclic match as lines between face centres to "
244  << str.name() << endl;
245 
246  forAll(cycPatch, facei)
247  {
248  const point& fc0 = mesh.faceCentres()[cycPatch.start()+facei];
249  meshTools::writeOBJ(str, fc0);
250  vertI++;
251  const point& fc1 = mesh.faceCentres()[nbrPatch.start()+facei];
252  meshTools::writeOBJ(str, fc1);
253  vertI++;
254 
255  str<< "l " << vertI-1 << ' ' << vertI << nl;
256  }
257  }
258  }
259 }
260 
261 
262 // Synchronise points on both sides of coupled boundaries.
263 template<class CombineOp>
264 void syncPoints
265 (
266  const polyMesh& mesh,
267  pointField& points,
268  const CombineOp& cop,
269  const point& nullValue
270 )
271 {
272  if (points.size() != mesh.nPoints())
273  {
275  << "Number of values " << points.size()
276  << " is not equal to the number of points in the mesh "
277  << mesh.nPoints() << abort(FatalError);
278  }
279 
280  const polyBoundaryMesh& patches = mesh.boundaryMesh();
281 
282  // Is there any coupled patch with transformation?
283  bool hasTransformation = false;
284 
285  if (Pstream::parRun())
286  {
287  // Send
288 
289  forAll(patches, patchi)
290  {
291  const polyPatch& pp = patches[patchi];
292 
293  if
294  (
295  isA<processorPolyPatch>(pp)
296  && pp.nPoints() > 0
297  && refCast<const processorPolyPatch>(pp).owner()
298  )
299  {
300  const processorPolyPatch& procPatch =
301  refCast<const processorPolyPatch>(pp);
302 
303  // Get data per patchPoint in neighbouring point numbers.
304  pointField patchInfo(procPatch.nPoints(), nullValue);
305 
306  const labelList& meshPts = procPatch.meshPoints();
307  const labelList& nbrPts = procPatch.nbrPoints();
308 
309  forAll(nbrPts, pointi)
310  {
311  label nbrPointi = nbrPts[pointi];
312  if (nbrPointi >= 0 && nbrPointi < patchInfo.size())
313  {
314  patchInfo[nbrPointi] = points[meshPts[pointi]];
315  }
316  }
317 
318  OPstream toNbr
319  (
321  procPatch.neighbProcNo()
322  );
323  toNbr << patchInfo;
324  }
325  }
326 
327 
328  // Receive and set.
329 
330  forAll(patches, patchi)
331  {
332  const polyPatch& pp = patches[patchi];
333 
334  if
335  (
336  isA<processorPolyPatch>(pp)
337  && pp.nPoints() > 0
338  && !refCast<const processorPolyPatch>(pp).owner()
339  )
340  {
341  const processorPolyPatch& procPatch =
342  refCast<const processorPolyPatch>(pp);
343 
344  pointField nbrPatchInfo(procPatch.nPoints());
345  {
346  // We do not know the number of points on the other side
347  // so cannot use Pstream::read.
348  IPstream fromNbr
349  (
351  procPatch.neighbProcNo()
352  );
353  fromNbr >> nbrPatchInfo;
354  }
355  // Null any value which is not on neighbouring processor
356  nbrPatchInfo.setSize(procPatch.nPoints(), nullValue);
357 
358  if (procPatch.transform().transformsPosition())
359  {
360  hasTransformation = true;
361  procPatch.transform().transformPosition
362  (
363  nbrPatchInfo,
364  nbrPatchInfo
365  );
366  }
367 
368  const labelList& meshPts = procPatch.meshPoints();
369 
370  forAll(meshPts, pointi)
371  {
372  label meshPointi = meshPts[pointi];
373  points[meshPointi] = nbrPatchInfo[pointi];
374  }
375  }
376  }
377  }
378 
379  // Do the cyclics.
380  forAll(patches, patchi)
381  {
382  const polyPatch& pp = patches[patchi];
383 
384  if
385  (
386  isA<cyclicPolyPatch>(pp)
387  && refCast<const cyclicPolyPatch>(pp).owner()
388  )
389  {
390  const cyclicPolyPatch& cycPatch =
391  refCast<const cyclicPolyPatch>(pp);
392 
393  const edgeList& coupledPoints = cycPatch.coupledPoints();
394  const labelList& meshPts = cycPatch.meshPoints();
395  const cyclicPolyPatch& nbrPatch = cycPatch.nbrPatch();
396  const labelList& nbrMeshPts = nbrPatch.meshPoints();
397 
398  pointField patchPoints(coupledPoints.size());
399 
400  forAll(coupledPoints, i)
401  {
402  const edge& e = coupledPoints[i];
403  label point0 = meshPts[e[0]];
404  patchPoints[i] = points[point0];
405  }
406 
407  if (cycPatch.transform().transformsPosition())
408  {
409  hasTransformation = true;
411  (
412  patchPoints,
413  patchPoints
414  );
415  }
416 
417  forAll(coupledPoints, i)
418  {
419  const edge& e = coupledPoints[i];
420  label point1 = nbrMeshPts[e[1]];
421  points[point1] = patchPoints[i];
422  }
423  }
424  }
425 
426  //- Note: hasTransformation is only used for warning messages so
427  // reduction not strictly necessary.
428  // reduce(hasTransformation, orOp<bool>());
429 
430  // Synchronise multiple shared points.
431  const globalMeshData& pd = mesh.globalData();
432 
433  if (pd.nGlobalPoints() > 0)
434  {
435  if (hasTransformation)
436  {
438  << "There are decomposed cyclics in this mesh with"
439  << " transformations." << endl
440  << "This is not supported. The result will be incorrect"
441  << endl;
442  }
443 
444 
445  // Values on shared points.
446  pointField sharedPts(pd.nGlobalPoints(), nullValue);
447 
448  forAll(pd.sharedPointLabels(), i)
449  {
450  label meshPointi = pd.sharedPointLabels()[i];
451  // Fill my entries in the shared points
452  sharedPts[pd.sharedPointAddr()[i]] = points[meshPointi];
453  }
454 
455  // Combine on master.
456  Pstream::listCombineGather(sharedPts, cop);
457  Pstream::listCombineScatter(sharedPts);
458 
459  // Now we will all have the same information. Merge it back with
460  // my local information.
461  forAll(pd.sharedPointLabels(), i)
462  {
463  label meshPointi = pd.sharedPointLabels()[i];
464  points[meshPointi] = sharedPts[pd.sharedPointAddr()[i]];
465  }
466  }
467 }
468 
469 
470 
471 int main(int argc, char *argv[])
472 {
473  #include "addOverwriteOption.H"
474  #include "addRegionOption.H"
475  #include "addDictOption.H"
476  #include "setRootCase.H"
477  #include "createTime.H"
478  runTime.functionObjects().off();
479 
480  Foam::word meshRegionName = polyMesh::defaultRegion;
481  args.optionReadIfPresent("region", meshRegionName);
482 
483  const bool overwrite = args.optionFound("overwrite");
484 
485  #include "createNamedPolyMesh.H"
486 
487  const word oldInstance = mesh.pointsInstance();
488 
489  const dictionary dict(systemDict("createPatchDict", args, mesh));
490 
491  // Whether to synchronise points
492  const Switch pointSync(dict.lookupOrDefault("pointSync", false));
493 
494  // Whether to write cyclic matches to .OBJ files
495  const Switch writeCyclicMatch
496  (
497  dict.lookupOrDefault("writeCyclicMatch", false)
498  );
499 
500 
501  const polyBoundaryMesh& patches = mesh.boundaryMesh();
502 
503  // If running parallel check same patches everywhere
504  patches.checkParallelSync(true);
505 
506 
507  if (writeCyclicMatch)
508  {
509  writeCyclicMatchObjs("initial_", mesh);
510  }
511 
512  // Read patch construct info from dictionary
513  PtrList<dictionary> patchSources(dict.lookup("patches"));
514 
515  HashSet<word> addedPatchNames;
516  forAll(patchSources, addedI)
517  {
518  const dictionary& dict = patchSources[addedI];
519  addedPatchNames.insert(dict.lookup("name"));
520  }
521 
522 
523  // 1. Add all new patches
524  // ~~~~~~~~~~~~~~~~~~~~~~
525 
526  if (patchSources.size())
527  {
528  // Old and new patches.
529  DynamicList<polyPatch*> allPatches(patches.size()+patchSources.size());
530 
531  label startFacei = mesh.nInternalFaces();
532 
533  // Copy old patches.
534  forAll(patches, patchi)
535  {
536  const polyPatch& pp = patches[patchi];
537 
538  if (!isA<processorPolyPatch>(pp))
539  {
540  allPatches.append
541  (
542  pp.clone
543  (
544  patches,
545  patchi,
546  pp.size(),
547  startFacei
548  ).ptr()
549  );
550  startFacei += pp.size();
551  }
552  }
553 
554  forAll(patchSources, addedI)
555  {
556  const dictionary& dict = patchSources[addedI];
557 
558  word patchName(dict.lookup("name"));
559 
560  label destPatchi = patches.findPatchID(patchName);
561 
562  if (destPatchi == -1)
563  {
564  dictionary patchDict(dict.subDict("patchInfo"));
565 
566  destPatchi = allPatches.size();
567 
568  Info<< "Adding new patch " << patchName
569  << " as patch " << destPatchi
570  << " from " << patchDict << endl;
571 
572  patchDict.set("nFaces", 0);
573  patchDict.set("startFace", startFacei);
574 
575  // Add an empty patch.
576  allPatches.append
577  (
579  (
580  patchName,
581  patchDict,
582  destPatchi,
583  patches
584  ).ptr()
585  );
586  }
587  else
588  {
589  Info<< "Patch '" << patchName << "' already exists. Only "
590  << "moving patch faces - type will remain the same" << endl;
591  }
592  }
593 
594  // Copy old patches.
595  forAll(patches, patchi)
596  {
597  const polyPatch& pp = patches[patchi];
598 
599  if (isA<processorPolyPatch>(pp))
600  {
601  allPatches.append
602  (
603  pp.clone
604  (
605  patches,
606  patchi,
607  pp.size(),
608  startFacei
609  ).ptr()
610  );
611  startFacei += pp.size();
612  }
613  }
614 
615  allPatches.shrink();
616  mesh.removeBoundary();
617  mesh.addPatches(allPatches);
618 
619  Info<< endl;
620  }
621 
622 
623 
624  // 2. Repatch faces
625  // ~~~~~~~~~~~~~~~~
626 
627  polyTopoChange meshMod(mesh);
628 
629 
630  forAll(patchSources, addedI)
631  {
632  const dictionary& dict = patchSources[addedI];
633 
634  const word patchName(dict.lookup("name"));
635  label destPatchi = patches.findPatchID(patchName);
636 
637  if (destPatchi == -1)
638  {
640  << "patch " << patchName << " not added. Problem."
641  << abort(FatalError);
642  }
643 
644  const word sourceType(dict.lookup("constructFrom"));
645 
646  if (sourceType == "patches")
647  {
648  labelHashSet patchSources
649  (
650  patches.patchSet
651  (
652  wordReList(dict.lookup("patches"))
653  )
654  );
655 
656  // Repatch faces of the patches.
657  forAllConstIter(labelHashSet, patchSources, iter)
658  {
659  const polyPatch& pp = patches[iter.key()];
660 
661  Info<< "Moving faces from patch " << pp.name()
662  << " to patch " << destPatchi << endl;
663 
664  forAll(pp, i)
665  {
666  changePatchID
667  (
668  mesh,
669  pp.start() + i,
670  destPatchi,
671  meshMod
672  );
673  }
674  }
675  }
676  else if (sourceType == "set")
677  {
678  const word setName(dict.lookup("set"));
679 
680  faceSet faces(mesh, setName);
681 
682  Info<< "Read " << returnReduce(faces.size(), sumOp<label>())
683  << " faces from faceSet " << faces.name() << endl;
684 
685  // Sort (since faceSet contains faces in arbitrary order)
686  labelList faceLabels(faces.toc());
687 
688  SortableList<label> patchFaces(faceLabels);
689 
690  forAll(patchFaces, i)
691  {
692  label facei = patchFaces[i];
693 
694  if (mesh.isInternalFace(facei))
695  {
697  << "Face " << facei << " specified in set "
698  << faces.name()
699  << " is not an external face of the mesh." << endl
700  << "This application can only repatch existing boundary"
701  << " faces." << exit(FatalError);
702  }
703 
704  changePatchID
705  (
706  mesh,
707  facei,
708  destPatchi,
709  meshMod
710  );
711  }
712  }
713  else
714  {
716  << "Invalid source type " << sourceType << endl
717  << "Valid source types are 'patches' 'set'" << exit(FatalError);
718  }
719  }
720  Info<< endl;
721 
722 
723  // Change mesh, use inflation to reforce calculation of transformation
724  // tensors.
725  Info<< "Doing topology modification to order faces." << nl << endl;
726  autoPtr<polyTopoChangeMap> map = meshMod.changeMesh(mesh, true);
727  mesh.setPoints(map().preMotionPoints());
728 
729  if (writeCyclicMatch)
730  {
731  writeCyclicMatchObjs("coupled_", mesh);
732  }
733 
734  // Synchronise points.
735  if (!pointSync)
736  {
737  Info<< "Not synchronising points." << nl << endl;
738  }
739  else
740  {
741  Info<< "Synchronising points." << nl << endl;
742 
743  pointField newPoints(mesh.points());
744 
745  syncPoints
746  (
747  mesh,
748  newPoints,
750  point(great, great, great)
751  );
752 
753  scalarField diff(mag(newPoints-mesh.points()));
754  Info<< "Points changed by average:" << gAverage(diff)
755  << " max:" << gMax(diff) << nl << endl;
756 
757  mesh.setPoints(newPoints);
758  }
759 
760  // 3. Remove zeros-sized patches
761  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
762 
763  Info<< "Removing patches with no faces in them." << nl<< endl;
764  filterPatches(mesh, addedPatchNames);
765 
766 
767  if (writeCyclicMatch)
768  {
769  writeCyclicMatchObjs("final_", mesh);
770  }
771 
772 
773  // Set the precision of the points data to 10
775 
776  if (!overwrite)
777  {
778  runTime++;
779  }
780  else
781  {
782  mesh.setInstance(oldInstance);
783  }
784 
785  // Write resulting mesh
786  Info<< "Writing repatched mesh to " << runTime.timeName() << nl << endl;
787  mesh.write();
788 
789  Info<< "End\n" << endl;
790 
791  return 0;
792 }
793 
794 
795 // ************************************************************************* //
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:453
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:82
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:403
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
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:434
autoPtr< polyTopoChangeMap > 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.
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
const word & name() const
Return name.
A class for handling file names.
Definition: fileName.H:79
A list of face labels.
Definition: faceSet.H:48
const labelList & nbrPoints() const
Return neighbour point labels. WIP.
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.
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
label nInternalFaces() const
A list that is sorted upon construction or when explicitly requested with the sort() method...
Definition: List.H:80
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:458
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
const boolList & flipMap() const
Return face flip map.
Definition: faceZone.H:249
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:325
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
bool optionFound(const word &opt) const
Return true if the named option is found.
Definition: argListI.H:114
A simple wrapper around bool so that it can be read as a word: true/false, on/off, yes/no, y/n, t/f, or none/any.
Definition: Switch.H:60
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:307
virtual const transformer & transform() const
Return transformation between the coupled patches.
const cyclicPolyPatch & nbrPatch() const
IOdictionary systemDict(const word &dictName, const argList &args, const objectRegistry &ob, const word &regionName=polyMesh::defaultRegion)
Definition: systemDict.C:92
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:204
const labelList & meshPoints() const
Return labelList of mesh points in patch. They are constructed.
Neighbour processor patch.
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1211
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:1002
Input inter-processor communications stream.
Definition: IPstream.H:50
static bool constraintType(const word &pt)
Return true if the given type is a constraint type.
Definition: polyPatch.C:265
virtual autoPtr< polyPatch > clone(const polyBoundaryMesh &bm) const
Construct and return a clone, resetting the boundary mesh.
Definition: polyPatch.H:225
label nGlobalPoints() const
Return number of globally shared points.
vector invTransformPosition(const vector &v) const
Inverse transform the given position.
Definition: transformerI.H:177
const fileName & pointsInstance() const
Return the current instance directory for points.
Definition: polyMesh.C:882
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:195
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1237
const Field< PointType > & points() const
Return reference to global points.
const globalMeshData & globalData() const
Return parallel info.
Definition: polyMesh.C:1435
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1224
int neighbProcNo() const
Return neighbour processor number.
errorManip< error > abort(error &err)
Definition: errorManip.H:131
void addPatches(const List< polyPatch *> &, const bool validBoundary=true)
Add boundary patches.
Definition: polyMesh.C:989
defineTemplateTypeNameAndDebug(IOPtrList< ensightPart >, 0)
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:260
Type gMax(const FieldField< Field, Type > &f)
Output inter-processor communications stream.
Definition: OPstream.H:50
void setInstance(const fileName &)
Set the instance for mesh files.
Definition: polyMeshIO.C:77
bool transformsPosition() const
Return true if the transformer transforms a point.
Definition: transformerI.H:146
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
const vectorField & faceCentres() const
label whichZone(const label objectIndex) const
Given a global object index, return the zone it is in.
Definition: MeshZones.C:221
A PtrList of objects of type <T> with automated input and output.
Definition: IOPtrList.H:50
virtual const transformer & transform() const
Return null transform between processor patches.
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:399
label patchi
vector point
Point is a vector.
Definition: point.H:41
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
#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:70
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
const meshFaceZones & faceZones() const
Return face zones.
Definition: polyMesh.H:495
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:306
messageStream Info
dimensioned< scalar > mag(const dimensioned< Type > &)
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)
vector transformPosition(const vector &v) const
Transform the given position.
Definition: transformerI.H:153
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
const doubleScalar e
Elementary charge.
Definition: doubleScalar.H:105
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:76
A subset of mesh faces organised as a primitive patch.
Definition: faceZone.H:65
virtual bool write(const bool write=true) const
Write using setting from DB.
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
virtual void setPoints(const pointField &)
Reset the points.
Definition: polyMesh.C:1299
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:864
label setAction(const topoAction &action)
For compatibility with polyTopoChange: set topological action.