perfectInterface.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2018 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 Description
28  Best thing is probably to look at attachDetach which does almost exactly
29  the same but for the geometric matching of points and face centres.
30 
31 \*---------------------------------------------------------------------------*/
32 
33 #include "perfectInterface.H"
34 #include "polyTopoChanger.H"
35 #include "polyMesh.H"
36 #include "polyTopoChange.H"
38 #include "mapPolyMesh.H"
39 #include "matchPoints.H"
40 #include "polyModifyFace.H"
41 #include "polyRemovePoint.H"
42 #include "polyRemoveFace.H"
43 #include "indirectPrimitivePatch.H"
44 
45 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
46 
47 namespace Foam
48 {
49  defineTypeNameAndDebug(perfectInterface, 0);
51  (
52  polyMeshModifier,
53  perfectInterface,
54  dictionary
55  );
56 }
57 
58 
59 // Tolerance used as fraction of minimum edge length.
60 const Foam::scalar Foam::perfectInterface::tol_ = 1e-3;
61 
62 
63 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
64 
65 Foam::pointField Foam::perfectInterface::calcFaceCentres
66 (
67  const indirectPrimitivePatch& pp
68 )
69 {
70  const pointField& points = pp.points();
71 
72  pointField ctrs(pp.size());
73 
74  forAll(ctrs, patchFacei)
75  {
76  ctrs[patchFacei] = pp[patchFacei].centre(points);
77  }
78 
79  return ctrs;
80 }
81 
82 
83 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
84 
85 Foam::perfectInterface::perfectInterface
86 (
87  const word& name,
88  const label index,
89  const polyTopoChanger& mme,
90  const word& faceZoneName,
91  const word& masterPatchName,
92  const word& slavePatchName
93 )
94 :
95  polyMeshModifier(name, index, mme, true),
96  faceZoneID_(faceZoneName, mme.mesh().faceZones()),
97  masterPatchID_(masterPatchName, mme.mesh().boundaryMesh()),
98  slavePatchID_(slavePatchName, mme.mesh().boundaryMesh())
99 {}
100 
101 
102 Foam::perfectInterface::perfectInterface
103 (
104  const word& name,
105  const dictionary& dict,
106  const label index,
107  const polyTopoChanger& mme
108 )
109 :
110  polyMeshModifier(name, index, mme, dict.get<bool>("active")),
111  faceZoneID_
112  (
113  dict.lookup("faceZoneName"),
114  mme.mesh().faceZones()
115  ),
116  masterPatchID_
117  (
118  dict.lookup("masterPatchName"),
119  mme.mesh().boundaryMesh()
120  ),
121  slavePatchID_
122  (
123  dict.lookup("slavePatchName"),
124  mme.mesh().boundaryMesh()
125  )
126 {}
127 
128 
129 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
130 
132 {}
133 
134 
135 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
136 
138 {
139  // If modifier is inactive, skip change
140  if (!active())
141  {
142  if (debug)
143  {
144  Pout<< "bool perfectInterface::changeTopology() const "
145  << "for object " << name() << " : "
146  << "Inactive" << endl;
147  }
148 
149  return false;
150  }
151  else
152  {
153  // I want topo change every time step.
154  return true;
155  }
156 }
157 
158 
160 (
161  const indirectPrimitivePatch& pp0,
162  const indirectPrimitivePatch& pp1,
164 ) const
165 {
166  const polyMesh& mesh = topoChanger().mesh();
167 
169 
170  // Some aliases
171  const edgeList& edges0 = pp0.edges();
172  const pointField& pts0 = pp0.localPoints();
173  const pointField& pts1 = pp1.localPoints();
174  const labelList& meshPts0 = pp0.meshPoints();
175  const labelList& meshPts1 = pp1.meshPoints();
176 
177 
178  // Get local dimension as fraction of minimum edge length
179 
180  scalar minLen = GREAT;
181 
182  forAll(edges0, edgeI)
183  {
184  minLen = min(minLen, edges0[edgeI].mag(pts0));
185  }
186  scalar typDim = tol_*minLen;
187 
188  if (debug)
189  {
190  Pout<< "typDim:" << typDim << " edges0:" << edges0.size()
191  << " pts0:" << pts0.size() << " pts1:" << pts1.size()
192  << " pp0:" << pp0.size() << " pp1:" << pp1.size() << endl;
193  }
194 
195 
196  // Determine pointMapping in mesh point labels. Uses geometric
197  // comparison to find correspondence between patch points.
198 
199  labelList renumberPoints(mesh.points().size());
200  forAll(renumberPoints, i)
201  {
202  renumberPoints[i] = i;
203  }
204  {
205  labelList from1To0Points(pts1.size());
206 
207  bool matchOk = matchPoints
208  (
209  pts1,
210  pts0,
211  scalarField(pts1.size(), typDim), // tolerance
212  true, // verbose
213  from1To0Points
214  );
215 
216  if (!matchOk)
217  {
219  << "Points on patch sides do not match to within tolerance "
220  << typDim << exit(FatalError);
221  }
222 
223  forAll(pts1, i)
224  {
225  renumberPoints[meshPts1[i]] = meshPts0[from1To0Points[i]];
226  }
227  }
228 
229 
230 
231  // Calculate correspondence between patch faces
232 
233  labelList from0To1Faces(pp1.size());
234 
235  bool matchOk = matchPoints
236  (
237  calcFaceCentres(pp0),
238  calcFaceCentres(pp1),
239  scalarField(pp0.size(), typDim), // tolerance
240  true, // verbose
241  from0To1Faces
242  );
243 
244  if (!matchOk)
245  {
247  << "Face centres of patch sides do not match to within tolerance "
248  << typDim << exit(FatalError);
249  }
250 
251 
252 
253  // Now
254  // - renumber faces using pts1 (except patch1 faces)
255  // - remove patch1 faces. Remember cell label on owner side.
256  // - modify patch0 faces to be internal.
257 
258  // 1. Get faces to be renumbered
259  labelHashSet affectedFaces(2*pp1.size());
260  forAll(meshPts1, i)
261  {
262  label meshPointi = meshPts1[i];
263 
264  if (meshPointi != renumberPoints[meshPointi])
265  {
266  const labelList& pFaces = mesh.pointFaces()[meshPointi];
267 
268  affectedFaces.insert(pFaces);
269  }
270  }
271  forAll(pp1, i)
272  {
273  affectedFaces.erase(pp1.addressing()[i]);
274  }
275  // Remove patch0 from renumbered faces. Should not be necessary since
276  // patch0 and 1 should not share any point (if created by mergeMeshing)
277  // so affectedFaces should not contain any patch0 faces but you can
278  // never be sure what the user is doing.
279  forAll(pp0, i)
280  {
281  label facei = pp0.addressing()[i];
282 
283  if (affectedFaces.erase(facei))
284  {
286  << "Found face " << facei << " vertices "
287  << mesh.faces()[facei] << " whose points are"
288  << " used both by master patch and slave patch" << endl;
289  }
290  }
291 
292 
293  // 2. Renumber (non patch0/1) faces.
294  for (const label facei : affectedFaces)
295  {
296  const face& f = mesh.faces()[facei];
297 
298  face newFace(f.size());
299 
300  forAll(newFace, fp)
301  {
302  newFace[fp] = renumberPoints[f[fp]];
303  }
304 
305  label nbr = -1;
306 
307  label patchi = -1;
308 
309  if (mesh.isInternalFace(facei))
310  {
311  nbr = mesh.faceNeighbour()[facei];
312  }
313  else
314  {
315  patchi = patches.whichPatch(facei);
316  }
317 
318  label zoneID = mesh.faceZones().whichZone(facei);
319 
320  bool zoneFlip = false;
321 
322  if (zoneID >= 0)
323  {
324  const faceZone& fZone = mesh.faceZones()[zoneID];
325 
326  zoneFlip = fZone.flipMap()[fZone.whichFace(facei)];
327  }
328 
329  ref.setAction
330  (
332  (
333  newFace, // modified face
334  facei, // label of face being modified
335  mesh.faceOwner()[facei], // owner
336  nbr, // neighbour
337  false, // face flip
338  patchi, // patch for face
339  false, // remove from zone
340  zoneID, // zone for face
341  zoneFlip // face flip in zone
342  )
343  );
344  }
345 
346 
347  // 3. Remove patch1 points
348  forAll(meshPts1, i)
349  {
350  label meshPointi = meshPts1[i];
351 
352  if (meshPointi != renumberPoints[meshPointi])
353  {
354  ref.setAction(polyRemovePoint(meshPointi));
355  }
356  }
357 
358 
359  // 4. Remove patch1 faces
360  forAll(pp1, i)
361  {
362  label facei = pp1.addressing()[i];
363  ref.setAction(polyRemoveFace(facei));
364  }
365 
366 
367  // 5. Modify patch0 faces for new points (not really necessary; see
368  // comment above about patch1 and patch0 never sharing points) and
369  // becoming internal.
370  const boolList& mfFlip =
371  mesh.faceZones()[faceZoneID_.index()].flipMap();
372 
373  forAll(pp0, i)
374  {
375  label facei = pp0.addressing()[i];
376 
377  const face& f = mesh.faces()[facei];
378 
379  face newFace(f.size());
380 
381  forAll(newFace, fp)
382  {
383  newFace[fp] = renumberPoints[f[fp]];
384  }
385 
386  label own = mesh.faceOwner()[facei];
387 
388  label pp1Facei = pp1.addressing()[from0To1Faces[i]];
389 
390  label nbr = mesh.faceOwner()[pp1Facei];
391 
392  if (own < nbr)
393  {
394  ref.setAction
395  (
397  (
398  newFace, // modified face
399  facei, // label of face being modified
400  own, // owner
401  nbr, // neighbour
402  false, // face flip
403  -1, // patch for face
404  false, // remove from zone
405  faceZoneID_.index(), // zone for face
406  mfFlip[i] // face flip in zone
407  )
408  );
409  }
410  else
411  {
412  ref.setAction
413  (
415  (
416  newFace.reverseFace(), // modified face
417  facei, // label of face being modified
418  nbr, // owner
419  own, // neighbour
420  true, // face flip
421  -1, // patch for face
422  false, // remove from zone
423  faceZoneID_.index(), // zone for face
424  !mfFlip[i] // face flip in zone
425  )
426  );
427  }
428  }
429 }
430 
431 
433 {
434  if (debug)
435  {
436  Pout<< "bool perfectInterface::setRefinement(polyTopoChange&) const : "
437  << "for object " << name() << " : "
438  << "masterPatchID_:" << masterPatchID_
439  << " slavePatchID_:" << slavePatchID_
440  << " faceZoneID_:" << faceZoneID_ << endl;
441  }
442 
443  if
444  (
445  masterPatchID_.active()
446  && slavePatchID_.active()
447  && faceZoneID_.active()
448  )
449  {
450  const polyMesh& mesh = topoChanger().mesh();
451 
453  const polyPatch& patch0 = patches[masterPatchID_.index()];
454  const polyPatch& patch1 = patches[slavePatchID_.index()];
455 
456 
457  labelList pp0Labels(identity(patch0.size(), patch0.start()));
459  (
460  IndirectList<face>(mesh.faces(), pp0Labels),
461  mesh.points()
462  );
463 
464  labelList pp1Labels(identity(patch1.size(), patch1.start()));
466  (
467  IndirectList<face>(mesh.faces(), pp1Labels),
468  mesh.points()
469  );
470 
471  setRefinement(pp0, pp1, ref);
472  }
473 }
474 
475 
477 {
478  // Update only my points. Nothing to be done here as points already
479  // shared by now.
480 }
481 
482 
484 {
485  // Mesh has changed topologically. Update local topological data
486  const polyMesh& mesh = topoChanger().mesh();
487 
488  faceZoneID_.update(mesh.faceZones());
489  masterPatchID_.update(mesh.boundaryMesh());
490  slavePatchID_.update(mesh.boundaryMesh());
491 }
492 
493 
495 {
496  os << nl << type() << nl
497  << name()<< nl
498  << faceZoneID_.name() << nl
499  << masterPatchID_.name() << nl
500  << slavePatchID_.name() << endl;
501 }
502 
503 
505 {
506  os << nl;
507 
508  os.beginBlock(name());
509  os.writeEntry("type", type());
510  os.writeEntry("active", active());
511  os.writeEntry("faceZoneName", faceZoneID_.name());
512  os.writeEntry("masterPatchName", masterPatchID_.name());
513  os.writeEntry("slavePatchName", slavePatchID_.name());
514  os.endBlock();
515 }
516 
517 
518 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
Foam::perfectInterface::~perfectInterface
virtual ~perfectInterface()
Destructor.
Definition: perfectInterface.C:131
Foam::pointField
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:44
Foam::polyMesh::points
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1038
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
ref
rDeltaT ref()
polyRemovePoint.H
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::polyBoundaryMesh
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
Definition: polyBoundaryMesh.H:62
Foam::PrimitivePatch::edges
const edgeList & edges() const
Return list of edges, address into LOCAL point list.
Definition: PrimitivePatch.C:238
Foam::polyRemovePoint
Class containing data for point removal.
Definition: polyRemovePoint.H:48
Foam::dynamicFvMesh::update
virtual bool update()=0
Update the mesh for both mesh motion and topology change.
Foam::PrimitivePatch::localPoints
const Field< PointType > & localPoints() const
Return pointField of points in patch.
Definition: PrimitivePatch.C:458
polyRemoveFace.H
mapPolyMesh.H
Foam::polyTopoChanger
List of mesh modifiers defining the mesh dynamics.
Definition: polyTopoChanger.H:67
Foam::perfectInterface::modifyMotionPoints
virtual void modifyMotionPoints(pointField &motionPoints) const
Modify motion points to comply with the topological change.
Definition: perfectInterface.C:476
polyTopoChanger.H
polyTopoChange.H
Foam::primitiveMesh::pointFaces
const labelListList & pointFaces() const
Definition: primitiveMeshPointFaces.C:34
Foam::polyTopoChange
Direct mesh changes based on v1.3 polyTopoChange syntax.
Definition: polyTopoChange.H:100
Foam::perfectInterface::updateMesh
virtual void updateMesh(const mapPolyMesh &)
Force recalculation of locally stored data on topological change.
Definition: perfectInterface.C:483
Foam::polyModifyFace
Class describing modification of a face.
Definition: polyModifyFace.H:50
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:435
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:337
Foam::Ostream::beginBlock
virtual Ostream & beginBlock(const keyType &kw)
Write begin block group with the given name.
Definition: Ostream.C:91
Foam::dictionary::get
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:81
Foam::Pout
prefixOSstream Pout
An Ostream wrapper for parallel output to std::cout.
Foam::perfectInterface::setRefinement
virtual void setRefinement(polyTopoChange &) const
Insert the layer addition/removal instructions.
Definition: perfectInterface.C:432
polyMesh.H
Foam::HashSet< label, Hash< label > >
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
Foam::perfectInterface::writeDict
virtual void writeDict(Ostream &) const
Write dictionary.
Definition: perfectInterface.C:504
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:77
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
Foam::polyMesh::faceZones
const faceZoneMesh & faceZones() const
Return face zone mesh.
Definition: polyMesh.H:477
matchPoints.H
Determine correspondence between points. See below.
Foam::label
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:62
pFaces
Info<< "Finished reading KIVA file"<< endl;cellShapeList cellShapes(nPoints);labelList cellZoning(nPoints, -1);const cellModel &hex=cellModel::ref(cellModel::HEX);labelList hexLabels(8);label activeCells=0;labelList pointMap(nPoints);forAll(pointMap, i){ pointMap[i]=i;}for(label i=0;i< nPoints;i++){ if(f[i] > 0.0) { hexLabels[0]=i;hexLabels[1]=i1tab[i];hexLabels[2]=i3tab[i1tab[i]];hexLabels[3]=i3tab[i];hexLabels[4]=i8tab[i];hexLabels[5]=i1tab[i8tab[i]];hexLabels[6]=i3tab[i1tab[i8tab[i]]];hexLabels[7]=i3tab[i8tab[i]];cellShapes[activeCells]=cellShape(hex, hexLabels);edgeList edges=cellShapes[activeCells].edges();forAll(edges, ei) { if(edges[ei].mag(points)< SMALL) { label start=pointMap[edges[ei].start()];while(start !=pointMap[start]) { start=pointMap[start];} label end=pointMap[edges[ei].end()];while(end !=pointMap[end]) { end=pointMap[end];} label minLabel=min(start, end);pointMap[start]=pointMap[end]=minLabel;} } cellZoning[activeCells]=idreg[i];activeCells++;}}cellShapes.setSize(activeCells);cellZoning.setSize(activeCells);forAll(cellShapes, celli){ cellShape &cs=cellShapes[celli];forAll(cs, i) { cs[i]=pointMap[cs[i]];} cs.collapse();}label bcIDs[11]={-1, 0, 2, 4, -1, 5, -1, 6, 7, 8, 9};const label nBCs=12;const word *kivaPatchTypes[nBCs]={ &wallPolyPatch::typeName, &wallPolyPatch::typeName, &wallPolyPatch::typeName, &wallPolyPatch::typeName, &symmetryPolyPatch::typeName, &wedgePolyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &symmetryPolyPatch::typeName, &oldCyclicPolyPatch::typeName};enum patchTypeNames{ PISTON, VALVE, LINER, CYLINDERHEAD, AXIS, WEDGE, INFLOW, OUTFLOW, PRESIN, PRESOUT, SYMMETRYPLANE, CYCLIC};const char *kivaPatchNames[nBCs]={ "piston", "valve", "liner", "cylinderHead", "axis", "wedge", "inflow", "outflow", "presin", "presout", "symmetryPlane", "cyclic"};List< SLList< face > > pFaces[nBCs]
Definition: readKivaGrid.H:235
Foam::Field< vector >
Foam::faceZone
A subset of mesh faces organised as a primitive patch.
Definition: faceZone.H:65
Foam::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
Foam::polyMesh::faceOwner
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1076
Foam::IndirectList
A List with indirect addressing.
Definition: IndirectList.H:56
Foam::perfectInterface::changeTopology
virtual bool changeTopology() const
Check for topology change.
Definition: perfectInterface.C:137
Foam::polyBoundaryMesh::whichPatch
label whichPatch(const label faceIndex) const
Return patch index for a given face label.
Definition: polyBoundaryMesh.C:805
indirectPrimitivePatch.H
Foam::dictionary::lookup
ITstream & lookup(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionary.C:419
Foam::ZoneMesh::whichZone
label whichZone(const label objectIndex) const
Given a global object index, return the zone it is in.
Definition: ZoneMesh.C:315
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::Ostream::endBlock
virtual Ostream & endBlock()
Write end block group.
Definition: Ostream.C:109
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
zoneID
const labelIOList & zoneID
Definition: interpolatedFaces.H:22
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::faceZone::whichFace
label whichFace(const label globalCellID) const
Helper function to re-direct to zone::localID(...)
Definition: faceZone.C:388
Foam::polyPatch::start
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:311
polyModifyFace.H
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::polyMeshModifier
Virtual base class for mesh modifiers.
Definition: polyMeshModifier.H:71
Foam::indirectPrimitivePatch
PrimitivePatch< face, IndirectList, const pointField & > indirectPrimitivePatch
A PrimitivePatch using an IndirectList for the faces.
Definition: indirectPrimitivePatch.H:47
Foam::polyMesh::faces
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1063
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:355
Foam::nl
constexpr char nl
Definition: Ostream.H:372
Foam::polyTopoChanger::mesh
const polyMesh & mesh() const
Return the mesh reference.
Definition: polyTopoChanger.H:116
f
labelList f(nPoints)
Foam::polyRemoveFace
Class containing data for face removal.
Definition: polyRemoveFace.H:48
Foam::List< edge >
Foam::primitiveMesh::isInternalFace
bool isInternalFace(const label faceIndex) const
Return true if given face label is internal to the mesh.
Definition: primitiveMeshI.H:102
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::type
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:590
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::identity
labelList identity(const label len, label start=0)
Create identity map of the given length with (map[i] == i)
Definition: labelList.C:38
patches
const polyBoundaryMesh & patches
Definition: convertProcessorPatches.H:65
Foam::Ostream::writeEntry
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:219
Foam::mapPolyMesh
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:160
perfectInterface.H
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:74
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::PrimitivePatch::meshPoints
const labelList & meshPoints() const
Return labelList of mesh points in patch.
Definition: PrimitivePatch.C:418
Foam::faceZone::flipMap
const boolList & flipMap() const
Return face flip map.
Definition: faceZone.H:271
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::matchPoints
bool matchPoints(const UList< point > &pts0, const UList< point > &pts1, const UList< scalar > &matchDistance, const bool verbose, labelList &from0To1, const point &origin=point::zero)
Determine correspondence between pointFields. Gets passed.
Definition: matchPoints.C:35
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:294
Foam::polyMesh::faceNeighbour
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1082
Foam::PrimitivePatch
A list of faces which address into the list of points.
Definition: PrimitivePatch.H:90
Foam::perfectInterface::write
virtual void write(Ostream &) const
Write.
Definition: perfectInterface.C:494