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-2021 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 Application
25  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 "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 
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 // Write current match for all patches the as OBJ files
199 void writeCyclicMatchObjs(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  // Write patches
215  {
216  OFstream str(prefix+cycPatch.name() + ".obj");
217  Pout<< "Writing " << cycPatch.name()
218  << " faces to " << str.name() << endl;
220  (
221  str,
222  cycPatch,
223  cycPatch.points()
224  );
225  }
226 
227  const cyclicPolyPatch& nbrPatch = cycPatch.nbrPatch();
228  {
229  OFstream str(prefix+nbrPatch.name()+".obj");
230  Pout<< "Writing " << 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<< "Writing 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 // Synchronise points on both sides of coupled boundaries.
265 template<class CombineOp>
266 void syncPoints
267 (
268  const polyMesh& mesh,
269  pointField& points,
270  const CombineOp& cop,
271  const point& nullValue
272 )
273 {
274  if (points.size() != mesh.nPoints())
275  {
277  << "Number of values " << points.size()
278  << " is not equal to the number of points in the mesh "
279  << mesh.nPoints() << abort(FatalError);
280  }
281 
282  const polyBoundaryMesh& patches = mesh.boundaryMesh();
283 
284  // Is there any coupled patch with transformation?
285  bool hasTransformation = false;
286 
287  if (Pstream::parRun())
288  {
289  // Send
290 
291  forAll(patches, patchi)
292  {
293  const polyPatch& pp = patches[patchi];
294 
295  if
296  (
297  isA<processorPolyPatch>(pp)
298  && pp.nPoints() > 0
299  && refCast<const processorPolyPatch>(pp).owner()
300  )
301  {
302  const processorPolyPatch& procPatch =
303  refCast<const processorPolyPatch>(pp);
304 
305  // Get data per patchPoint in neighbouring point numbers.
306  pointField patchInfo(procPatch.nPoints(), nullValue);
307 
308  const labelList& meshPts = procPatch.meshPoints();
309  const labelList& nbrPts = procPatch.nbrPoints();
310 
311  forAll(nbrPts, pointi)
312  {
313  label nbrPointi = nbrPts[pointi];
314  if (nbrPointi >= 0 && nbrPointi < patchInfo.size())
315  {
316  patchInfo[nbrPointi] = points[meshPts[pointi]];
317  }
318  }
319 
320  OPstream toNbr
321  (
323  procPatch.neighbProcNo()
324  );
325  toNbr << patchInfo;
326  }
327  }
328 
329 
330  // Receive and set.
331 
332  forAll(patches, patchi)
333  {
334  const polyPatch& pp = patches[patchi];
335 
336  if
337  (
338  isA<processorPolyPatch>(pp)
339  && pp.nPoints() > 0
340  && !refCast<const processorPolyPatch>(pp).owner()
341  )
342  {
343  const processorPolyPatch& procPatch =
344  refCast<const processorPolyPatch>(pp);
345 
346  pointField nbrPatchInfo(procPatch.nPoints());
347  {
348  // We do not know the number of points on the other side
349  // so cannot use Pstream::read.
350  IPstream fromNbr
351  (
353  procPatch.neighbProcNo()
354  );
355  fromNbr >> nbrPatchInfo;
356  }
357  // Null any value which is not on neighbouring processor
358  nbrPatchInfo.setSize(procPatch.nPoints(), nullValue);
359 
360  if (procPatch.transform().transformsPosition())
361  {
362  hasTransformation = true;
363  procPatch.transform().transformPosition
364  (
365  nbrPatchInfo,
366  nbrPatchInfo
367  );
368  }
369 
370  const labelList& meshPts = procPatch.meshPoints();
371 
372  forAll(meshPts, pointi)
373  {
374  label meshPointi = meshPts[pointi];
375  points[meshPointi] = nbrPatchInfo[pointi];
376  }
377  }
378  }
379  }
380 
381  // Do the cyclics.
382  forAll(patches, patchi)
383  {
384  const polyPatch& pp = patches[patchi];
385 
386  if
387  (
388  isA<cyclicPolyPatch>(pp)
389  && refCast<const cyclicPolyPatch>(pp).owner()
390  )
391  {
392  const cyclicPolyPatch& cycPatch =
393  refCast<const cyclicPolyPatch>(pp);
394 
395  const edgeList& coupledPoints = cycPatch.coupledPoints();
396  const labelList& meshPts = cycPatch.meshPoints();
397  const cyclicPolyPatch& nbrPatch = cycPatch.nbrPatch();
398  const labelList& nbrMeshPts = nbrPatch.meshPoints();
399 
400  pointField patchPoints(coupledPoints.size());
401 
402  forAll(coupledPoints, i)
403  {
404  const edge& e = coupledPoints[i];
405  label point0 = meshPts[e[0]];
406  patchPoints[i] = points[point0];
407  }
408 
409  if (cycPatch.transform().transformsPosition())
410  {
411  hasTransformation = true;
413  (
414  patchPoints,
415  patchPoints
416  );
417  }
418 
419  forAll(coupledPoints, i)
420  {
421  const edge& e = coupledPoints[i];
422  label point1 = nbrMeshPts[e[1]];
423  points[point1] = patchPoints[i];
424  }
425  }
426  }
427 
428  //- Note: hasTransformation is only used for warning messages so
429  // reduction not strictly necessary.
430  // reduce(hasTransformation, orOp<bool>());
431 
432  // Synchronise multiple shared points.
433  const globalMeshData& pd = mesh.globalData();
434 
435  if (pd.nGlobalPoints() > 0)
436  {
437  if (hasTransformation)
438  {
440  << "There are decomposed cyclics in this mesh with"
441  << " transformations." << endl
442  << "This is not supported. The result will be incorrect"
443  << endl;
444  }
445 
446 
447  // Values on shared points.
448  pointField sharedPts(pd.nGlobalPoints(), nullValue);
449 
450  forAll(pd.sharedPointLabels(), i)
451  {
452  label meshPointi = pd.sharedPointLabels()[i];
453  // Fill my entries in the shared points
454  sharedPts[pd.sharedPointAddr()[i]] = points[meshPointi];
455  }
456 
457  // Combine on master.
458  Pstream::listCombineGather(sharedPts, cop);
459  Pstream::listCombineScatter(sharedPts);
460 
461  // Now we will all have the same information. Merge it back with
462  // my local information.
463  forAll(pd.sharedPointLabels(), i)
464  {
465  label meshPointi = pd.sharedPointLabels()[i];
466  points[meshPointi] = sharedPts[pd.sharedPointAddr()[i]];
467  }
468  }
469 }
470 
471 
472 
473 int main(int argc, char *argv[])
474 {
475  #include "addOverwriteOption.H"
476  #include "addRegionOption.H"
477  #include "addDictOption.H"
478  #include "setRootCase.H"
479  #include "createTime.H"
481 
482  Foam::word meshRegionName = polyMesh::defaultRegion;
483  args.optionReadIfPresent("region", meshRegionName);
484 
485  const bool overwrite = args.optionFound("overwrite");
486 
487  #include "createNamedPolyMesh.H"
488 
489  const word oldInstance = mesh.pointsInstance();
490 
491  const dictionary dict(systemDict("createPatchDict", args, mesh));
492 
493  // Whether to synchronise points
494  const Switch pointSync(dict.lookup("pointSync"));
495 
496  // Whether to write cyclic matches to .OBJ files
497  const Switch writeCyclicMatch
498  (
499  dict.lookupOrDefault("writeCyclicMatch", false)
500  );
501 
502 
503  const polyBoundaryMesh& patches = mesh.boundaryMesh();
504 
505  // If running parallel check same patches everywhere
506  patches.checkParallelSync(true);
507 
508 
509  if (writeCyclicMatch)
510  {
511  writeCyclicMatchObjs("initial_", mesh);
512  }
513 
514  // Read patch construct info from dictionary
515  PtrList<dictionary> patchSources(dict.lookup("patches"));
516 
517  HashSet<word> addedPatchNames;
518  forAll(patchSources, addedI)
519  {
520  const dictionary& dict = patchSources[addedI];
521  addedPatchNames.insert(dict.lookup("name"));
522  }
523 
524 
525  // 1. Add all new patches
526  // ~~~~~~~~~~~~~~~~~~~~~~
527 
528  if (patchSources.size())
529  {
530  // Old and new patches.
531  DynamicList<polyPatch*> allPatches(patches.size()+patchSources.size());
532 
533  label startFacei = mesh.nInternalFaces();
534 
535  // Copy old patches.
536  forAll(patches, patchi)
537  {
538  const polyPatch& pp = patches[patchi];
539 
540  if (!isA<processorPolyPatch>(pp))
541  {
542  allPatches.append
543  (
544  pp.clone
545  (
546  patches,
547  patchi,
548  pp.size(),
549  startFacei
550  ).ptr()
551  );
552  startFacei += pp.size();
553  }
554  }
555 
556  forAll(patchSources, addedI)
557  {
558  const dictionary& dict = patchSources[addedI];
559 
560  word patchName(dict.lookup("name"));
561 
562  label destPatchi = patches.findPatchID(patchName);
563 
564  if (destPatchi == -1)
565  {
566  dictionary patchDict(dict.subDict("patchInfo"));
567 
568  destPatchi = allPatches.size();
569 
570  Info<< "Adding new patch " << patchName
571  << " as patch " << destPatchi
572  << " from " << patchDict << endl;
573 
574  patchDict.set("nFaces", 0);
575  patchDict.set("startFace", startFacei);
576 
577  // Add an empty patch.
578  allPatches.append
579  (
581  (
582  patchName,
583  patchDict,
584  destPatchi,
585  patches
586  ).ptr()
587  );
588  }
589  else
590  {
591  Info<< "Patch '" << patchName << "' already exists. Only "
592  << "moving patch faces - type will remain the same" << endl;
593  }
594  }
595 
596  // Copy old patches.
597  forAll(patches, patchi)
598  {
599  const polyPatch& pp = patches[patchi];
600 
601  if (isA<processorPolyPatch>(pp))
602  {
603  allPatches.append
604  (
605  pp.clone
606  (
607  patches,
608  patchi,
609  pp.size(),
610  startFacei
611  ).ptr()
612  );
613  startFacei += pp.size();
614  }
615  }
616 
617  allPatches.shrink();
618  mesh.removeBoundary();
619  mesh.addPatches(allPatches);
620 
621  Info<< endl;
622  }
623 
624 
625 
626  // 2. Repatch faces
627  // ~~~~~~~~~~~~~~~~
628 
629  polyTopoChange meshMod(mesh);
630 
631 
632  forAll(patchSources, addedI)
633  {
634  const dictionary& dict = patchSources[addedI];
635 
636  const word patchName(dict.lookup("name"));
637  label destPatchi = patches.findPatchID(patchName);
638 
639  if (destPatchi == -1)
640  {
642  << "patch " << patchName << " not added. Problem."
643  << abort(FatalError);
644  }
645 
646  const word sourceType(dict.lookup("constructFrom"));
647 
648  if (sourceType == "patches")
649  {
650  labelHashSet patchSources
651  (
652  patches.patchSet
653  (
654  wordReList(dict.lookup("patches"))
655  )
656  );
657 
658  // Repatch faces of the patches.
659  forAllConstIter(labelHashSet, patchSources, iter)
660  {
661  const polyPatch& pp = patches[iter.key()];
662 
663  Info<< "Moving faces from patch " << pp.name()
664  << " to patch " << destPatchi << endl;
665 
666  forAll(pp, i)
667  {
668  changePatchID
669  (
670  mesh,
671  pp.start() + i,
672  destPatchi,
673  meshMod
674  );
675  }
676  }
677  }
678  else if (sourceType == "set")
679  {
680  const word setName(dict.lookup("set"));
681 
682  faceSet faces(mesh, setName);
683 
684  Info<< "Read " << returnReduce(faces.size(), sumOp<label>())
685  << " faces from faceSet " << faces.name() << endl;
686 
687  // Sort (since faceSet contains faces in arbitrary order)
688  labelList faceLabels(faces.toc());
689 
690  SortableList<label> patchFaces(faceLabels);
691 
692  forAll(patchFaces, i)
693  {
694  label facei = patchFaces[i];
695 
696  if (mesh.isInternalFace(facei))
697  {
699  << "Face " << facei << " specified in set "
700  << faces.name()
701  << " is not an external face of the mesh." << endl
702  << "This application can only repatch existing boundary"
703  << " faces." << exit(FatalError);
704  }
705 
706  changePatchID
707  (
708  mesh,
709  facei,
710  destPatchi,
711  meshMod
712  );
713  }
714  }
715  else
716  {
718  << "Invalid source type " << sourceType << endl
719  << "Valid source types are 'patches' 'set'" << exit(FatalError);
720  }
721  }
722  Info<< endl;
723 
724 
725  // Change mesh, use inflation to reforce calculation of transformation
726  // tensors.
727  Info<< "Doing topology modification to order faces." << nl << endl;
728  autoPtr<mapPolyMesh> map = meshMod.changeMesh(mesh, true);
729  mesh.movePoints(map().preMotionPoints());
730 
731  if (writeCyclicMatch)
732  {
733  writeCyclicMatchObjs("coupled_", mesh);
734  }
735 
736  // Synchronise points.
737  if (!pointSync)
738  {
739  Info<< "Not synchronising points." << nl << endl;
740  }
741  else
742  {
743  Info<< "Synchronising points." << nl << endl;
744 
745  pointField newPoints(mesh.points());
746 
747  syncPoints
748  (
749  mesh,
750  newPoints,
752  point(great, great, great)
753  );
754 
755  scalarField diff(mag(newPoints-mesh.points()));
756  Info<< "Points changed by average:" << gAverage(diff)
757  << " max:" << gMax(diff) << nl << endl;
758 
759  mesh.movePoints(newPoints);
760  }
761 
762  // 3. Remove zeros-sized patches
763  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
764 
765  Info<< "Removing patches with no faces in them." << nl<< endl;
766  filterPatches(mesh, addedPatchNames);
767 
768 
769  if (writeCyclicMatch)
770  {
771  writeCyclicMatchObjs("final_", mesh);
772  }
773 
774 
775  // Set the precision of the points data to 10
777 
778  if (!overwrite)
779  {
780  runTime++;
781  }
782  else
783  {
784  mesh.setInstance(oldInstance);
785  }
786 
787  // Write resulting mesh
788  Info<< "Writing repatched mesh to " << runTime.timeName() << nl << endl;
789  mesh.write();
790 
791  Info<< "End\n" << endl;
792 
793  return 0;
794 }
795 
796 
797 // ************************************************************************* //
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:434
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
virtual tmp< scalarField > movePoints(const pointField &)
Move points, returns volumes swept by faces in motion.
Definition: polyMesh.C:1287
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
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: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.
void off()
Switch the function objects off.
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:323
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:248
engineTime & runTime
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:309
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:306
virtual const transformer & transform() const
Return transformation between the coupled patches.
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:636
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:1131
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:982
Input inter-processor communications stream.
Definition: IPstream.H:50
virtual autoPtr< polyPatch > clone(const polyBoundaryMesh &bm) const
Construct and return a clone, resetting the boundary mesh.
Definition: polyPatch.H:228
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:802
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:1169
const Field< PointType > & points() const
Return reference to global points.
const globalMeshData & globalData() const
Return parallel info.
Definition: polyMesh.C:1394
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1156
int neighbProcNo() const
Return neighbour 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:909
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:32
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,.
const functionObjectList & functionObjects() const
Return the list of function objects.
Definition: Time.H:410
#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:476
Direct mesh changes based on v1.3 polyTopoChange syntax.
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:309
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: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.
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
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:844
label setAction(const topoAction &action)
For compatibility with polyTopoChange: set topological action.