cyclicAMIPolyPatch.H
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-2023 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 Class
28  Foam::cyclicAMIPolyPatch
29 
30 Description
31  Cyclic patch for Arbitrary Mesh Interface (AMI)
32 
33  Includes provision for updating the patch topology to enforce a 1-to-1
34  face match across the interface, based on the \c createAMIFaces flag.
35 
36  The manipulations are based on the reference:
37 
38  \verbatim
39  H.J. Aguerre, S. Márquez Damián, J.M. Gimenez, N.M.Nigro, Conservative
40  handling of arbitrary non-conformal interfaces using an efficient
41  supermesh, Journal of Computational Physics 335(15) 21-49. 2017.
42  https://doi.org/10.1016/j.jcp.2017.01.018.
43  \endverbatim
44 
45 SourceFiles
46  cyclicAMIPolyPatch.C
47 
48 \*---------------------------------------------------------------------------*/
49 
50 #ifndef Foam_cyclicAMIPolyPatch_H
51 #define Foam_cyclicAMIPolyPatch_H
52 
53 #include "AMIInterpolation.H"
54 #include "coupledPolyPatch.H"
56 #include "polyBoundaryMesh.H"
57 #include "coupleGroupIdentifier.H"
58 #include "faceAreaWeightAMI.H"
59 #include "cylindricalCS.H"
60 
61 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
62 
63 namespace Foam
64 {
65 
66 /*---------------------------------------------------------------------------*\
67  Class cyclicAMIPolyPatch Declaration
68 \*---------------------------------------------------------------------------*/
69 
71 :
72  public coupledPolyPatch
73 {
74  // Private Member Functions
75 
76  //- Return normal of face at max distance from rotation axis
77  vector findFaceNormalMaxRadius(const pointField& faceCentres) const;
78 
79  void calcTransforms
80  (
81  const primitivePatch& half0,
82  const pointField& half0Ctrs,
83  const vectorField& half0Areas,
84  const pointField& half1Ctrs,
85  const vectorField& half1Areas
86  );
87 
88 
89 protected:
90 
91  // Protected Data
92 
93  //- Name of other half
94  mutable word nbrPatchName_;
95 
96  //- Optional patchGroup to find neighbPatch
98 
99  //- Index of other half
100  mutable label nbrPatchID_;
101 
102  //- Particle displacement fraction accross AMI
103  const scalar fraction_;
104 
105 
106  // Transformations
107 
108  // For rotation
109 
110  //- Axis of rotation for rotational cyclics
112 
113  //- Point on axis of rotation for rotational cyclics
115 
116  //- Flag to show whether the rotation angle is defined
119  //- Rotation angle
120  scalar rotationAngle_;
121 
122 
123  // For translation
124 
125  //- Translation vector
127 
129  // Triggering periodic interpolation
130 
131  //- Periodic patch name
134  //- Periodic patch
135  mutable label periodicPatchID_;
136 
137 
138  //- AMI interpolation class
140 
141  //- Dictionary used during projection surface construction
142  const dictionary surfDict_;
143 
144  //- Projection surface
146 
147 
148  // Change of topology as AMI is updated
150  //- Flag to indicate that new AMI faces will created
151  // Set by the call to changeTopology
152  mutable bool createAMIFaces_;
153 
154  //- Move face centres (default = no)
155  bool moveFaceCentres_;
156 
157  mutable bool updatingAMI_;
158 
162 
163  //- Temporary storage for AMI face areas
164  mutable vectorField faceAreas0_;
166  //- Temporary storage for AMI face centres
167  mutable vectorField faceCentres0_;
168 
169 
170  // Protected Member Functions
171 
172  // Topology change
173 
174  //- Collect faces to remove in the topoChange container
175  virtual bool removeAMIFaces(polyTopoChange& topoChange);
176 
177  //- Collect faces to add in the topoChange container
178  virtual bool addAMIFaces(polyTopoChange& topoChange);
179 
180  //- Set properties of newly inserted faces after topological changes
181  virtual void setAMIFaces();
182 
183  //- Helper to re-apply the geometric scaling lost during mesh
184  //- updates
185  virtual void restoreScaledGeometry();
186 
188  //- Create and return pointer to the projection surface
190 
191  //- Create a coordinate system from the periodic patch (or nullptr)
193 
194  //- Reset the AMI interpolator, supply patch points
195  virtual void resetAMI(const UList<point>& points) const;
197  //- Reset the AMI interpolator, use current patch points
198  virtual void resetAMI() const;
199 
200  //- Recalculate the transformation tensors
201  virtual void calcTransforms();
202 
203  //- Initialise the calculation of the patch geometry
204  virtual void initGeometry(PstreamBuffers&);
205 
206  //- Calculate the patch geometry
207  virtual void calcGeometry(PstreamBuffers&);
208 
209  //- Initialise the patches for moving points
210  virtual void initMovePoints(PstreamBuffers& pBufs, const pointField&);
211 
212  //- Correct patches after moving points
213  virtual void movePoints(PstreamBuffers& pBufs, const pointField&);
214 
215  //- Initialise the update of the patch topology
216  virtual void initUpdateMesh(PstreamBuffers&);
217 
218  //- Update of the patch topology
219  virtual void updateMesh(PstreamBuffers&);
220 
221  //- Clear geometry
222  virtual void clearGeom();
223 
224 
225 public:
226 
227  //- Runtime type information
228  TypeName("cyclicAMI");
229 
230 
231  // Constructors
232 
233  //- Construct from (base coupled patch) components
235  (
236  const word& name,
237  const label size,
238  const label start,
239  const label index,
240  const polyBoundaryMesh& bm,
241  const word& patchType,
243  const word& defaultAMIMethod = faceAreaWeightAMI::typeName
244  );
245 
246  //- Construct from dictionary
248  (
249  const word& name,
250  const dictionary& dict,
251  const label index,
252  const polyBoundaryMesh& bm,
253  const word& patchType,
254  const word& defaultAMIMethod = faceAreaWeightAMI::typeName
255  );
256 
257  //- Construct as copy, resetting the boundary mesh
259 
260  //- Construct given the original patch and resetting the
261  // face list and boundary mesh information
263  (
264  const cyclicAMIPolyPatch& pp,
265  const polyBoundaryMesh& bm,
266  const label index,
267  const label newSize,
268  const label newStart,
269  const word& nbrPatchName
270  );
271 
272  //- Construct given the original patch and a map
274  (
275  const cyclicAMIPolyPatch& pp,
276  const polyBoundaryMesh& bm,
277  const label index,
278  const labelUList& mapAddressing,
279  const label newStart
280  );
281 
282 
283  //- Construct and return a clone, resetting the boundary mesh
284  virtual autoPtr<polyPatch> clone(const polyBoundaryMesh& bm) const
285  {
287  }
288 
289  //- Construct and return a clone, resetting the face list
290  //- and boundary mesh
291  virtual autoPtr<polyPatch> clone
292  (
293  const polyBoundaryMesh& bm,
294  const label index,
295  const label newSize,
296  const label newStart
297  ) const
298  {
299  return autoPtr<polyPatch>
300  (
302  (
303  *this,
304  bm,
305  index,
306  newSize,
307  newStart,
309  )
310  );
311  }
312 
313  //- Construct and return a clone, resetting the face list
314  // and boundary mesh
315  virtual autoPtr<polyPatch> clone
316  (
317  const polyBoundaryMesh& bm,
318  const label index,
319  const labelUList& mapAddressing,
320  const label newStart
321  ) const
322  {
323  return autoPtr<polyPatch>
324  (
326  (
327  *this,
328  bm,
329  index,
330  mapAddressing,
331  newStart
332  )
333  );
334  }
335 
336 
337  //- Destructor
338  virtual ~cyclicAMIPolyPatch() = default;
339 
340 
341  // Member Functions
342 
343  // Implicit treatment functions
344 
345  //- Return number of new internal of this polyPatch faces
346  virtual void newInternalProcFaces(label&, label&) const;
347 
348 
349  //- Return nbrCells
350  virtual const labelUList& nbrCells() const
351  {
352  return neighbPatch().faceCells();
353  }
354 
355  //- Return nbr patch ID
356  virtual label neighbPolyPatchID() const
357  {
358  return neighbPatch().index();
359  }
360 
361  //- Return collocated faces map
362  virtual refPtr<labelListList> mapCollocatedFaces() const
363  {
364  const labelListList& sourceFaces = AMI().srcAddress();
365  return refPtr<labelListList>(new labelListList(sourceFaces));
366  }
367 
368  //- Return implicit master
369  virtual bool masterImplicit() const
370  {
371  return owner();
372  }
373 
374 
375  // Access
376 
377  //- Tolerance used e.g. for area calculations/limits
378  static const scalar tolerance_;
379 
380  //- Flag to indicate whether the AMI can be reset
381  inline bool canResetAMI() const;
382 
383  //- Return access to the createAMIFaces flag
384  inline bool createAMIFaces() const;
385 
386  //- Return access to the updated flag
387  inline bool updatingAMI() const;
388 
389  //- Return true if this patch changes the mesh topology
390  // True when createAMIFaces is true
391  virtual bool changeTopology() const;
392 
393  //- Set topology changes in the polyTopoChange object
394  virtual bool setTopology(polyTopoChange& topoChange);
395 
396  //- Is patch 'coupled'. Note that on AMI the geometry is not
397  //- coupled but the fields are!
398  virtual bool coupled() const
399  {
400  return false;
401  }
403  //- Neighbour patch name
404  inline const word& neighbPatchName() const;
405 
406  //- Neighbour patch ID
407  virtual label neighbPatchID() const;
408 
409  //- Particle fraction increase between AMI pathces
410  inline scalar fraction() const;
411 
412  //- Does this side own the patch?
413  virtual bool owner() const;
414 
415  //- Return a reference to the neighbour patch
416  virtual const cyclicAMIPolyPatch& neighbPatch() const;
417 
418  //- Periodic patch ID (or -1)
419  label periodicPatchID() const;
420 
421  //- Return a reference to the AMI interpolator
422  const AMIPatchToPatchInterpolation& AMI() const;
423 
424  //- Helper function to return the weights
425  inline const scalarListList& weights() const;
426 
427  //- Helper function to return the weights sum
428  inline const scalarField& weightsSum() const;
429 
430  //- Return true if applying the low weight correction
431  bool applyLowWeightCorrection() const;
432 
433  //- Return access to the initial face areas
434  // Used for topology change
435  inline vectorField& faceAreas0() const;
436 
437  //- Return access to the initial face centres
438  // Used for topology change
439  inline vectorField& faceCentres0() const;
440 
441 
442  // Transformations
443 
444  //- Axis of rotation for rotational cyclic AMI
445  inline const vector& rotationAxis() const;
446 
447  //- Point on axis of rotation for rotational cyclic AMI
448  inline const point& rotationCentre() const;
449 
450  //- Translation vector for translational cyclic AMI
451  inline const vector& separationVector() const;
452 
453  //- Transform patch-based positions from nbr side to this side
454  virtual void transformPosition(pointField&) const;
455 
456  //- Transform a patch-based position from nbr side to this side
457  virtual void transformPosition
458  (
459  point& l,
460  const label facei
461  ) const;
462 
463  //- Transform a patch-based position from this side to nbr side
464  virtual void reverseTransformPosition
465  (
466  point& l,
467  const label facei
468  ) const;
469 
470  //- Transform a patch-based direction from this side to
471  //- nbr side
472  virtual void reverseTransformDirection
473  (
474  vector& d,
475  const label facei
476  ) const;
477 
479  // Interpolations
480 
481  //- Interpolate field
482  template<class Type>
484  (
485  const Field<Type>& fld,
486  const UList<Type>& defaultValues = UList<Type>()
487  ) const;
488 
489  //- Interpolate tmp field
490  template<class Type>
492  (
493  const tmp<Field<Type>>& tFld,
494  const UList<Type>& defaultValues = UList<Type>()
495  ) const;
496 
497  //- Interpolate without periodic
498  template<class Type>
500  (
501  const Field<Type>& fld,
502  const UList<Type>& defaultValues
503  ) const;
504 
505  //- Low-level interpolate List
506  template<class Type, class CombineOp>
507  void interpolate
508  (
509  const UList<Type>& fld,
510  const CombineOp& cop,
511  List<Type>& result,
512  const UList<Type>& defaultValues = UList<Type>()
513  ) const;
514 
515 
516  // Interpolations (non-blocking). Subject to change
517 
518  template<class Type>
520  (
521  const Field<Type>& fld,
522  labelRange& sendRequests,
523  PtrList<List<Type>>& sendBuffers,
524  labelRange& recvRequests,
525  PtrList<List<Type>>& recvBuffers
526  ) const;
527 
528  template<class Type>
529  void initInterpolate
530  (
531  const Field<Type>& fld,
532  labelRange& sendRequests,
533  PtrList<List<Type>>& sendBuffers,
534  labelRange& recvRequests,
535  PtrList<List<Type>>& recvBuffers
536  ) const;
537 
538  template<class Type>
540  (
541  const Field<Type>& localFld,
542  const labelRange& requests, // The receive requests
543  const PtrList<List<Type>>& recvBuffers,
544  const UList<Type>& defaultValues
545  ) const;
546 
547 
548  //- Calculate the patch geometry
549  virtual void calcGeometry
550  (
551  const primitivePatch& referPatch,
552  const pointField& thisCtrs,
553  const vectorField& thisAreas,
554  const pointField& thisCc,
555  const pointField& nbrCtrs,
556  const vectorField& nbrAreas,
557  const pointField& nbrCc
558  );
559 
560  //- Initialize ordering for primitivePatch. Does not
561  //- refer to *this (except for name() and type() etc.)
562  virtual void initOrder
563  (
565  const primitivePatch&
566  ) const;
567 
568  //- Return new ordering for primitivePatch.
569  // Ordering is -faceMap: for every face
570  // index of the new face -rotation:for every new face the clockwise
571  // shift of the original face. Return false if nothing changes
572  // (faceMap is identity, rotation is 0), true otherwise.
573  virtual bool order
574  (
576  const primitivePatch&,
578  labelList& rotation
579  ) const;
580 
581  //- Return face index on neighbour patch which shares point p
582  //- following trajectory vector n
583  label pointFace
584  (
585  const label facei,
586  const vector& n,
587  point& p
588  ) const;
589 
590  //- Write the polyPatch data as a dictionary
591  virtual void write(Ostream&) const;
592 };
593 
594 
595 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
596 
597 } // End namespace Foam
598 
599 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
600 
601 #include "cyclicAMIPolyPatchI.H"
602 
603 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
604 
605 #ifdef NoRepository
607 #endif
608 
609 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
610 
611 #endif
612 
613 // ************************************************************************* //
virtual void clearGeom()
Clear geometry.
virtual void resetAMI() const
Reset the AMI interpolator, use current patch points.
vector separationVector_
Translation vector.
dictionary dict
virtual refPtr< labelListList > mapCollocatedFaces() const
Return collocated faces map.
const AMIPatchToPatchInterpolation & AMI() const
Return a reference to the AMI interpolator.
virtual label neighbPatchID() const
Neighbour patch ID.
virtual bool order(PstreamBuffers &, const primitivePatch &, labelList &faceMap, labelList &rotation) const
Return new ordering for primitivePatch.
cyclicAMIPolyPatch(const word &name, const label size, const label start, const label index, const polyBoundaryMesh &bm, const word &patchType, const transformType transform=UNKNOWN, const word &defaultAMIMethod=faceAreaWeightAMI::typeName)
Construct from (base coupled patch) components.
autoPtr< coordSystem::cylindrical > cylindricalCS() const
Create a coordinate system from the periodic patch (or nullptr)
friend class polyBoundaryMesh
Definition: polyPatch.H:107
virtual void initUpdateMesh(PstreamBuffers &)
Initialise the update of the patch topology.
virtual void initGeometry(PstreamBuffers &)
Initialise the calculation of the patch geometry.
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
label start() const noexcept
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:441
virtual bool coupled() const
Is patch &#39;coupled&#39;. Note that on AMI the geometry is not coupled but the fields are! ...
bool moveFaceCentres_
Move face centres (default = no)
virtual bool addAMIFaces(polyTopoChange &topoChange)
Collect faces to add in the topoChange container.
A range or interval of labels defined by a start and a size.
Definition: labelRange.H:52
virtual void write(Ostream &) const
Write the polyPatch data as a dictionary.
static const scalar tolerance_
Tolerance used e.g. for area calculations/limits.
virtual void calcTransforms()
Recalculate the transformation tensors.
virtual autoPtr< polyPatch > clone(const polyBoundaryMesh &bm) const
Construct and return a clone, resetting the boundary mesh.
virtual void transformPosition(pointField &) const
Transform patch-based positions from nbr side to this side.
bool applyLowWeightCorrection() const
Return true if applying the low weight correction.
vectorField & faceCentres0() const
Return access to the initial face centres.
The coupledPolyPatch is an abstract base class for patches that couple regions of the computational d...
vector rotationAxis_
Axis of rotation for rotational cyclics.
bool rotationAngleDefined_
Flag to show whether the rotation angle is defined.
label periodicPatchID_
Periodic patch.
virtual void reverseTransformDirection(vector &d, const label facei) const
Transform a patch-based direction from this side to nbr side.
tmp< Field< Type > > interpolate(const Field< Type > &fld, const UList< Type > &defaultValues=UList< Type >()) const
Interpolate field.
List< labelList > labelListList
List of labelList.
Definition: labelList.H:38
virtual bool setTopology(polyTopoChange &topoChange)
Set topology changes in the polyTopoChange object.
scalar fraction() const
Particle fraction increase between AMI pathces.
A class for managing references or pointers (no reference counting)
Definition: HashPtrTable.H:49
UList< label > labelUList
A UList of labels.
Definition: UList.H:78
virtual void newInternalProcFaces(label &, label &) const
Return number of new internal of this polyPatch faces.
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
bool createAMIFaces_
Flag to indicate that new AMI faces will created.
void initInterpolate(const Field< Type > &fld, labelRange &sendRequests, PtrList< List< Type >> &sendBuffers, labelRange &recvRequests, PtrList< List< Type >> &recvBuffers) const
const dictionary surfDict_
Dictionary used during projection surface construction.
virtual void calcGeometry(PstreamBuffers &)
Calculate the patch geometry.
label nbrPatchID_
Index of other half.
virtual void setAMIFaces()
Set properties of newly inserted faces after topological changes.
vectorField faceAreas0_
Temporary storage for AMI face areas.
const vector & separationVector() const
Translation vector for translational cyclic AMI.
const coupleGroupIdentifier coupleGroup_
Optional patchGroup to find neighbPatch.
virtual ~cyclicAMIPolyPatch()=default
Destructor.
A list of faces which address into the list of points.
const vector & rotationAxis() const
Axis of rotation for rotational cyclic AMI.
const scalar fraction_
Particle displacement fraction accross AMI.
const word & neighbPatchName() const
Neighbour patch name.
const point & rotationCentre() const
Point on axis of rotation for rotational cyclic AMI.
const labelUList & faceCells() const
Return face-cell addressing.
Definition: polyPatch.C:365
bool updatingAMI() const
Return access to the updated flag.
A class for handling words, derived from Foam::string.
Definition: word.H:63
const autoPtr< searchableSurface > & surfPtr() const
Create and return pointer to the projection surface.
autoPtr< AMIPatchToPatchInterpolation > AMIPtr_
AMI interpolation class.
point rotationCentre_
Point on axis of rotation for rotational cyclics.
virtual void updateMesh(PstreamBuffers &)
Update of the patch topology.
bool createAMIFaces() const
Return access to the createAMIFaces flag.
word periodicPatchName_
Periodic patch name.
Cyclic patch for Arbitrary Mesh Interface (AMI)
const Field< point_type > & points() const noexcept
Return reference to global points.
label periodicPatchID() const
Periodic patch ID (or -1)
bool canResetAMI() const
Flag to indicate whether the AMI can be reset.
void initInterpolateUntransformed(const Field< Type > &fld, labelRange &sendRequests, PtrList< List< Type >> &sendBuffers, labelRange &recvRequests, PtrList< List< Type >> &recvBuffers) const
vectorField & faceAreas0() const
Return access to the initial face areas.
virtual const labelUList & nbrCells() const
Return nbrCells.
const scalarField & weightsSum() const
Helper function to return the weights sum.
Encapsulates using "patchGroups" to specify coupled patch.
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: HashTable.H:105
vectorField faceCentres0_
Temporary storage for AMI face centres.
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO...
const word & name() const noexcept
The patch name.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
Buffers for inter-processor communications streams (UOPstream, UIPstream).
virtual void restoreScaledGeometry()
Helper to re-apply the geometric scaling lost during mesh updates.
virtual bool removeAMIFaces(polyTopoChange &topoChange)
Collect faces to remove in the topoChange container.
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< ' ';}gmvFile<< nl;for(const word &name :lagrangianScalarNames){ IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
tmp< Field< Type > > interpolateUntransformed(const Field< Type > &fld, const UList< Type > &defaultValues) const
Interpolate without periodic.
TypeName("cyclicAMI")
Runtime type information.
virtual void initOrder(PstreamBuffers &, const primitivePatch &) const
Initialize ordering for primitivePatch. Does not refer to *this (except for name() and type() etc...
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers...
Definition: List.H:55
virtual const cyclicAMIPolyPatch & neighbPatch() const
Return a reference to the neighbour patch.
virtual bool owner() const
Does this side own the patch?
Direct mesh changes based on v1.3 polyTopoChange syntax.
virtual bool changeTopology() const
Return true if this patch changes the mesh topology.
Interpolation class dealing with transfer of data between two primitive patches with an arbitrary mes...
label n
virtual transformType transform() const
Type of transform.
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
label index() const noexcept
The index of this patch in the boundaryMesh.
volScalarField & p
A class for managing temporary objects.
Definition: HashPtrTable.H:50
const labelListList & srcAddress() const
Return const access to source patch addressing.
word nbrPatchName_
Name of other half.
virtual bool masterImplicit() const
Return implicit master.
virtual label neighbPolyPatchID() const
Return nbr patch ID.
virtual void reverseTransformPosition(point &l, const label facei) const
Transform a patch-based position from this side to nbr side.
scalar rotationAngle_
Rotation angle.
label pointFace(const label facei, const vector &n, point &p) const
Return face index on neighbour patch which shares point p following trajectory vector n...
const scalarListList & weights() const
Helper function to return the weights.
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
Namespace for OpenFOAM.
virtual void movePoints(PstreamBuffers &pBufs, const pointField &)
Correct patches after moving points.
virtual void initMovePoints(PstreamBuffers &pBufs, const pointField &)
Initialise the patches for moving points.
autoPtr< searchableSurface > surfPtr_
Projection surface.
const vectorField::subField faceCentres() const
Return face centres.
Definition: polyPatch.C:315