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-2020 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 "cyclicTransform.H"
40 #include "AMIInterpolation.H"
41 #include "polyBoundaryMesh.H"
42 #include "coupleGroupIdentifier.H"
43 
44 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
45 
46 namespace Foam
47 {
48 
49 /*---------------------------------------------------------------------------*\
50  Class cyclicAMIPolyPatch Declaration
51 \*---------------------------------------------------------------------------*/
52 
54 :
55  public coupledPolyPatch,
56  public cyclicTransform
57 {
58 protected:
59 
60  // Protected data
61 
62  //- Name of cyclic neighbour patch
63  mutable word nbrPatchName_;
64 
65  //- Optional patchGroup to find nbrPatch
67 
68  //- Index of cyclic neighbour patch
69  mutable label nbrPatchID_;
70 
71  //- AMI interpolation classes
73 
74  //- AMI transforms (from source to target)
76 
77  //- Flag to indicate that slave patch should be reversed for AMI
78  const bool AMIReverse_;
79 
80  //- Flag to indicate that patches should match/overlap
81  const bool AMIRequireMatch_;
82 
83  //- Low weight correction threshold for AMI
84  const scalar AMILowWeightCorrection_;
85 
86  //- AMI Method
88 
89  //- Projection surface
91 
92  //- Dictionary used during projection surface construction
93  const dictionary surfDict_;
94 
95 
96  // Protected Member Functions
97 
98  //- Reset the AMI interpolator
99  virtual void resetAMI() const;
100 
101  //- Initialise the calculation of the patch geometry
102  virtual void initCalcGeometry(PstreamBuffers&);
103 
104  //- Calculate the patch geometry
105  virtual void calcGeometry(PstreamBuffers&);
106 
107  //- Initialise the patches for moving points
108  virtual void initMovePoints(PstreamBuffers& pBufs, const pointField&);
109 
110  //- Correct patches after moving points
111  virtual void movePoints(PstreamBuffers& pBufs, const pointField&);
112 
113  //- Initialise the update of the patch topology
114  virtual void initUpdateMesh(PstreamBuffers&);
115 
116  //- Update of the patch topology
117  virtual void updateMesh(PstreamBuffers&);
118 
119  //- Clear geometry
120  virtual void clearGeom();
121 
122 
123 public:
124 
125  //- Runtime type information
126  TypeName("cyclicAMI");
127 
128 
129  // Constructors
130 
131  //- Construct from (base couped patch) components
133  (
134  const word& name,
135  const label size,
136  const label start,
137  const label index,
138  const polyBoundaryMesh& bm,
139  const word& patchType,
140  const bool AMIRequireMatch = true,
143  );
144 
145  //- Construct from dictionary
147  (
148  const word& name,
149  const dictionary& dict,
150  const label index,
151  const polyBoundaryMesh& bm,
152  const word& patchType,
153  const bool AMIRequireMatch = true,
156  );
157 
158  //- Construct as copy, resetting the boundary mesh
160 
161  //- Construct given the original patch and resetting the
162  // face list and boundary mesh information
164  (
165  const cyclicAMIPolyPatch& pp,
166  const polyBoundaryMesh& bm,
167  const label index,
168  const label newSize,
169  const label newStart,
170  const word& nbrPatchName
171  );
172 
173  //- Construct given the original patch and a map
175  (
176  const cyclicAMIPolyPatch& pp,
177  const polyBoundaryMesh& bm,
178  const label index,
179  const labelUList& mapAddressing,
180  const label newStart
181  );
182 
183 
184  //- Construct and return a clone, resetting the boundary mesh
185  virtual autoPtr<polyPatch> clone(const polyBoundaryMesh& bm) const
186  {
187  return autoPtr<polyPatch>(new cyclicAMIPolyPatch(*this, bm));
188  }
189 
190  //- Construct and return a clone, resetting the face list
191  // and boundary mesh
193  (
194  const polyBoundaryMesh& bm,
195  const label index,
196  const label newSize,
197  const label newStart
198  ) const
199  {
200  return autoPtr<polyPatch>
201  (
203  (
204  *this,
205  bm,
206  index,
207  newSize,
208  newStart,
209  nbrPatchName_
210  )
211  );
212  }
213 
214  //- Construct and return a clone, resetting the face list
215  // and boundary mesh
217  (
218  const polyBoundaryMesh& bm,
219  const label index,
220  const labelUList& mapAddressing,
221  const label newStart
222  ) const
223  {
224  return autoPtr<polyPatch>
225  (
227  (
228  *this,
229  bm,
230  index,
231  mapAddressing,
232  newStart
233  )
234  );
235  }
236 
237 
238  //- Destructor
239  virtual ~cyclicAMIPolyPatch();
240 
241 
242  // Member Functions
243 
244  // Access
245 
246  //- Is patch 'coupled'. Note that on AMI the geometry is not
247  // coupled but the fields are!
248  virtual bool coupled() const
249  {
250  return false;
251  }
252 
253  //- Neighbour patch name
254  inline const word& nbrPatchName() const;
255 
256  //- Neighbour patch ID
257  virtual label nbrPatchID() const;
258 
259  //- Does this side own the patch?
260  virtual bool owner() const;
261 
262  //- Return a reference to the neighbour patch
263  virtual const cyclicAMIPolyPatch& nbrPatch() const;
264 
265  //- Return a reference to the projection surface
266  const autoPtr<searchableSurface>& surfPtr() const;
267 
268  //- Return a reference to the AMI interpolators
269  const PtrList<AMIInterpolation>& AMIs() const;
270 
271  //- Return a reference to the AMI transforms
272  const List<transformer>& AMITransforms() const;
273 
274  //- Return true if applying the low weight correction
275  bool applyLowWeightCorrection() const;
276 
277  //- Return the weights sum for this patch
278  virtual const scalarField& weightsSum() const;
279 
280  //- Return the weights sum for the neighbour patch
281  virtual const scalarField& nbrWeightsSum() const;
282 
283 
284  // Transformations
285 
286  //- Return transformation between the coupled patches
287  virtual const transformer& transform() const
288  {
290  }
291 
292 
293  // Interpolations
294 
295  //- Interpolate field
296  template<class Type>
298  (
299  const Field<Type>& fld,
300  const UList<Type>& defaultValues = UList<Type>()
301  ) const;
302 
303  //- Interpolate tmp field
304  template<class Type>
306  (
307  const tmp<Field<Type>>& tFld,
308  const UList<Type>& defaultValues = UList<Type>()
309  ) const;
310 
311  //- Interpolate field component
313  (
314  const scalarField& field,
315  const direction cmpt,
316  const direction rank,
317  const scalarUList& defaultValues = scalarUList()
318  ) const;
319 
320 
321  //- Initialize ordering for primitivePatch. Does not
322  // refer to *this (except for name() and type() etc.)
323  virtual void initOrder
324  (
326  const primitivePatch&
327  ) const;
328 
329  //- Return new ordering for primitivePatch.
330  // Ordering is -faceMap: for every face
331  // index of the new face -rotation:for every new face the clockwise
332  // shift of the original face. Return false if nothing changes
333  // (faceMap is identity, rotation is 0), true otherwise.
334  virtual bool order
335  (
337  const primitivePatch&,
339  labelList& rotation
340  ) const;
341 
342  //- Return the transform and face indices on neighbour patch which
343  // shares point p following trajectory vector n
345  (
346  const label facei,
347  const vector& n,
348  point& p
349  ) const;
350 
351  //- Index of processor that holds all of both sides, or -1 if
352  // distributed
353  label singlePatchProc() const;
354 
355  //- Write the polyPatch data as a dictionary
356  virtual void write(Ostream&) const;
357 };
358 
359 
360 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
361 
362 } // End namespace Foam
363 
364 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
365 
366 #include "cyclicAMIPolyPatchI.H"
367 
368 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
369 
370 #ifdef NoRepository
372 #endif
373 
374 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
375 
376 #endif
377 
378 // ************************************************************************* //
virtual void clearGeom()
Clear geometry.
virtual void resetAMI() const
Reset the AMI interpolator.
virtual void initCalcGeometry(PstreamBuffers &)
Initialise the calculation of the patch geometry.
dictionary dict
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.
Vector-tensor class used to perform translations and rotations in 3D space.
Definition: transformer.H:83
virtual void initUpdateMesh(PstreamBuffers &)
Initialise the update of the patch topology.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
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
labelPair pointAMIAndFace(const label facei, const vector &n, point &p) const
Return the transform and face indices on neighbour patch which.
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...
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 cyclic neighbour patch.
virtual const scalarField & nbrWeightsSum() const
Return the weights sum for the neighbour patch.
const List< transformer > & AMITransforms() const
Return a reference to the AMI transforms.
virtual label nbrPatchID() const
Neighbour patch ID.
virtual const cyclicAMIPolyPatch & nbrPatch() const
Return a reference to the neighbour patch.
const coupleGroupIdentifier coupleGroup_
Optional patchGroup to find nbrPatch.
cyclicAMIPolyPatch(const word &name, const label size, const label start, const label index, const polyBoundaryMesh &bm, const word &patchType, const bool AMIRequireMatch=true, const AMIInterpolation::interpolationMethod AMIMethod=AMIInterpolation::imFaceAreaWeight)
Construct from (base couped patch) components.
A list of faces which address into the list of points.
const bool AMIReverse_
Flag to indicate that slave patch should be reversed for AMI.
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 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.
virtual void updateMesh(PstreamBuffers &)
Update of the patch topology.
virtual const transformer & transform() const
Return transformation between the coupled patches.
PtrList< AMIInterpolation > AMIs_
AMI interpolation classes.
Cyclic patch for Arbitrary Mesh Interface (AMI)
const transformer & transform() const
Return transformation between the coupled patches.
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:54
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
const word & nbrPatchName() const
Neighbour patch name.
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
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
Cyclic plane tranformation.
rDeltaTY field()
volScalarField & p
A class for managing temporary objects.
Definition: PtrList.H:53
word nbrPatchName_
Name of cyclic neighbour patch.
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.
List< transformer > AMITransforms_
AMI transforms (from source to target)