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