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