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-2025 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 "zoneGenerator.H"
48 #include "faceSet.H"
49 #include "polyTopoChange.H"
50 #include "systemDict.H"
51 
52 using namespace Foam;
53 
54 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
55 
56 void changePatchID
57 (
58  const polyMesh& mesh,
59  const label faceID,
60  const label patchID,
61  polyTopoChange& meshMod
62 )
63 {
64  meshMod.modifyFace
65  (
66  mesh.faces()[faceID], // face
67  faceID, // face ID
68  mesh.faceOwner()[faceID], // owner
69  -1, // neighbour
70  false, // flip flux
71  patchID // patch ID
72  );
73 }
74 
75 
76 // Filter out the empty patches.
77 void filterPatches(polyMesh& mesh, const HashSet<word>& addedPatchNames)
78 {
80 
81  // Patches to keep
82  DynamicList<polyPatch*> allPatches(patches.size());
83 
84  label nOldPatches = returnReduce(patches.size(), sumOp<label>());
85 
86  // Copy old patches.
88  {
89  const polyPatch& pp = patches[patchi];
90 
91  // Note: reduce possible since non-proc patches guaranteed in same order
92  if (!isA<processorPolyPatch>(pp))
93  {
94  // Add if
95  // - listed in createPatchDict
96  // - or constraint type
97  // - or non zero size
98  if
99  (
100  addedPatchNames.found(pp.name())
101  || polyPatch::constraintType(pp.type())
102  || returnReduce(pp.size(), sumOp<label>()) > 0
103  )
104  {
105  allPatches.append
106  (
107  pp.clone
108  (
109  patches,
110  allPatches.size(),
111  pp.size(),
112  pp.start()
113  ).ptr()
114  );
115  }
116  else
117  {
118  Info<< "Removing zero-sized patch " << pp.name()
119  << " type " << pp.type()
120  << " at position " << patchi << endl;
121  }
122  }
123  }
124 
125  // Copy non-empty processor patches
127  {
128  const polyPatch& pp = patches[patchi];
129 
130  if (isA<processorPolyPatch>(pp))
131  {
132  if (pp.size())
133  {
134  allPatches.append
135  (
136  pp.clone
137  (
138  patches,
139  allPatches.size(),
140  pp.size(),
141  pp.start()
142  ).ptr()
143  );
144  }
145  else
146  {
147  Info<< "Removing empty processor patch " << pp.name()
148  << " at position " << patchi << endl;
149  }
150  }
151  }
152 
153  label nAllPatches = returnReduce(allPatches.size(), sumOp<label>());
154  if (nAllPatches != nOldPatches)
155  {
156  Info<< "Removing patches." << endl;
157  allPatches.shrink();
159  mesh.addPatches(allPatches);
160  }
161  else
162  {
163  Info<< "No patches removed." << endl;
164  forAll(allPatches, i)
165  {
166  delete allPatches[i];
167  }
168  }
169 }
170 
171 
172 // Write current match for all patches the as OBJ files
173 void writeCyclicMatchObjs(const fileName& prefix, const polyMesh& mesh)
174 {
176 
178  {
179  if
180  (
181  isA<cyclicPolyPatch>(patches[patchi])
182  && refCast<const cyclicPolyPatch>(patches[patchi]).owner()
183  )
184  {
185  const cyclicPolyPatch& cycPatch =
186  refCast<const cyclicPolyPatch>(patches[patchi]);
187 
188  // Write patches
189  {
190  OFstream str(prefix+cycPatch.name() + ".obj");
191  Pout<< "Writing " << cycPatch.name()
192  << " faces to " << str.name() << endl;
194  (
195  str,
196  cycPatch,
197  cycPatch.points()
198  );
199  }
200 
201  const cyclicPolyPatch& nbrPatch = cycPatch.nbrPatch();
202  {
203  OFstream str(prefix+nbrPatch.name()+".obj");
204  Pout<< "Writing " << nbrPatch.name()
205  << " faces to " << str.name() << endl;
207  (
208  str,
209  nbrPatch,
210  nbrPatch.points()
211  );
212  }
213 
214 
215  // Lines between corresponding face centres
216  OFstream str(prefix+cycPatch.name()+nbrPatch.name()+"_match.obj");
217  label vertI = 0;
218 
219  Pout<< "Writing cyclic match as lines between face centres to "
220  << str.name() << endl;
221 
222  forAll(cycPatch, facei)
223  {
224  const point& fc0 = mesh.faceCentres()[cycPatch.start()+facei];
225  meshTools::writeOBJ(str, fc0);
226  vertI++;
227  const point& fc1 = mesh.faceCentres()[nbrPatch.start()+facei];
228  meshTools::writeOBJ(str, fc1);
229  vertI++;
230 
231  str<< "l " << vertI-1 << ' ' << vertI << nl;
232  }
233  }
234  }
235 }
236 
237 
238 // Synchronise points on both sides of coupled boundaries.
239 template<class CombineOp>
240 void syncPoints
241 (
242  const polyMesh& mesh,
244  const CombineOp& cop,
245  const point& nullValue
246 )
247 {
248  if (points.size() != mesh.nPoints())
249  {
251  << "Number of values " << points.size()
252  << " is not equal to the number of points in the mesh "
253  << mesh.nPoints() << abort(FatalError);
254  }
255 
257 
258  // Is there any coupled patch with transformation?
259  bool hasTransformation = false;
260 
261  if (Pstream::parRun())
262  {
263  // Send
264 
266  {
267  const polyPatch& pp = patches[patchi];
268 
269  if
270  (
271  isA<processorPolyPatch>(pp)
272  && pp.nPoints() > 0
273  && refCast<const processorPolyPatch>(pp).owner()
274  )
275  {
276  const processorPolyPatch& procPatch =
277  refCast<const processorPolyPatch>(pp);
278 
279  // Get data per patchPoint in neighbouring point numbers.
280  pointField patchInfo(procPatch.nPoints(), nullValue);
281 
282  const labelList& meshPts = procPatch.meshPoints();
283  const labelList& nbrPts = procPatch.nbrPoints();
284 
285  forAll(nbrPts, pointi)
286  {
287  label nbrPointi = nbrPts[pointi];
288  if (nbrPointi >= 0 && nbrPointi < patchInfo.size())
289  {
290  patchInfo[nbrPointi] = points[meshPts[pointi]];
291  }
292  }
293 
294  OPstream toNbr
295  (
297  procPatch.neighbProcNo()
298  );
299  toNbr << patchInfo;
300  }
301  }
302 
303 
304  // Receive and set.
305 
307  {
308  const polyPatch& pp = patches[patchi];
309 
310  if
311  (
312  isA<processorPolyPatch>(pp)
313  && pp.nPoints() > 0
314  && !refCast<const processorPolyPatch>(pp).owner()
315  )
316  {
317  const processorPolyPatch& procPatch =
318  refCast<const processorPolyPatch>(pp);
319 
320  pointField nbrPatchInfo(procPatch.nPoints());
321  {
322  // We do not know the number of points on the other side
323  // so cannot use Pstream::read.
324  IPstream fromNbr
325  (
327  procPatch.neighbProcNo()
328  );
329  fromNbr >> nbrPatchInfo;
330  }
331  // Null any value which is not on neighbouring processor
332  nbrPatchInfo.setSize(procPatch.nPoints(), nullValue);
333 
334  if (procPatch.transform().transformsPosition())
335  {
336  hasTransformation = true;
337  procPatch.transform().transformPosition
338  (
339  nbrPatchInfo,
340  nbrPatchInfo
341  );
342  }
343 
344  const labelList& meshPts = procPatch.meshPoints();
345 
346  forAll(meshPts, pointi)
347  {
348  label meshPointi = meshPts[pointi];
349  points[meshPointi] = nbrPatchInfo[pointi];
350  }
351  }
352  }
353  }
354 
355  // Do the cyclics.
357  {
358  const polyPatch& pp = patches[patchi];
359 
360  if
361  (
362  isA<cyclicPolyPatch>(pp)
363  && refCast<const cyclicPolyPatch>(pp).owner()
364  )
365  {
366  const cyclicPolyPatch& cycPatch =
367  refCast<const cyclicPolyPatch>(pp);
368 
369  const edgeList& coupledPoints = cycPatch.coupledPoints();
370  const labelList& meshPts = cycPatch.meshPoints();
371  const cyclicPolyPatch& nbrPatch = cycPatch.nbrPatch();
372  const labelList& nbrMeshPts = nbrPatch.meshPoints();
373 
374  pointField patchPoints(coupledPoints.size());
375 
376  forAll(coupledPoints, i)
377  {
378  const edge& e = coupledPoints[i];
379  label point0 = meshPts[e[0]];
380  patchPoints[i] = points[point0];
381  }
382 
383  if (cycPatch.transform().transformsPosition())
384  {
385  hasTransformation = true;
387  (
388  patchPoints,
389  patchPoints
390  );
391  }
392 
393  forAll(coupledPoints, i)
394  {
395  const edge& e = coupledPoints[i];
396  label point1 = nbrMeshPts[e[1]];
397  points[point1] = patchPoints[i];
398  }
399  }
400  }
401 
402  //- Note: hasTransformation is only used for warning messages so
403  // reduction not strictly necessary.
404  // reduce(hasTransformation, orOp<bool>());
405 
406  // Synchronise multiple shared points.
407  const globalMeshData& pd = mesh.globalData();
408 
409  if (pd.nGlobalPoints() > 0)
410  {
411  if (hasTransformation)
412  {
414  << "There are decomposed cyclics in this mesh with"
415  << " transformations." << endl
416  << "This is not supported. The result will be incorrect"
417  << endl;
418  }
419 
420 
421  // Values on shared points.
422  pointField sharedPts(pd.nGlobalPoints(), nullValue);
423 
424  forAll(pd.sharedPointLabels(), i)
425  {
426  label meshPointi = pd.sharedPointLabels()[i];
427  // Fill my entries in the shared points
428  sharedPts[pd.sharedPointAddr()[i]] = points[meshPointi];
429  }
430 
431  // Combine on master.
432  Pstream::listCombineGather(sharedPts, cop);
433  Pstream::listCombineScatter(sharedPts);
434 
435  // Now we will all have the same information. Merge it back with
436  // my local information.
437  forAll(pd.sharedPointLabels(), i)
438  {
439  label meshPointi = pd.sharedPointLabels()[i];
440  points[meshPointi] = sharedPts[pd.sharedPointAddr()[i]];
441  }
442  }
443 }
444 
445 
446 
447 int main(int argc, char *argv[])
448 {
449  #include "addNoOverwriteOption.H"
450  #include "addMeshOption.H"
451  #include "addRegionOption.H"
452  #include "addDictOption.H"
453  #include "setRootCase.H"
455 
456  Foam::word meshRegionName = polyMesh::defaultRegion;
457  args.optionReadIfPresent("region", meshRegionName);
458 
459  #include "setNoOverwrite.H"
460 
461  #include "createSpecifiedPolyMesh.H"
462 
463  const word oldInstance = mesh.pointsInstance();
464 
465  const dictionary dict(systemDict("createPatchDict", args, mesh));
466 
467  // Whether to write cyclic matches to .OBJ files
468  const Switch writeCyclicMatch
469  (
470  dict.lookupOrDefault("writeCyclicMatch", false)
471  );
472 
474 
475  // If running parallel check same patches everywhere
476  patches.checkParallelSync(true);
477 
478 
479  if (writeCyclicMatch)
480  {
481  writeCyclicMatchObjs("initial_", mesh);
482  }
483 
484  // For backwards-compatibility read patches as list of dictionaries
485  // if not a dictionary
486  dictionary patchesListDict("patches", dict);
487 
488  if (!dict.isDict("patches"))
489  {
490  // Read patch construct info from dictionary
491  PtrList<dictionary> patchSources(dict.lookup("patches"));
492 
493  forAll(patchSources, psi)
494  {
495  const dictionary& dict = patchSources[psi];
496  patchesListDict.add(dict.lookup<word>("name"), dict);
497  }
498  }
499 
500  const dictionary& patchesDict =
501  dict.isDict("patches")
502  ? dict.subDict("patches")
503  : patchesListDict;
504 
505 
506  HashSet<word> addedPatchNames;
507  forAllConstIter(dictionary, patchesDict, iter)
508  {
509  addedPatchNames.insert(iter().keyword());
510  }
511 
512 
513  // 1. Add all new patches
514  // ~~~~~~~~~~~~~~~~~~~~~~
515 
516  if (patchesDict.size())
517  {
518  // Old and new patches.
519  DynamicList<polyPatch*> allPatches(patches.size()+patchesDict.size());
520 
521  label startFacei = mesh.nInternalFaces();
522 
523  // Copy old patches.
525  {
526  const polyPatch& pp = patches[patchi];
527 
528  if (!isA<processorPolyPatch>(pp))
529  {
530  allPatches.append
531  (
532  pp.clone
533  (
534  patches,
535  patchi,
536  pp.size(),
537  pp.start()
538  ).ptr()
539  );
540  startFacei += pp.size();
541  }
542  }
543 
544  forAllConstIter(dictionary, patchesDict, iter)
545  {
546  const word& patchName = iter().keyword();
547  const dictionary& dict = iter().dict();
548 
549  label destPatchi = patches.findIndex(patchName);
550 
551  if (destPatchi == -1)
552  {
553  dictionary patchDict(dict.subDict("patchInfo"));
554 
555  destPatchi = allPatches.size();
556 
557  Info<< "Adding new patch " << patchName
558  << " as patch " << destPatchi
559  << " from " << patchDict << endl;
560 
561  patchDict.set("nFaces", 0);
562  patchDict.set("startFace", startFacei);
563 
564  // Add an empty patch.
565  allPatches.append
566  (
568  (
569  patchName,
570  patchDict,
571  destPatchi,
572  patches
573  ).ptr()
574  );
575  }
576  else
577  {
578  Info<< "Patch '" << patchName << "' already exists. Only "
579  << "moving patch faces - type will remain the same" << endl;
580  }
581  }
582 
583  // Copy old patches.
585  {
586  const polyPatch& pp = patches[patchi];
587 
588  if (isA<processorPolyPatch>(pp))
589  {
590  allPatches.append
591  (
592  pp.clone
593  (
594  patches,
595  patchi,
596  pp.size(),
597  pp.start()
598  ).ptr()
599  );
600  startFacei += pp.size();
601  }
602  }
603 
604  allPatches.shrink();
606  mesh.addPatches(allPatches);
607 
608  Info<< endl;
609  }
610 
611 
612 
613  // 2. Repatch faces
614  // ~~~~~~~~~~~~~~~~
615 
616  polyTopoChange meshMod(mesh);
617 
618  forAllConstIter(dictionary, patchesDict, iter)
619  {
620  const word& patchName = iter().keyword();
621  const dictionary& dict = iter().dict();
622 
623  label destPatchi = patches.findIndex(patchName);
624 
625  if (destPatchi == -1)
626  {
628  << "patch " << patchName << " not added. Problem."
629  << abort(FatalError);
630  }
631 
632  const word sourceType(dict.lookup("constructFrom"));
633 
634  if (sourceType == "patches")
635  {
636  labelHashSet patchesDict
637  (
638  patches.patchSet
639  (
640  wordReList(dict.lookup("patches"))
641  )
642  );
643 
644  // Repatch faces of the patches.
645  forAllConstIter(labelHashSet, patchesDict, iter)
646  {
647  const polyPatch& pp = patches[iter.key()];
648 
649  Info<< "Moving faces from patch " << pp.name()
650  << " to patch " << destPatchi << endl;
651 
652  forAll(pp, i)
653  {
654  changePatchID
655  (
656  mesh,
657  pp.start() + i,
658  destPatchi,
659  meshMod
660  );
661  }
662  }
663  }
664  else if (sourceType == "zone")
665  {
666  SortableList<label> patchFaces;
667 
668  if (dict.isDict("zone"))
669  {
671  (
673  (
674  "zone",
676  mesh,
677  dict.subDict("zone")
678  )
679  );
680 
681  patchFaces = zg->generate().fZone();
682 
683  Info<< "Set "
684  << returnReduce(patchFaces.size(), sumOp<label>())
685  << " faces from zoneGenerator " << zg->name()
686  << " of type " << zg->type() << endl;
687  }
688  else
689  {
690  const word zoneName(dict.lookup("zone"));
691 
692  patchFaces = mesh.faceZones()[zoneName];
693 
694  Info<< "Read "
695  << returnReduce(patchFaces.size(), sumOp<label>())
696  << " faces from faceZone " << zoneName << endl;
697  }
698 
699  patchFaces.sort();
700 
701  forAll(patchFaces, i)
702  {
703  label facei = patchFaces[i];
704 
705  if (mesh.isInternalFace(facei))
706  {
708  << "Face " << facei << " specified in faceZone "
709  << " is not an external face of the mesh." << endl
710  << "This application can only repatch existing boundary"
711  << " faces." << exit(FatalError);
712  }
713 
714  changePatchID
715  (
716  mesh,
717  facei,
718  destPatchi,
719  meshMod
720  );
721  }
722  }
723  else if (sourceType == "set")
724  {
725  const word setName(dict.lookup("set"));
726 
727  faceSet faces(mesh, setName);
728 
729  Info<< "Read " << returnReduce(faces.size(), sumOp<label>())
730  << " faces from faceSet " << faces.name() << endl;
731 
732  // Sort (since faceSet contains faces in arbitrary order)
733  labelList faceLabels(faces.toc());
734 
735  SortableList<label> patchFaces(faceLabels);
736 
737  forAll(patchFaces, i)
738  {
739  label facei = patchFaces[i];
740 
741  if (mesh.isInternalFace(facei))
742  {
744  << "Face " << facei << " specified in set "
745  << faces.name()
746  << " is not an external face of the mesh." << endl
747  << "This application can only repatch existing boundary"
748  << " faces." << exit(FatalError);
749  }
750 
751  changePatchID
752  (
753  mesh,
754  facei,
755  destPatchi,
756  meshMod
757  );
758  }
759  }
760  else
761  {
763  << "Invalid source type " << sourceType << endl
764  << "Valid source types are 'patches' 'set'" << exit(FatalError);
765  }
766  }
767  Info<< endl;
768 
769 
770  // Change mesh, use motion to reforce calculation of transformation
771  // tensors.
772  Info<< "Doing topology modification to order faces." << nl << endl;
774 
775  // 3. Remove zeros-sized patches
776  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
777 
778  Info<< "Removing patches with no faces in them." << nl<< endl;
779  filterPatches(mesh, addedPatchNames);
780 
781 
782  if (writeCyclicMatch)
783  {
784  writeCyclicMatchObjs("final_", mesh);
785  }
786 
787 
788  // Set the precision of the points data to 10
790 
791  if (!overwrite)
792  {
793  runTime++;
794  }
795  else
796  {
797  mesh.setInstance(oldInstance);
798  }
799 
800  // Write resulting mesh
801  Info<< "Writing repatched mesh to " << runTime.name() << nl << endl;
802  mesh.write();
803 
804  Info<< "End\n" << endl;
805 
806  return 0;
807 }
808 
809 
810 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:476
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:109
bool found(const Key &) const
Return true if hashedEntry is found in table.
Definition: HashTable.C:138
static unsigned int defaultPrecision()
Return the default precision.
Definition: IOstream.H:473
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
Output to file stream.
Definition: OFstream.H:87
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
void sort()
(stable) sort the list (if changed after construction time)
Definition: SortableList.C:112
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition: Switch.H:61
label size() const
Return the number of elements in the UList.
Definition: UListI.H:311
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 optionReadIfPresent(const word &opt, T &) const
Read a value from the named option if present.
Definition: argListI.H:255
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 keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:740
bool isDict(const word &) const
Check if entry is a sub-dictionary.
Definition: dictionary.C:803
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:849
bool add(entry *, bool mergeEntry=false)
Add a new entry.
Definition: dictionary.C:1020
T lookupOrDefault(const word &, const T &, const bool writeDefault=writeOptionalEntries > 0) const
Find and return a T, if not found return the given default.
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 class for handling file names.
Definition: fileName.H:82
virtual bool write(const bool write=true) const
Write mesh using IO settings from time.
Definition: fvMesh.C:1785
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:270
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1344
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1357
const globalMeshData & globalData() const
Return parallel info.
Definition: polyMesh.C:1521
const fileName & pointsInstance() const
Return the current instance directory for points.
Definition: polyMesh.C:988
const faceZoneList & faceZones() const
Return face zones.
Definition: polyMesh.H:443
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:401
void addPatches(const List< polyPatch * > &, const bool validBoundary=true)
Add boundary patches.
Definition: polyMesh.C:1107
void removeBoundary()
Remove boundary patches.
Definition: polyMeshClear.C:36
void setInstance(const fileName &)
Set the instance for mesh files.
Definition: polyMeshIO.C:91
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 syncParallel=true, const bool orderCells=false, const bool orderPoints=false)
Inplace changes mesh without change of patches.
void modifyFace(const face &f, const label facei, const label own, const label nei, const bool flipFaceFlux, const label patchID)
Modify vertices or cell of face.
const vectorField & faceCentres() const
label nInternalFaces() const
label nPoints() 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.
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
static autoPtr< zoneGenerator > New(const word &name, const polyMesh &mesh, const dictionary &dict)
Select constructed from name, mesh and dictionary.
Definition: zoneGenerator.C:66
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
int main(int argc, char *argv[])
Definition: financialFoam.C:44
label patchi
const pointField & points
const fvPatchList & patches
const volScalarField & psi
#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:106
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
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:258
errorManip< error > abort(error &err)
Definition: errorManip.H:131
messageStream Info
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
error FatalError
IOdictionary systemDict(const word &dictName, const argList &args, const objectRegistry &ob, const word &regionName=polyMesh::defaultRegion, const fileName &path=fileName::null)
Definition: systemDict.C:93
static const char nl
Definition: Ostream.H:267
dictionary dict
const bool overwrite
Definition: setNoOverwrite.H:1
Foam::argList args(argc, argv)