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 | Copyright (C) 2011-2016 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"
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 class
110 
111  //- Flag to indicate that slave patch should be reversed for AMI
112  const bool AMIReverse_;
113 
114  //- Flag to indicate that patches should match/overlap
115  bool AMIRequireMatch_;
116 
117  //- Low weight correction threshold for AMI
118  const scalar AMILowWeightCorrection_;
119 
120  //- Projection surface
122 
123  //- Dictionary used during projection surface construction
124  const dictionary surfDict_;
125 
126 
127  // Protected Member Functions
128 
129  //- Reset the AMI interpolator
130  virtual void resetAMI
131  (
134  ) const;
135 
136  //- Recalculate the transformation tensors
137  virtual void calcTransforms();
138 
139  //- Initialise the calculation of the patch geometry
140  virtual void initGeometry(PstreamBuffers&);
141 
142  //- Calculate the patch geometry
143  virtual void calcGeometry(PstreamBuffers&);
144 
145  //- Initialise the patches for moving points
146  virtual void initMovePoints(PstreamBuffers& pBufs, const pointField&);
147 
148  //- Correct patches after moving points
149  virtual void movePoints(PstreamBuffers& pBufs, const pointField&);
150 
151  //- Initialise the update of the patch topology
152  virtual void initUpdateMesh(PstreamBuffers&);
153 
154  //- Update of the patch topology
155  virtual void updateMesh(PstreamBuffers&);
156 
157  //- Clear geometry
158  virtual void clearGeom();
159 
160 
161 public:
162 
163  //- Runtime type information
164  TypeName("cyclicAMI");
165 
166 
167  // Constructors
168 
169  //- Construct from (base couped patch) components
171  (
172  const word& name,
173  const label size,
174  const label start,
175  const label index,
176  const polyBoundaryMesh& bm,
177  const word& patchType,
179  );
180 
181  //- Construct from dictionary
183  (
184  const word& name,
185  const dictionary& dict,
186  const label index,
187  const polyBoundaryMesh& bm,
188  const word& patchType
189  );
190 
191  //- Construct as copy, resetting the boundary mesh
193 
194  //- Construct given the original patch and resetting the
195  // face list and boundary mesh information
197  (
198  const cyclicAMIPolyPatch& pp,
199  const polyBoundaryMesh& bm,
200  const label index,
201  const label newSize,
202  const label newStart,
203  const word& nbrPatchName
204  );
205 
206  //- Construct given the original patch and a map
208  (
209  const cyclicAMIPolyPatch& pp,
210  const polyBoundaryMesh& bm,
211  const label index,
212  const labelUList& mapAddressing,
213  const label newStart
214  );
215 
216 
217  //- Construct and return a clone, resetting the boundary mesh
218  virtual autoPtr<polyPatch> clone(const polyBoundaryMesh& bm) const
219  {
220  return autoPtr<polyPatch>(new cyclicAMIPolyPatch(*this, bm));
221  }
222 
223  //- Construct and return a clone, resetting the face list
224  // and boundary mesh
226  (
227  const polyBoundaryMesh& bm,
228  const label index,
229  const label newSize,
230  const label newStart
231  ) const
232  {
233  return autoPtr<polyPatch>
234  (
236  (
237  *this,
238  bm,
239  index,
240  newSize,
241  newStart,
242  nbrPatchName_
243  )
244  );
245  }
246 
247  //- Construct and return a clone, resetting the face list
248  // and boundary mesh
250  (
251  const polyBoundaryMesh& bm,
252  const label index,
253  const labelUList& mapAddressing,
254  const label newStart
255  ) const
256  {
257  return autoPtr<polyPatch>
258  (
260  (
261  *this,
262  bm,
263  index,
264  mapAddressing,
265  newStart
266  )
267  );
268  }
269 
270 
271  //- Destructor
272  virtual ~cyclicAMIPolyPatch();
273 
274 
275  // Member Functions
276 
277  // Access
278 
279  //- Is patch 'coupled'. Note that on AMI the geometry is not
280  // coupled but the fields are!
281  virtual bool coupled() const
282  {
283  return false;
284  }
285 
286  //- Neighbour patch name
287  inline const word& neighbPatchName() const;
288 
289  //- Neighbour patch ID
290  virtual label neighbPatchID() const;
291 
292  //- Does this side own the patch?
293  virtual bool owner() const;
294 
295  //- Return a reference to the neighbour patch
296  virtual const cyclicAMIPolyPatch& neighbPatch() const;
297 
298  //- Return a reference to the projection surface
299  const autoPtr<searchableSurface>& surfPtr() const;
300 
301  //- Return a reference to the AMI interpolator
302  const AMIPatchToPatchInterpolation& AMI() const;
303 
304  //- Return true if applying the low weight correction
305  bool applyLowWeightCorrection() const;
306 
307 
308  // Transformations
309 
310  //- Axis of rotation for rotational cyclic AMI
311  inline const vector& rotationAxis() const;
312 
313  //- Point on axis of rotation for rotational cyclic AMI
314  inline const point& rotationCentre() const;
315 
316  //- Translation vector for translational cyclic AMI
317  inline const vector& separationVector() const;
318 
319  //- Transform patch-based positions from nbr side to this side
320  virtual void transformPosition(pointField&) const;
321 
322  //- Transform a patch-based position from nbr side to this side
323  virtual void transformPosition
324  (
325  point& l,
326  const label facei
327  ) const;
328 
329  //- Transform a patch-based position from this side to nbr side
330  virtual void reverseTransformPosition
331  (
332  point& l,
333  const label facei
334  ) const;
335 
336  //- Transform a patch-based direction from this side to nbr side
337  virtual void reverseTransformDirection
338  (
339  vector& d,
340  const label facei
341  ) const;
342 
343 
344  // Interpolations
345 
346  //- Interpolate field
347  template<class Type>
349  (
350  const Field<Type>& fld,
351  const UList<Type>& defaultValues = UList<Type>()
352  ) const;
353 
354  //- Interpolate tmp field
355  template<class Type>
357  (
358  const tmp<Field<Type>>& tFld,
359  const UList<Type>& defaultValues = UList<Type>()
360  ) const;
361 
362  //- Low-level interpolate List
363  template<class Type, class CombineOp>
364  void interpolate
365  (
366  const UList<Type>& fld,
367  const CombineOp& cop,
368  List<Type>& result,
369  const UList<Type>& defaultValues = UList<Type>()
370  ) const;
371 
372 
373  //- Calculate the patch geometry
374  virtual void calcGeometry
375  (
376  const primitivePatch& referPatch,
377  const pointField& thisCtrs,
378  const vectorField& thisAreas,
379  const pointField& thisCc,
380  const pointField& nbrCtrs,
381  const vectorField& nbrAreas,
382  const pointField& nbrCc
383  );
384 
385  //- Initialize ordering for primitivePatch. Does not
386  // refer to *this (except for name() and type() etc.)
387  virtual void initOrder
388  (
390  const primitivePatch&
391  ) const;
392 
393  //- Return new ordering for primitivePatch.
394  // Ordering is -faceMap: for every face
395  // index of the new face -rotation:for every new face the clockwise
396  // shift of the original face. Return false if nothing changes
397  // (faceMap is identity, rotation is 0), true otherwise.
398  virtual bool order
399  (
401  const primitivePatch&,
403  labelList& rotation
404  ) const;
405 
406  //- Return face index on neighbour patch which shares point p
407  // following trajectory vector n
409  (
410  const label facei,
411  const vector& n,
412  point& p
413  ) const;
414 
415  //- Write the polyPatch data as a dictionary
416  virtual void write(Ostream&) const;
417 };
418 
419 
420 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
421 
422 } // End namespace Foam
423 
424 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
425 
426 #include "cyclicAMIPolyPatchI.H"
427 
428 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
429 
430 #ifdef NoRepository
432 #endif
433 
434 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
435 
436 #endif
437 
438 // ************************************************************************* //
virtual void clearGeom()
Clear geometry.
vector separationVector_
Translation vector.
dictionary dict
const AMIPatchToPatchInterpolation & AMI() const
Return a reference to the AMI interpolator.
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.
bool AMIRequireMatch_
Flag to indicate that patches should match/overlap.
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:137
virtual bool coupled() const
Is patch &#39;coupled&#39;. Note that on AMI the geometry is not.
virtual void write(Ostream &) const
Write the polyPatch data as a dictionary.
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.
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.
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
interpolationMethod
Enumeration specifying interpolation method.
const dictionary surfDict_
Dictionary used during projection surface construction.
virtual void calcGeometry(PstreamBuffers &)
Calculate the patch geometry.
label nbrPatchID_
Index of other half.
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.
cyclicAMIPolyPatch(const word &name, const label size, const label start, const label index, const polyBoundaryMesh &bm, const word &patchType, const transformType transform=UNKNOWN)
Construct from (base couped patch) components.
const bool AMIReverse_
Flag to indicate that slave patch should be reversed for AMI.
const word & neighbPatchName() const
Neighbour patch name.
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.
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.
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.
virtual void resetAMI(const AMIPatchToPatchInterpolation::interpolationMethod &AMIMethod=AMIPatchToPatchInterpolation::imFaceAreaWeight) const
Reset the AMI interpolator.
Cyclic patch for Arbitrary Mesh Interface (AMI)
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:61
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.
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:300
Interpolation class dealing with transfer of data between two primitive patches with an arbitrary mes...
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.
label pointFace(const label facei, const vector &n, point &p) const
Return face index on neighbour patch which shares point p.
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.