cyclicAMIPolyPatch.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2018 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 Class
25  Foam::cyclicAMIPolyPatch
26 
27 Description
28  Cyclic patch for Arbitrary Mesh Interface (AMI)
29 
30 SourceFiles
31  cyclicAMIPolyPatch.C
32 
33 \*---------------------------------------------------------------------------*/
34 
35 #ifndef cyclicAMIPolyPatch_H
36 #define cyclicAMIPolyPatch_H
37 
38 #include "coupledPolyPatch.H"
39 #include "AMIInterpolation.H"
40 #include "polyBoundaryMesh.H"
41 #include "coupleGroupIdentifier.H"
42 
43 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
44 
45 namespace Foam
46 {
47 
48 /*---------------------------------------------------------------------------*\
49  Class cyclicAMIPolyPatch Declaration
50 \*---------------------------------------------------------------------------*/
51 
53 :
54  public coupledPolyPatch
55 {
56  // Private Member Functions
57 
58  //- Return normal of face at max distance from rotation axis
59  vector findFaceNormalMaxRadius(const pointField& faceCentres) const;
60 
61  void calcTransforms
62  (
63  const primitivePatch& half0,
64  const pointField& half0Ctrs,
65  const vectorField& half0Areas,
66  const pointField& half1Ctrs,
67  const vectorField& half1Areas
68  );
69 
70 
71 protected:
72 
73  // Protected data
74 
75  //- Name of other half
76  mutable word nbrPatchName_;
77 
78  //- Optional patchGroup to find neighbPatch
80 
81  //- Index of other half
82  mutable label nbrPatchID_;
83 
84 
85  // Transformations
86 
87  // For rotation
88 
89  //- Axis of rotation for rotational cyclics
91 
92  //- Point on axis of rotation for rotational cyclics
94 
95  //- Flag to show whether the rotation angle is defined
97 
98  //- Rotation angle
99  scalar rotationAngle_;
100 
101 
102  // For translation
103 
104  //- Translation vector
106 
107 
108  //- AMI interpolation classes
110 
111  //- AMI transforms (from source to target)
113 
114  //- Flag to indicate that slave patch should be reversed for AMI
115  const bool AMIReverse_;
116 
117  //- Flag to indicate that patches should match/overlap
118  const bool AMIRequireMatch_;
119 
120  //- Low weight correction threshold for AMI
121  const scalar AMILowWeightCorrection_;
122 
123  //- AMI Method
125 
126  //- Projection surface
128 
129  //- Dictionary used during projection surface construction
130  const dictionary surfDict_;
131 
132 
133  // Protected Member Functions
134 
135  //- Reset the AMI interpolator
136  virtual void resetAMI() const;
137 
138  //- Recalculate the transformation tensors
139  virtual void calcTransforms();
140 
141  //- Initialise the calculation of the patch geometry
142  virtual void initGeometry(PstreamBuffers&);
143 
144  //- Calculate the patch geometry
145  virtual void calcGeometry(PstreamBuffers&);
146 
147  //- Initialise the patches for moving points
148  virtual void initMovePoints(PstreamBuffers& pBufs, const pointField&);
149 
150  //- Correct patches after moving points
151  virtual void movePoints(PstreamBuffers& pBufs, const pointField&);
152 
153  //- Initialise the update of the patch topology
154  virtual void initUpdateMesh(PstreamBuffers&);
155 
156  //- Update of the patch topology
157  virtual void updateMesh(PstreamBuffers&);
158 
159  //- Clear geometry
160  virtual void clearGeom();
161 
162 
163 public:
164 
165  //- Runtime type information
166  TypeName("cyclicAMI");
167 
168 
169  // Constructors
170 
171  //- Construct from (base couped patch) components
173  (
174  const word& name,
175  const label size,
176  const label start,
177  const label index,
178  const polyBoundaryMesh& bm,
179  const word& patchType,
181  const bool AMIRequireMatch = true,
184  );
185 
186  //- Construct from dictionary
188  (
189  const word& name,
190  const dictionary& dict,
191  const label index,
192  const polyBoundaryMesh& bm,
193  const word& patchType,
194  const bool AMIRequireMatch = true,
197  );
198 
199  //- Construct as copy, resetting the boundary mesh
201 
202  //- Construct given the original patch and resetting the
203  // face list and boundary mesh information
205  (
206  const cyclicAMIPolyPatch& pp,
207  const polyBoundaryMesh& bm,
208  const label index,
209  const label newSize,
210  const label newStart,
211  const word& nbrPatchName
212  );
213 
214  //- Construct given the original patch and a map
216  (
217  const cyclicAMIPolyPatch& pp,
218  const polyBoundaryMesh& bm,
219  const label index,
220  const labelUList& mapAddressing,
221  const label newStart
222  );
223 
224 
225  //- Construct and return a clone, resetting the boundary mesh
226  virtual autoPtr<polyPatch> clone(const polyBoundaryMesh& bm) const
227  {
228  return autoPtr<polyPatch>(new cyclicAMIPolyPatch(*this, bm));
229  }
230 
231  //- Construct and return a clone, resetting the face list
232  // and boundary mesh
234  (
235  const polyBoundaryMesh& bm,
236  const label index,
237  const label newSize,
238  const label newStart
239  ) const
240  {
241  return autoPtr<polyPatch>
242  (
244  (
245  *this,
246  bm,
247  index,
248  newSize,
249  newStart,
250  nbrPatchName_
251  )
252  );
253  }
254 
255  //- Construct and return a clone, resetting the face list
256  // and boundary mesh
258  (
259  const polyBoundaryMesh& bm,
260  const label index,
261  const labelUList& mapAddressing,
262  const label newStart
263  ) const
264  {
265  return autoPtr<polyPatch>
266  (
268  (
269  *this,
270  bm,
271  index,
272  mapAddressing,
273  newStart
274  )
275  );
276  }
277 
278 
279  //- Destructor
280  virtual ~cyclicAMIPolyPatch();
281 
282 
283  // Member Functions
284 
285  // Access
286 
287  //- Is patch 'coupled'. Note that on AMI the geometry is not
288  // coupled but the fields are!
289  virtual bool coupled() const
290  {
291  return false;
292  }
293 
294  //- Neighbour patch name
295  inline const word& neighbPatchName() const;
296 
297  //- Neighbour patch ID
298  virtual label neighbPatchID() const;
299 
300  //- Does this side own the patch?
301  virtual bool owner() const;
302 
303  //- Return a reference to the neighbour patch
304  virtual const cyclicAMIPolyPatch& neighbPatch() const;
305 
306  //- Return a reference to the projection surface
307  const autoPtr<searchableSurface>& surfPtr() const;
308 
309  //- Return a reference to the AMI interpolators
310  const PtrList<AMIInterpolation>& AMIs() const;
311 
312  //- Return a reference to the AMI transforms
314 
315  //- Return true if applying the low weight correction
316  bool applyLowWeightCorrection() const;
317 
318  //- Return the weights sum for this patch
319  virtual const scalarField& weightsSum() const;
320 
321  //- Return the weights sum for the neighbour patch
322  virtual const scalarField& neighbWeightsSum() const;
323 
324 
325 
326  // Transformations
327 
328  //- Axis of rotation for rotational cyclic AMI
329  inline const vector& rotationAxis() const;
330 
331  //- Point on axis of rotation for rotational cyclic AMI
332  inline const point& rotationCentre() const;
333 
334  //- Translation vector for translational cyclic AMI
335  inline const vector& separationVector() const;
336 
337  //- Transform patch-based positions from nbr side to this side
338  virtual void transformPosition(pointField&) const;
339 
340  //- Transform a patch-based position from nbr side to this side
341  virtual void transformPosition
342  (
343  point& l,
344  const label facei
345  ) const;
346 
347  //- Transform a patch-based direction from nbr side to this side
348  virtual void transformDirection
349  (
350  vector& d,
351  const label facei
352  ) const;
353 
354  //- Transform a patch-based position from this side to nbr side
355  virtual void reverseTransformPosition
356  (
357  point& l,
358  const label facei
359  ) const;
360 
361  //- Transform a patch-based direction from this side to nbr side
362  virtual void reverseTransformDirection
363  (
364  vector& d,
365  const label facei
366  ) const;
367 
368 
369  // Interpolations
370 
371  //- Interpolate field
372  template<class Type>
374  (
375  const Field<Type>& fld,
376  const UList<Type>& defaultValues = UList<Type>()
377  ) const;
378 
379  //- Interpolate tmp field
380  template<class Type>
382  (
383  const tmp<Field<Type>>& tFld,
384  const UList<Type>& defaultValues = UList<Type>()
385  ) const;
386 
387  //- Interpolate field component
389  (
390  const scalarField& field,
391  const direction cmpt,
392  const direction rank,
393  const scalarUList& defaultValues = scalarUList()
394  ) const;
395 
396 
397  //- Calculate the patch geometry
398  virtual void calcGeometry
399  (
400  const primitivePatch& referPatch,
401  const pointField& thisCtrs,
402  const vectorField& thisAreas,
403  const pointField& thisCc,
404  const pointField& nbrCtrs,
405  const vectorField& nbrAreas,
406  const pointField& nbrCc
407  );
408 
409  //- Initialize ordering for primitivePatch. Does not
410  // refer to *this (except for name() and type() etc.)
411  virtual void initOrder
412  (
414  const primitivePatch&
415  ) const;
416 
417  //- Return new ordering for primitivePatch.
418  // Ordering is -faceMap: for every face
419  // index of the new face -rotation:for every new face the clockwise
420  // shift of the original face. Return false if nothing changes
421  // (faceMap is identity, rotation is 0), true otherwise.
422  virtual bool order
423  (
425  const primitivePatch&,
427  labelList& rotation
428  ) const;
429 
430  //- Return the transform and face indices on neighbour patch which
431  // shares point p following trajectory vector n
433  (
434  const label facei,
435  const vector& n,
436  point& p
437  ) const;
438 
439  //- Index of processor that holds all of both sides, or -1 if
440  // distributed
441  label singlePatchProc() const;
442 
443  //- Write the polyPatch data as a dictionary
444  virtual void write(Ostream&) const;
445 };
446 
447 
448 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
449 
450 } // End namespace Foam
451 
452 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
453 
454 #include "cyclicAMIPolyPatchI.H"
455 
456 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
457 
458 #ifdef NoRepository
460 #endif
461 
462 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
463 
464 #endif
465 
466 // ************************************************************************* //
virtual void clearGeom()
Clear geometry.
virtual void resetAMI() const
Reset the AMI interpolator.
vector separationVector_
Translation vector.
dictionary dict
virtual label neighbPatchID() const
Neighbour patch ID.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
const word & name() const
Return name.
virtual bool order(PstreamBuffers &, const primitivePatch &, labelList &faceMap, labelList &rotation) const
Return new ordering for primitivePatch.
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 any number of values (e...
Definition: dictionary.H:158
List< vectorTensorTransform > AMITransforms_
AMI transforms (from source to target)
virtual bool coupled() const
Is patch &#39;coupled&#39;. Note that on AMI the geometry is not.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
virtual void write(Ostream &) const
Write the polyPatch data as a dictionary.
uint8_t direction
Definition: direction.H:45
virtual void calcTransforms()
Recalculate the transformation tensors.
virtual const scalarField & neighbWeightsSum() const
Return the weights sum for the neighbour patch.
labelPair pointAMIAndFace(const label facei, const vector &n, point &p) const
Return the transform and face indices on neighbour patch which.
cyclicAMIPolyPatch(const word &name, const label size, const label start, const label index, const polyBoundaryMesh &bm, const word &patchType, const transformType transform=UNKNOWN, const bool AMIRequireMatch=true, const AMIInterpolation::interpolationMethod AMIMethod=AMIInterpolation::imFaceAreaWeight)
Construct from (base couped patch) components.
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.
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.
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.
virtual autoPtr< PrimitivePatch< FaceList, PointField > > clone() const
Construct and return a clone.
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
interpolationMethod
Enumeration specifying interpolation method.
const bool AMIRequireMatch_
Flag to indicate that patches should match/overlap.
const dictionary surfDict_
Dictionary used during projection surface construction.
virtual void calcGeometry(PstreamBuffers &)
Calculate the patch geometry.
label nbrPatchID_
Index of other half.
const List< vectorTensorTransform > & AMITransforms() const
Return a reference to the AMI transforms.
const vector & separationVector() const
Translation vector for translational cyclic AMI.
const coupleGroupIdentifier coupleGroup_
Optional patchGroup to find neighbPatch.
A list of faces which address into the list of points.
const vector & rotationAxis() const
Axis of rotation for rotational cyclic AMI.
const bool AMIReverse_
Flag to indicate that slave patch should be reversed for AMI.
const word & neighbPatchName() const
Neighbour patch name.
An ordered pair of two objects of type <T> with first() and second() elements.
Definition: contiguous.H:49
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< ' ';}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< ' ';}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< ' ';}gmvFile<< nl;forAll(lagrangianScalarNames, i){ const word &name=lagrangianScalarNames[i];IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
const point & rotationCentre() const
Point on axis of rotation for rotational cyclic AMI.
const AMIInterpolation::interpolationMethod AMIMethod_
AMI Method.
const PtrList< AMIInterpolation > & AMIs() const
Return a reference to the AMI interpolators.
A class for handling words, derived from string.
Definition: word.H:59
const autoPtr< searchableSurface > & surfPtr() const
Return a reference to the projection surface.
const scalar AMILowWeightCorrection_
Low weight correction threshold for AMI.
label singlePatchProc() const
Index of processor that holds all of both sides, or -1 if.
point rotationCentre_
Point on axis of rotation for rotational cyclics.
virtual void updateMesh(PstreamBuffers &)
Update of the patch topology.
PtrList< AMIInterpolation > AMIs_
AMI interpolation classes.
Cyclic patch for Arbitrary Mesh Interface (AMI)
virtual const scalarField & weightsSum() const
Return the weights sum for this patch.
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:60
Base class for Arbitrary Mesh Interface (AMI) methods.
Definition: AMIMethod.H:55
virtual ~cyclicAMIPolyPatch()
Destructor.
Foam::polyBoundaryMesh.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
Buffers for inter-processor communications streams (UOPstream, UIPstream).
TypeName("cyclicAMI")
Runtime type information.
virtual void initOrder(PstreamBuffers &, const primitivePatch &) const
Initialize ordering for primitivePatch. Does not.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:70
virtual const cyclicAMIPolyPatch & neighbPatch() const
Return a reference to the neighbour patch.
virtual bool owner() const
Does this side own the patch?
label index() const
Return the index of this patch in the boundaryMesh.
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:303
UList< scalar > scalarUList
A UList of scalars.
Definition: scalarList.H:48
label n
virtual transformType transform() const
Type of transform.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
volScalarField & p
A class for managing temporary objects.
Definition: PtrList.H:53
word nbrPatchName_
Name of other half.
virtual void reverseTransformPosition(point &l, const label facei) const
Transform a patch-based position from this side to nbr side.
scalar rotationAngle_
Rotation angle.
Namespace for OpenFOAM.
virtual void transformDirection(vector &d, const label facei) const
Transform a patch-based direction from nbr side to this side.
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.