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-2024 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 "wordReList.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 {
79  const polyBoundaryMesh& patches = mesh.boundaryMesh();
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();
158  mesh.removeBoundary();
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 {
175  const polyBoundaryMesh& patches = mesh.boundaryMesh();
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 
256  const polyBoundaryMesh& patches = mesh.boundaryMesh();
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 "addOverwriteOption.H"
450  #include "addRegionOption.H"
451  #include "addDictOption.H"
452  #include "setRootCase.H"
454 
455  Foam::word meshRegionName = polyMesh::defaultRegion;
456  args.optionReadIfPresent("region", meshRegionName);
457 
458  const bool overwrite = args.optionFound("overwrite");
459 
460  #include "createNamedPolyMesh.H"
461 
462  const word oldInstance = mesh.pointsInstance();
463 
464  const dictionary dict(systemDict("createPatchDict", args, mesh));
465 
466  // Whether to write cyclic matches to .OBJ files
467  const Switch writeCyclicMatch
468  (
469  dict.lookupOrDefault("writeCyclicMatch", false)
470  );
471 
472  const polyBoundaryMesh& patches = mesh.boundaryMesh();
473 
474  // If running parallel check same patches everywhere
475  patches.checkParallelSync(true);
476 
477 
478  if (writeCyclicMatch)
479  {
480  writeCyclicMatchObjs("initial_", mesh);
481  }
482 
483  // Read patch construct info from dictionary
484  PtrList<dictionary> patchSources(dict.lookup("patches"));
485 
486  HashSet<word> addedPatchNames;
487  forAll(patchSources, addedI)
488  {
489  const dictionary& dict = patchSources[addedI];
490  addedPatchNames.insert(dict.lookup<word>("name"));
491  }
492 
493 
494  // 1. Add all new patches
495  // ~~~~~~~~~~~~~~~~~~~~~~
496 
497  if (patchSources.size())
498  {
499  // Old and new patches.
500  DynamicList<polyPatch*> allPatches(patches.size()+patchSources.size());
501 
502  label startFacei = mesh.nInternalFaces();
503 
504  // Copy old patches.
506  {
507  const polyPatch& pp = patches[patchi];
508 
509  if (!isA<processorPolyPatch>(pp))
510  {
511  allPatches.append
512  (
513  pp.clone
514  (
515  patches,
516  patchi,
517  pp.size(),
518  startFacei
519  ).ptr()
520  );
521  startFacei += pp.size();
522  }
523  }
524 
525  forAll(patchSources, addedI)
526  {
527  const dictionary& dict = patchSources[addedI];
528 
529  word patchName(dict.lookup("name"));
530 
531  label destPatchi = patches.findIndex(patchName);
532 
533  if (destPatchi == -1)
534  {
535  dictionary patchDict(dict.subDict("patchInfo"));
536 
537  destPatchi = allPatches.size();
538 
539  Info<< "Adding new patch " << patchName
540  << " as patch " << destPatchi
541  << " from " << patchDict << endl;
542 
543  patchDict.set("nFaces", 0);
544  patchDict.set("startFace", startFacei);
545 
546  // Add an empty patch.
547  allPatches.append
548  (
550  (
551  patchName,
552  patchDict,
553  destPatchi,
554  patches
555  ).ptr()
556  );
557  }
558  else
559  {
560  Info<< "Patch '" << patchName << "' already exists. Only "
561  << "moving patch faces - type will remain the same" << endl;
562  }
563  }
564 
565  // Copy old patches.
567  {
568  const polyPatch& pp = patches[patchi];
569 
570  if (isA<processorPolyPatch>(pp))
571  {
572  allPatches.append
573  (
574  pp.clone
575  (
576  patches,
577  patchi,
578  pp.size(),
579  startFacei
580  ).ptr()
581  );
582  startFacei += pp.size();
583  }
584  }
585 
586  allPatches.shrink();
587  mesh.removeBoundary();
588  mesh.addPatches(allPatches);
589 
590  Info<< endl;
591  }
592 
593 
594 
595  // 2. Repatch faces
596  // ~~~~~~~~~~~~~~~~
597 
598  polyTopoChange meshMod(mesh);
599 
600 
601  forAll(patchSources, addedI)
602  {
603  const dictionary& dict = patchSources[addedI];
604 
605  const word patchName(dict.lookup("name"));
606  label destPatchi = patches.findIndex(patchName);
607 
608  if (destPatchi == -1)
609  {
611  << "patch " << patchName << " not added. Problem."
612  << abort(FatalError);
613  }
614 
615  const word sourceType(dict.lookup("constructFrom"));
616 
617  if (sourceType == "patches")
618  {
619  labelHashSet patchSources
620  (
621  patches.patchSet
622  (
623  wordReList(dict.lookup("patches"))
624  )
625  );
626 
627  // Repatch faces of the patches.
628  forAllConstIter(labelHashSet, patchSources, iter)
629  {
630  const polyPatch& pp = patches[iter.key()];
631 
632  Info<< "Moving faces from patch " << pp.name()
633  << " to patch " << destPatchi << endl;
634 
635  forAll(pp, i)
636  {
637  changePatchID
638  (
639  mesh,
640  pp.start() + i,
641  destPatchi,
642  meshMod
643  );
644  }
645  }
646  }
647  else if (sourceType == "set")
648  {
649  const word setName(dict.lookup("set"));
650 
651  faceSet faces(mesh, setName);
652 
653  Info<< "Read " << returnReduce(faces.size(), sumOp<label>())
654  << " faces from faceSet " << faces.name() << endl;
655 
656  // Sort (since faceSet contains faces in arbitrary order)
657  labelList faceLabels(faces.toc());
658 
659  SortableList<label> patchFaces(faceLabels);
660 
661  forAll(patchFaces, i)
662  {
663  label facei = patchFaces[i];
664 
665  if (mesh.isInternalFace(facei))
666  {
668  << "Face " << facei << " specified in set "
669  << faces.name()
670  << " is not an external face of the mesh." << endl
671  << "This application can only repatch existing boundary"
672  << " faces." << exit(FatalError);
673  }
674 
675  changePatchID
676  (
677  mesh,
678  facei,
679  destPatchi,
680  meshMod
681  );
682  }
683  }
684  else
685  {
687  << "Invalid source type " << sourceType << endl
688  << "Valid source types are 'patches' 'set'" << exit(FatalError);
689  }
690  }
691  Info<< endl;
692 
693 
694  // Change mesh, use motion to reforce calculation of transformation
695  // tensors.
696  Info<< "Doing topology modification to order faces." << nl << endl;
697  autoPtr<polyTopoChangeMap> map = meshMod.changeMesh(mesh);
698 
699  // 3. Remove zeros-sized patches
700  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
701 
702  Info<< "Removing patches with no faces in them." << nl<< endl;
703  filterPatches(mesh, addedPatchNames);
704 
705 
706  if (writeCyclicMatch)
707  {
708  writeCyclicMatchObjs("final_", mesh);
709  }
710 
711 
712  // Set the precision of the points data to 10
714 
715  if (!overwrite)
716  {
717  runTime++;
718  }
719  else
720  {
721  mesh.setInstance(oldInstance);
722  }
723 
724  // Write resulting mesh
725  Info<< "Writing repatched mesh to " << runTime.name() << nl << endl;
726  mesh.write();
727 
728  Info<< "End\n" << endl;
729 
730  return 0;
731 }
732 
733 
734 // ************************************************************************* //
#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: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: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
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
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 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:162
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T, if not found return the given default.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:710
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:843
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
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:267
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1326
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1339
const globalMeshData & globalData() const
Return parallel info.
Definition: polyMesh.C:1515
const fileName & pointsInstance() const
Return the current instance directory for points.
Definition: polyMesh.C:965
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:404
void addPatches(const List< polyPatch * > &, const bool validBoundary=true)
Add boundary patches.
Definition: polyMesh.C:1084
void removeBoundary()
Remove boundary patches.
Definition: polyMeshClear.C:36
void setInstance(const fileName &)
Set the instance for mesh files.
Definition: polyMeshIO.C:78
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.
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:334
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: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
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:257
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
static const char nl
Definition: Ostream.H:266
dictionary dict
Foam::argList args(argc, argv)