Loading...
Searching...
No Matches
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-------------------------------------------------------------------------------
11License
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
27Class
28 Foam::cyclicAMIPolyPatch
29
30Description
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
45SourceFiles
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"
58#include "faceAreaWeightAMI.H"
59#include "cylindricalCS.H"
60
61// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
62
63namespace Foam
64{
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
80 (
81 const primitivePatch& half0,
82 const pointField& half0Ctrs,
83 const vectorField& half0Areas,
84 const pointField& half1Ctrs,
85 const vectorField& half1Areas
86 );
87
88
89protected:
90
91 // Protected Data
92
93 //- Name of other half
94 mutable word nbrPatchName_;
95
96 //- Optional patchGroup to find neighbPatch
99 //- Index of other half
100 mutable label nbrPatchID_;
101
102 //- Particle displacement fraction across 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
168
169 // Cache
171
172
173 // Protected Member Functions
174
175 // Topology change
176
177 //- Collect faces to remove in the topoChange container
178 virtual bool removeAMIFaces(polyTopoChange& topoChange);
179
180 //- Collect faces to add in the topoChange container
181 virtual bool addAMIFaces(polyTopoChange& topoChange);
182
183 //- Set properties of newly inserted faces after topological changes
184 virtual void setAMIFaces();
186 //- Helper to re-apply the geometric scaling lost during mesh
187 //- updates
188 virtual void restoreScaledGeometry();
190
191 //- Create and return pointer to the projection surface
192 const autoPtr<searchableSurface>& surfPtr() const;
193
194 //- Create a coordinate system from the periodic patch (or nullptr)
197 //- Reset the AMI interpolator, supply patch points
198 virtual void resetAMI(const UList<point>& points) const;
199
200 //- Reset the AMI interpolator, use current patch points
201 virtual void resetAMI() const;
202
203 //- Recalculate the transformation tensors
204 virtual void calcTransforms();
205
206 //- Initialise the calculation of the patch geometry
207 virtual void initGeometry(PstreamBuffers&);
208
209 //- Calculate the patch geometry
210 virtual void calcGeometry(PstreamBuffers&);
211
212 //- Initialise the patches for moving points
213 virtual void initMovePoints(PstreamBuffers& pBufs, const pointField&);
214
215 //- Correct patches after moving points
216 virtual void movePoints(PstreamBuffers& pBufs, const pointField&);
217
218 //- Initialise the update of the patch topology
219 virtual void initUpdateMesh(PstreamBuffers&);
220
221 //- Update of the patch topology
222 virtual void updateMesh(PstreamBuffers&);
223
224 //- Clear geometry
225 virtual void clearGeom();
226
227
228public:
229
230 //- Runtime type information
231 TypeName("cyclicAMI");
232
233
234 // Constructors
235
236 //- Construct from (base coupled patch) components
238 (
239 const word& name,
240 const label size,
241 const label start,
242 const label index,
243 const polyBoundaryMesh& bm,
244 const word& patchType,
246 const word& defaultAMIMethod = faceAreaWeightAMI::typeName
247 );
248
249 //- Construct from dictionary
251 (
252 const word& name,
253 const dictionary& dict,
254 const label index,
255 const polyBoundaryMesh& bm,
256 const word& patchType,
257 const word& defaultAMIMethod = faceAreaWeightAMI::typeName
258 );
259
260 //- Construct as copy, resetting the boundary mesh
262
263 //- Construct given the original patch and resetting the
264 // face list and boundary mesh information
266 (
267 const cyclicAMIPolyPatch& pp,
268 const polyBoundaryMesh& bm,
269 const label index,
270 const label newSize,
271 const label newStart,
272 const word& nbrPatchName
273 );
274
275 //- Construct given the original patch and a map
277 (
278 const cyclicAMIPolyPatch& pp,
279 const polyBoundaryMesh& bm,
280 const label index,
281 const labelUList& mapAddressing,
282 const label newStart
283 );
284
285
286 //- Construct and return a clone, resetting the boundary mesh
287 virtual autoPtr<polyPatch> clone(const polyBoundaryMesh& bm) const
288 {
290 }
291
292 //- Construct and return a clone, resetting the face list
293 //- and boundary mesh
295 (
296 const polyBoundaryMesh& bm,
297 const label index,
298 const label newSize,
299 const label newStart
300 ) const
301 {
302 return autoPtr<polyPatch>
303 (
305 (
306 *this,
307 bm,
308 index,
309 newSize,
310 newStart,
312 )
313 );
314 }
315
316 //- Construct and return a clone, resetting the face list
317 // and boundary mesh
319 (
320 const polyBoundaryMesh& bm,
321 const label index,
322 const labelUList& mapAddressing,
323 const label newStart
324 ) const
325 {
326 return autoPtr<polyPatch>
327 (
329 (
330 *this,
331 bm,
332 index,
333 mapAddressing,
334 newStart
335 )
336 );
337 }
338
339
340 //- Destructor
341 virtual ~cyclicAMIPolyPatch() = default;
342
343
344 // Member Functions
345
346 // Implicit treatment functions
347
348 //- Return number of new internal of this polyPatch faces
349 virtual void newInternalProcFaces(label&, label&) const;
350
351
352 //- Return nbrCells
353 virtual const labelUList& nbrCells() const
354 {
355 return neighbPatch().faceCells();
356 }
357
358 //- Return nbr patch ID
359 virtual label neighbPolyPatchID() const
360 {
361 return neighbPatch().index();
362 }
363
364 //- Return collocated faces map
365 virtual refPtr<labelListList> mapCollocatedFaces() const
366 {
367 const labelListList& sourceFaces = AMI().srcAddress();
368 return refPtr<labelListList>(new labelListList(sourceFaces));
369 }
370
371 //- Return implicit master
372 virtual bool masterImplicit() const
373 {
374 return owner();
375 }
376
378 // Access
379
380 //- Tolerance used e.g. for area calculations/limits
381 static const scalar tolerance_;
382
383 //- Flag to indicate whether the AMI can be reset
384 inline bool canResetAMI() const;
385
386 //- Return access to the createAMIFaces flag
387 inline bool createAMIFaces() const;
388
389 //- Return access to the updated flag
390 inline bool updatingAMI() const;
391
392 //- Return true if this patch changes the mesh topology
393 // True when createAMIFaces is true
394 virtual bool changeTopology() const;
395
396 //- Set topology changes in the polyTopoChange object
397 virtual bool setTopology(polyTopoChange& topoChange);
398
399 //- Is patch 'coupled'. Note that on AMI the geometry is not
400 //- coupled but the fields are!
401 virtual bool coupled() const
402 {
403 return false;
405
406 //- Neighbour patch name
407 inline const word& neighbPatchName() const;
408
409 //- Neighbour patch ID
410 virtual label neighbPatchID() const;
411
412 //- Particle fraction increase between AMI pathces
413 inline scalar fraction() const;
414
415 //- Does this side own the patch?
416 virtual bool owner() const;
417
418 //- Return a reference to the neighbour patch
419 virtual const cyclicAMIPolyPatch& neighbPatch() const;
420
421 //- Periodic patch ID (or -1)
422 label periodicPatchID() const;
423
424 //- Return a reference to the AMI interpolator
425 const AMIPatchToPatchInterpolation& AMI() const;
426
427 //- Helper function to return the weights
428 inline const scalarListList& weights() const;
430 //- Helper function to return the weights sum
431 inline const scalarField& weightsSum() const;
432
433 //- Return true if applying the low weight correction
434 bool applyLowWeightCorrection() const;
435
436 //- Return access to the initial face areas
437 // Used for topology change
438 inline vectorField& faceAreas0() const;
439
440 //- Return access to the initial face centres
441 // Used for topology change
442 inline vectorField& faceCentres0() const;
443
444
445 // Transformations
446
447 //- Axis of rotation for rotational cyclic AMI
448 inline const vector& rotationAxis() const;
449
450 //- Point on axis of rotation for rotational cyclic AMI
451 inline const point& rotationCentre() const;
452
453 //- Translation vector for translational cyclic AMI
454 inline const vector& separationVector() const;
455
456 //- Transform patch-based positions from nbr side to this side
457 virtual void transformPosition(pointField&) const;
458
459 //- Transform a patch-based position from nbr side to this side
460 virtual void transformPosition
462 point& l,
463 const label facei
464 ) const;
465
466 //- Transform a patch-based position from this side to nbr side
467 virtual void reverseTransformPosition
468 (
469 point& l,
470 const label facei
471 ) const;
472
473 //- Transform a patch-based direction from this side to
474 //- nbr side
475 virtual void reverseTransformDirection
476 (
477 vector& d,
478 const label facei
479 ) const;
480
482 // Interpolations
483
484 //- Interpolate field
485 template<class Type>
487 (
488 const Field<Type>& fld,
489 const UList<Type>& defaultValues = UList<Type>()
490 ) const;
491
492 //- Interpolate tmp field
493 template<class Type>
495 (
496 const tmp<Field<Type>>& tFld,
497 const UList<Type>& defaultValues = UList<Type>()
498 ) const;
499
500 //- Interpolate without periodic
501 template<class Type>
503 (
504 const Field<Type>& fld,
505 const UList<Type>& defaultValues
506 ) const;
507
508 //- Low-level interpolate List
509 template<class Type, class CombineOp>
510 void interpolate
511 (
512 const UList<Type>& fld,
513 const CombineOp& cop,
514 List<Type>& result,
515 const UList<Type>& defaultValues = UList<Type>()
516 ) const;
517
518
519 // Interpolations (non-blocking). Subject to change
520
521 template<class Type>
523 (
524 const Field<Type>& fld,
525 labelRange& sendRequests,
526 labelRange& recvRequests,
527 PtrList<List<Type>>& sendBuffers,
528 PtrList<List<Type>>& recvBuffers,
529
530 labelRange& sendRequests1,
531 labelRange& recvRequests1,
532 PtrList<List<Type>>& sendBuffers1,
533 PtrList<List<Type>>& recvBuffers1
534 ) const;
535
536 template<class Type>
537 void initInterpolate
538 (
539 const Field<Type>& fld,
540 labelRange& sendRequests,
541 labelRange& recvRequests,
542 PtrList<List<Type>>& sendBuffers,
543 PtrList<List<Type>>& recvBuffers,
544
545 labelRange& sendRequests1,
546 labelRange& recvRequests1,
547 PtrList<List<Type>>& sendBuffers1,
548 PtrList<List<Type>>& recvBuffers1
549 ) const;
550
551 template<class Type>
553 (
554 const Field<Type>& localFld,
555 const labelRange& requests, // The receive requests
556 const PtrList<List<Type>>& recvBuffers,
557 const labelRange& requests1, // The receive requests
558 const PtrList<List<Type>>& recvBuffers1,
559 const UList<Type>& defaultValues
560 ) const;
561
562
563 //- Calculate the patch geometry
564 virtual void calcGeometry
565 (
566 const primitivePatch& referPatch,
567 const pointField& thisCtrs,
568 const vectorField& thisAreas,
569 const pointField& thisCc,
570 const pointField& nbrCtrs,
571 const vectorField& nbrAreas,
572 const pointField& nbrCc
573 );
574
575 //- Initialize ordering for primitivePatch. Does not
576 //- refer to *this (except for name() and type() etc.)
577 virtual void initOrder
578 (
580 const primitivePatch&
581 ) const;
582
583 //- Return new ordering for primitivePatch.
584 // Ordering is -faceMap: for every face
585 // index of the new face -rotation:for every new face the clockwise
586 // shift of the original face. Return false if nothing changes
587 // (faceMap is identity, rotation is 0), true otherwise.
588 virtual bool order
589 (
591 const primitivePatch&,
593 labelList& rotation
594 ) const;
595
596 //- Return face index on neighbour patch which shares point p
597 //- following trajectory vector n
598 label pointFace
599 (
600 const label facei,
601 const vector& n,
602 point& p
603 ) const;
604
605 //- Write the polyPatch data as a dictionary
606 virtual void write(Ostream&) const;
607};
608
609
610// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
611
612} // End namespace Foam
613
614// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
615
616#include "cyclicAMIPolyPatchI.H"
617
618// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
619
620#ifdef NoRepository
622#endif
623
624// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
625
626#endif
627
628// ************************************************************************* //
label n
Info<< nl;Info<< "Write faMesh in vtk format:"<< nl;{ vtk::uindirectPatchWriter writer(aMesh.patch(), fileName(aMesh.time().globalPath()/vtkBaseFileName));writer.writeGeometry();globalIndex procAddr(aMesh.nFaces());labelList cellIDs;if(UPstream::master()) { cellIDs.resize(procAddr.totalSize());for(const labelRange &range :procAddr.ranges()) { auto slice=cellIDs.slice(range);slice=identity(range);} } writer.beginCellData(4);writer.writeProcIDs();writer.write("cellID", cellIDs);writer.write("area", aMesh.S().field());writer.write("normal", aMesh.faceAreaNormals());writer.beginPointData(1);writer.write("normal", aMesh.pointAreaNormals());Info<< " "<< writer.output().name()<< nl;}{ vtk::lineWriter writer(aMesh.points(), aMesh.edges(), fileName(aMesh.time().globalPath()/(vtkBaseFileName+"-edges")));writer.writeGeometry();writer.beginCellData(4);writer.writeProcIDs();{ Field< scalar > fld(faMeshTools::flattenEdgeField(aMesh.magLe(), true))
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
const labelListList & srcAddress() const
Return const access to source patch addressing.
Generic templated field type that is much like a Foam::List except that it is expected to hold numeri...
Definition Field.H:172
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
const Field< point_type > & points() const noexcept
Buffers for inter-processor communications streams (UOPstream, UIPstream).
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition PtrList.H:67
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition UList.H:89
void size(const label n)
Definition UList.H:118
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
static autoPtr< T > NewFrom(Args &&... args)
Construct autoPtr from derived type with forwarding arguments.
Definition autoPtr.H:193
Encapsulates using "patchGroups" to specify coupled patch.
virtual transformType transform() const
Type of transform.
coupledPolyPatch(const word &name, const label size, const label start, const label index, const polyBoundaryMesh &bm, const word &patchType, const transformType transform)
Construct from components.
Cyclic patch for Arbitrary Mesh Interface (AMI).
scalar rotationAngle_
Rotation angle.
virtual autoPtr< polyPatch > clone(const polyBoundaryMesh &bm) const
Construct and return a clone, resetting the boundary mesh.
vectorField & faceAreas0() const
Return access to the initial face areas.
virtual void updateMesh(PstreamBuffers &)
Update of the patch topology.
vector separationVector_
Translation vector.
word nbrPatchName_
Name of other half.
autoPtr< coordSystem::cylindrical > cylindricalCS() const
Create a coordinate system from the periodic patch (or nullptr).
virtual void resetAMI() const
Reset the AMI interpolator, use current patch points.
bool createAMIFaces_
Flag to indicate that new AMI faces will created.
virtual void initMovePoints(PstreamBuffers &pBufs, const pointField &)
Initialise the patches for moving points.
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 scalarField & weightsSum() const
Helper function to return the weights sum.
const word & neighbPatchName() const
Neighbour patch name.
const scalarListList & weights() const
Helper function to return the weights.
virtual bool owner() const
Does this side own the patch?
autoPtr< searchableSurface > surfPtr_
Projection surface.
virtual bool addAMIFaces(polyTopoChange &topoChange)
Collect faces to add in the topoChange container.
const scalar fraction_
Particle displacement fraction across AMI.
virtual label neighbPolyPatchID() const
Return nbr patch ID.
bool canResetAMI() const
Flag to indicate whether the AMI can be reset.
vectorField faceCentres0_
Temporary storage for AMI face centres.
virtual bool coupled() const
Is patch 'coupled'. Note that on AMI the geometry is not coupled but the fields are!
const vector & separationVector() const
Translation vector for translational cyclic AMI.
virtual bool masterImplicit() const
Return implicit master.
virtual void restoreScaledGeometry()
Helper to re-apply the geometric scaling lost during mesh updates.
label periodicPatchID_
Periodic patch.
virtual void clearGeom()
Clear geometry.
scalar fraction() const
Particle fraction increase between AMI pathces.
virtual bool removeAMIFaces(polyTopoChange &topoChange)
Collect faces to remove in the topoChange container.
virtual bool setTopology(polyTopoChange &topoChange)
Set topology changes in the polyTopoChange object.
virtual void calcGeometry(PstreamBuffers &)
Calculate the patch geometry.
const dictionary surfDict_
Dictionary used during projection surface construction.
virtual void movePoints(PstreamBuffers &pBufs, const pointField &)
Correct patches after moving points.
virtual void newInternalProcFaces(label &, label &) const
Return number of new internal of this polyPatch faces.
virtual void initGeometry(PstreamBuffers &)
Initialise the calculation of the patch geometry.
virtual ~cyclicAMIPolyPatch()=default
Destructor.
vectorField faceAreas0_
Temporary storage for AMI face areas.
virtual void initOrder(PstreamBuffers &, const primitivePatch &) const
Initialize ordering for primitivePatch. Does not refer to *this (except for name() and type() etc....
virtual void reverseTransformDirection(vector &d, const label facei) const
Transform a patch-based direction from this side to nbr side.
tmp< Field< Type > > interpolateUntransformed(const Field< Type > &fld, const UList< Type > &defaultValues) const
Interpolate without periodic.
virtual autoPtr< polyPatch > clone(const polyBoundaryMesh &bm, const label index, const label newSize, const label newStart) const
Construct and return a clone, resetting the face list and boundary mesh.
virtual const labelUList & nbrCells() const
Return nbrCells.
const point & rotationCentre() const
Point on axis of rotation for rotational cyclic AMI.
const vector & rotationAxis() const
Axis of rotation for rotational cyclic AMI.
const AMIPatchToPatchInterpolation & AMI() const
Return a reference to the AMI interpolator.
point rotationCentre_
Point on axis of rotation for rotational cyclics.
virtual bool changeTopology() const
Return true if this patch changes the mesh topology.
virtual const cyclicAMIPolyPatch & neighbPatch() const
Return a reference to the neighbour patch.
virtual void setAMIFaces()
Set properties of newly inserted faces after topological changes.
word periodicPatchName_
Periodic patch name.
virtual void transformPosition(pointField &) const
Transform patch-based positions from nbr side to this side.
label nbrPatchID_
Index of other half.
virtual refPtr< labelListList > mapCollocatedFaces() const
Return collocated faces map.
static const scalar tolerance_
Tolerance used e.g. for area calculations/limits.
TypeName("cyclicAMI")
Runtime type information.
tmp< Field< Type > > interpolate(const Field< Type > &fld, const UList< Type > &defaultValues=UList< Type >()) const
Interpolate field.
const autoPtr< searchableSurface > & surfPtr() const
Create and return pointer to the projection surface.
virtual autoPtr< polyPatch > clone(const polyBoundaryMesh &bm, const label index, const labelUList &mapAddressing, const label newStart) const
Construct and return a clone, resetting the face list.
bool rotationAngleDefined_
Flag to show whether the rotation angle is defined.
virtual bool order(PstreamBuffers &, const primitivePatch &, labelList &faceMap, labelList &rotation) const
Return new ordering for primitivePatch.
vectorField & faceCentres0() const
Return access to the initial face centres.
void initInterpolateUntransformed(const Field< Type > &fld, labelRange &sendRequests, labelRange &recvRequests, PtrList< List< Type > > &sendBuffers, PtrList< List< Type > > &recvBuffers, labelRange &sendRequests1, labelRange &recvRequests1, PtrList< List< Type > > &sendBuffers1, PtrList< List< Type > > &recvBuffers1) const
void initInterpolate(const Field< Type > &fld, labelRange &sendRequests, labelRange &recvRequests, PtrList< List< Type > > &sendBuffers, PtrList< List< Type > > &recvBuffers, labelRange &sendRequests1, labelRange &recvRequests1, PtrList< List< Type > > &sendBuffers1, PtrList< List< Type > > &recvBuffers1) const
bool createAMIFaces() const
Return access to the createAMIFaces flag.
tmp< Field< Type > > interpolate(const Field< Type > &localFld, const labelRange &requests, const PtrList< List< Type > > &recvBuffers, const labelRange &requests1, const PtrList< List< Type > > &recvBuffers1, const UList< Type > &defaultValues) const
virtual void initUpdateMesh(PstreamBuffers &)
Initialise the update of the patch topology.
bool applyLowWeightCorrection() const
Return true if applying the low weight correction.
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.
tmp< Field< Type > > interpolate(const tmp< Field< Type > > &tFld, const UList< Type > &defaultValues=UList< Type >()) const
Interpolate tmp field.
const coupleGroupIdentifier coupleGroup_
Optional patchGroup to find neighbPatch.
bool updatingAMI() const
Return access to the updated flag.
vector rotationAxis_
Axis of rotation for rotational cyclics.
virtual void reverseTransformPosition(point &l, const label facei) const
Transform a patch-based position from this side to nbr side.
autoPtr< AMIPatchToPatchInterpolation > AMIPtr_
AMI interpolation class.
bool moveFaceCentres_
Move face centres (default = no).
label periodicPatchID() const
Periodic patch ID (or -1).
virtual void calcTransforms()
Recalculate the transformation tensors.
virtual label neighbPatchID() const
Neighbour patch ID.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
A range or interval of labels defined by a start and a size.
Definition labelRange.H:66
label index() const noexcept
The index of this patch in the boundaryMesh.
A polyBoundaryMesh is a polyPatch list with registered IO, a reference to the associated polyMesh,...
const vectorField::subField faceCentres() const
Return face centres.
Definition polyPatch.C:320
label start() const noexcept
Return start label of this patch in the polyMesh face list.
Definition polyPatch.H:446
const labelUList & faceCells() const
Return face-cell addressing.
Definition polyPatch.C:401
Direct mesh changes based on v1.3 polyTopoChange syntax.
A class for managing references or pointers (no reference counting).
Definition refPtr.H:54
A class for managing temporary objects.
Definition tmp.H:75
A class for handling words, derived from Foam::string.
Definition word.H:66
volScalarField & p
Namespace for OpenFOAM.
List< scalarList > scalarListList
List of scalarList.
Definition scalarList.H:35
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
List< labelList > labelListList
List of labelList.
Definition labelList.H:38
List< label > labelList
A List of labels.
Definition List.H:62
PrimitivePatch< SubList< face >, const pointField & > primitivePatch
A PrimitivePatch with a SubList addressing for the faces, const reference for the point field.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
AMIInterpolation AMIPatchToPatchInterpolation
Patch-to-patch interpolation == Foam::AMIInterpolation.
Field< vector > vectorField
Specialisation of Field<T> for vector.
vector point
Point is a vector.
Definition point.H:37
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
vectorField pointField
pointField is a vectorField.
Vector< scalar > vector
Definition vector.H:57
UList< label > labelUList
A UList of labels.
Definition UList.H:75
runTime write()
dictionary dict
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition typeInfo.H:68