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-2021 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 cyclicAMIPolyPatch_H
51#define cyclicAMIPolyPatch_H
52
53#include "coupledPolyPatch.H"
55#include "polyBoundaryMesh.H"
57#include "faceAreaWeightAMI.H"
58#include "cylindricalCS.H"
59
60// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
61
62namespace Foam
63{
64
65/*---------------------------------------------------------------------------*\
66 Class cyclicAMIPolyPatch Declaration
67\*---------------------------------------------------------------------------*/
70:
71 public coupledPolyPatch
72{
73 // Private Member Functions
74
75 //- Return normal of face at max distance from rotation axis
76 vector findFaceNormalMaxRadius(const pointField& faceCentres) const;
77
79 (
80 const primitivePatch& half0,
81 const pointField& half0Ctrs,
82 const vectorField& half0Areas,
83 const pointField& half1Ctrs,
84 const vectorField& half1Areas
85 );
86
87
88protected:
89
90 // Protected data
91
92 //- Name of other half
93 mutable word nbrPatchName_;
94
95 //- Optional patchGroup to find neighbPatch
97
98 //- Index of other half
99 mutable label nbrPatchID_;
100
101 //- Particle displacement fraction accross AMI
102 const scalar fraction_;
103
104
105 // Transformations
106
107 // For rotation
108
109 //- Axis of rotation for rotational cyclics
111
112 //- Point on axis of rotation for rotational cyclics
114
115 //- Flag to show whether the rotation angle is defined
117
118 //- Rotation angle
119 scalar rotationAngle_;
120
121
122 // For translation
123
124 //- Translation vector
126
127
128 // Triggering periodic interpolation
129
130 //- Periodic patch name
132
133 //- Periodic patch
134 mutable label periodicPatchID_;
135
136
137 //- AMI interpolation class
139
140 //- Dictionary used during projection surface construction
141 const dictionary surfDict_;
142
143 //- Projection surface
145
146
147 // Change of topology as AMI is updated
148
149 //- Flag to indicate that new AMI faces will created
150 // Set by the call to changeTopology
151 mutable bool createAMIFaces_;
152
153 //- Move face centres (default = no)
154 bool moveFaceCentres_;
156 mutable bool updatingAMI_;
161
162 //- Temporary storage for AMI face areas
163 mutable vectorField faceAreas0_;
164
165 //- Temporary storage for AMI face centres
167
168
169 // Protected Member Functions
170
171 // Topology change
172
173 //- Collect faces to remove in the topoChange container
174 virtual bool removeAMIFaces(polyTopoChange& topoChange);
175
176 //- Collect faces to add in the topoChange container
177 virtual bool addAMIFaces(polyTopoChange& topoChange);
178
179 //- Set properties of newly inserted faces after topological changes
180 virtual void setAMIFaces();
181
182 //- Helper to re-apply the geometric scaling lost during mesh
183 //- updates
184 virtual void restoreScaledGeometry();
185
186
187 //- Create and return pointer to the projection surface
188 const autoPtr<searchableSurface>& surfPtr() const;
189
190 //- Create a coordinate system from the periodic patch (or nullptr)
192
193 //- Reset the AMI interpolator, supply patch points
194 virtual void resetAMI(const UList<point>& points) const;
195
196 //- Reset the AMI interpolator, use current patch points
197 virtual void resetAMI() const;
198
199 //- Recalculate the transformation tensors
200 virtual void calcTransforms();
201
202 //- Initialise the calculation of the patch geometry
203 virtual void initGeometry(PstreamBuffers&);
204
205 //- Calculate the patch geometry
206 virtual void calcGeometry(PstreamBuffers&);
207
208 //- Initialise the patches for moving points
209 virtual void initMovePoints(PstreamBuffers& pBufs, const pointField&);
210
211 //- Correct patches after moving points
212 virtual void movePoints(PstreamBuffers& pBufs, const pointField&);
213
214 //- Initialise the update of the patch topology
215 virtual void initUpdateMesh(PstreamBuffers&);
216
217 //- Update of the patch topology
218 virtual void updateMesh(PstreamBuffers&);
219
220 //- Clear geometry
221 virtual void clearGeom();
222
223
224public:
225
226 //- Runtime type information
227 TypeName("cyclicAMI");
228
229
230 // Constructors
231
232 //- Construct from (base coupled patch) components
234 (
235 const word& name,
236 const label size,
237 const label start,
238 const label index,
239 const polyBoundaryMesh& bm,
240 const word& patchType,
242 const word& defaultAMIMethod = faceAreaWeightAMI::typeName
243 );
244
245 //- Construct from dictionary
247 (
248 const word& name,
249 const dictionary& dict,
250 const label index,
251 const polyBoundaryMesh& bm,
252 const word& patchType,
253 const word& defaultAMIMethod = faceAreaWeightAMI::typeName
254 );
255
256 //- Construct as copy, resetting the boundary mesh
258
259 //- Construct given the original patch and resetting the
260 // face list and boundary mesh information
262 (
263 const cyclicAMIPolyPatch& pp,
264 const polyBoundaryMesh& bm,
265 const label index,
266 const label newSize,
267 const label newStart,
268 const word& nbrPatchName
269 );
270
271 //- Construct given the original patch and a map
273 (
274 const cyclicAMIPolyPatch& pp,
275 const polyBoundaryMesh& bm,
276 const label index,
277 const labelUList& mapAddressing,
278 const label newStart
279 );
280
281
282 //- Construct and return a clone, resetting the boundary mesh
283 virtual autoPtr<polyPatch> clone(const polyBoundaryMesh& bm) const
284 {
286 }
287
288 //- Construct and return a clone, resetting the face list
289 // and boundary mesh
291 (
292 const polyBoundaryMesh& bm,
293 const label index,
294 const label newSize,
295 const label newStart
296 ) const
297 {
298 return autoPtr<polyPatch>
299 (
301 (
302 *this,
303 bm,
304 index,
305 newSize,
306 newStart,
308 )
309 );
310 }
311
312 //- Construct and return a clone, resetting the face list
313 // and boundary mesh
315 (
316 const polyBoundaryMesh& bm,
317 const label index,
318 const labelUList& mapAddressing,
319 const label newStart
320 ) const
321 {
322 return autoPtr<polyPatch>
323 (
325 (
326 *this,
327 bm,
328 index,
329 mapAddressing,
330 newStart
331 )
332 );
333 }
334
335
336 //- Destructor
337 virtual ~cyclicAMIPolyPatch() = default;
338
339
340 // Member Functions
341
342 // Implicit treatment functions
343
344 //- Return number of new internal of this polyPatch faces
345 virtual void newInternalProcFaces(label&, label&) const;
346
347
348 //- Return nbrCells
349 virtual const labelUList& nbrCells() const
350 {
351 return neighbPatch().faceCells();
352 }
353
354 //- Return nbr patch ID
355 virtual label neighbPolyPatchID() const
356 {
357 return neighbPatch().index();
358 }
359
360 //- Return collocated faces map
362 {
363 const labelListList& sourceFaces = AMI().srcAddress();
364 return refPtr<labelListList>(new labelListList(sourceFaces));
365 }
366
367 //- Return implicit master
368 virtual bool masterImplicit() const
369 {
370 return owner();
371 }
372
373
374 // Access
375
376 //- Tolerance used e.g. for area calculations/limits
377 static const scalar tolerance_;
378
379 //- Flag to indicate whether the AMI can be reset
380 inline bool canResetAMI() const;
381
382 //- Return access to the createAMIFaces flag
383 inline bool createAMIFaces() const;
384
385 //- Return access to the updated flag
386 inline bool updatingAMI() const;
387
388 //- Return true if this patch changes the mesh topology
389 // True when createAMIFaces is true
390 virtual bool changeTopology() const;
391
392 //- Set topology changes in the polyTopoChange object
393 virtual bool setTopology(polyTopoChange& topoChange);
394
395 //- Is patch 'coupled'. Note that on AMI the geometry is not
396 //- coupled but the fields are!
397 virtual bool coupled() const
398 {
399 return false;
400 }
401
402 //- Neighbour patch name
403 inline const word& neighbPatchName() const;
404
405 //- Neighbour patch ID
406 virtual label neighbPatchID() const;
407
408 //- Particle fraction increase between AMI pathces
409 inline scalar fraction() const;
410
411 //- Does this side own the patch?
412 virtual bool owner() const;
413
414 //- Return a reference to the neighbour patch
415 virtual const cyclicAMIPolyPatch& neighbPatch() const;
416
417 //- Periodic patch ID (or -1)
418 label periodicPatchID() const;
419
420 //- Return a reference to the AMI interpolator
421 const AMIPatchToPatchInterpolation& AMI() const;
422
423 //- Helper function to return the weights
424 inline const scalarListList& weights() const;
425
426 //- Helper function to return the weights sum
427 inline const scalarField& weightsSum() const;
428
429 //- Return true if applying the low weight correction
430 bool applyLowWeightCorrection() const;
431
432 //- Return access to the initial face areas
433 // Used for topology change
434 inline vectorField& faceAreas0() const;
435
436 //- Return access to the initial face centres
437 // Used for topology change
438 inline vectorField& faceCentres0() const;
439
440
441 // Transformations
442
443 //- Axis of rotation for rotational cyclic AMI
444 inline const vector& rotationAxis() const;
445
446 //- Point on axis of rotation for rotational cyclic AMI
447 inline const point& rotationCentre() const;
448
449 //- Translation vector for translational cyclic AMI
450 inline const vector& separationVector() const;
451
452 //- Transform patch-based positions from nbr side to this side
453 virtual void transformPosition(pointField&) const;
454
455 //- Transform a patch-based position from nbr side to this side
456 virtual void transformPosition
457 (
458 point& l,
459 const label facei
460 ) const;
461
462 //- Transform a patch-based position from this side to nbr side
463 virtual void reverseTransformPosition
464 (
465 point& l,
466 const label facei
467 ) const;
468
469 //- Transform a patch-based direction from this side to
470 // nbr side
471 virtual void reverseTransformDirection
472 (
473 vector& d,
474 const label facei
475 ) const;
476
477
478 // Interpolations
479
480 //- Interpolate field
481 template<class Type>
483 (
484 const Field<Type>& fld,
485 const UList<Type>& defaultValues = UList<Type>()
486 ) const;
487
488 //- Interpolate tmp field
489 template<class Type>
491 (
492 const tmp<Field<Type>>& tFld,
493 const UList<Type>& defaultValues = UList<Type>()
494 ) const;
495
496 //- Interpolate without periodic
497 template<class Type>
499 (
500 const Field<Type>& fld,
501 const UList<Type>& defaultValues
502 ) const;
503
504 //- Low-level interpolate List
505 template<class Type, class CombineOp>
506 void interpolate
507 (
508 const UList<Type>& fld,
509 const CombineOp& cop,
510 List<Type>& result,
511 const UList<Type>& defaultValues = UList<Type>()
512 ) const;
513
514
515 //- Calculate the patch geometry
516 virtual void calcGeometry
517 (
518 const primitivePatch& referPatch,
519 const pointField& thisCtrs,
520 const vectorField& thisAreas,
521 const pointField& thisCc,
522 const pointField& nbrCtrs,
523 const vectorField& nbrAreas,
524 const pointField& nbrCc
525 );
526
527 //- Initialize ordering for primitivePatch. Does not
528 //- refer to *this (except for name() and type() etc.)
529 virtual void initOrder
530 (
532 const primitivePatch&
533 ) const;
534
535 //- Return new ordering for primitivePatch.
536 // Ordering is -faceMap: for every face
537 // index of the new face -rotation:for every new face the clockwise
538 // shift of the original face. Return false if nothing changes
539 // (faceMap is identity, rotation is 0), true otherwise.
540 virtual bool order
541 (
543 const primitivePatch&,
545 labelList& rotation
546 ) const;
547
548 //- Return face index on neighbour patch which shares point p
549 //- following trajectory vector n
550 label pointFace
551 (
552 const label facei,
553 const vector& n,
554 point& p
555 ) const;
556
557 //- Write the polyPatch data as a dictionary
558 virtual void write(Ostream&) const;
559};
560
561
562// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
563
564} // End namespace Foam
565
566// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
567
568#include "cyclicAMIPolyPatchI.H"
569
570// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
571
572#ifdef NoRepository
574#endif
575
576// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
577
578#endif
579
580// ************************************************************************* //
label n
Info<< nl<< "Wrote faMesh in vtk format: "<< writer.output().name()<< nl;}{ vtk::lineWriter writer(aMesh.points(), aMesh.edges(), fileName(aMesh.mesh().time().globalPath()/"finiteArea-edges"));writer.writeGeometry();writer.beginCellData(4);writer.writeProcIDs();{ Field< scalar > fld(faMeshTools::flattenEdgeField(aMesh.magLe(), true))
Interpolation class dealing with transfer of data between two primitive patches with an arbitrary mes...
const labelListList & srcAddress() const
Return const access to source patch addressing.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
A list of faces which address into the list of points.
const Field< point_type > & points() const noexcept
Return reference to global points.
Buffers for inter-processor communications streams (UOPstream, UIPstream).
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:94
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
Encapsulates using "patchGroups" to specify coupled patch.
The coupledPolyPatch is an abstract base class for patches that couple regions of the computational d...
virtual transformType transform() const
Type of transform.
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
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 accross 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
const vector & separationVector() const
Translation vector for translational cyclic AMI.
virtual bool masterImplicit() const
Return implicit master.
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
virtual void reverseTransformDirection(vector &d, const label facei) const
Transform a patch-based direction from this side to.
tmp< Field< Type > > interpolateUntransformed(const Field< Type > &fld, const UList< Type > &defaultValues) const
Interpolate without periodic.
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.
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.
bool createAMIFaces() const
Return access to the createAMIFaces flag.
virtual void initUpdateMesh(PstreamBuffers &)
Initialise the update of the patch topology.
bool applyLowWeightCorrection() const
Return true if applying the low weight correction.
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:126
label index() const noexcept
The index of this patch in the boundaryMesh.
const word & name() const noexcept
The patch name.
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:364
const vectorField::subField faceCentres() const
Return face centres.
Definition: polyPatch.C:321
const labelUList & faceCells() const
Return face-cell addressing.
Definition: polyPatch.C:371
Direct mesh changes based on v1.3 polyTopoChange syntax.
A class for managing references or pointers (no reference counting)
Definition: refPtr.H:58
A class for managing temporary objects.
Definition: tmp.H:65
A class for handling words, derived from Foam::string.
Definition: word.H:68
volScalarField & p
Namespace for OpenFOAM.
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:56
UList< label > labelUList
A UList of labels.
Definition: UList.H:85
runTime write()
dictionary dict
static const char *const typeName
The type name used in ensight case files.
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition: typeInfo.H:73