cyclicAMIFvPatch.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) 2019 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
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 
27 Class
28  Foam::cyclicAMIFvPatch
29 
30 Description
31  Cyclic patch for Arbitrary Mesh Interface (AMI)
32 
33 SourceFiles
34  cyclicAMIFvPatch.C
35 
36 \*---------------------------------------------------------------------------*/
37 
38 #ifndef cyclicAMIFvPatch_H
39 #define cyclicAMIFvPatch_H
40 
41 #include "coupledFvPatch.H"
42 #include "cyclicAMILduInterface.H"
43 #include "cyclicAMIPolyPatch.H"
44 #include "fvBoundaryMesh.H"
45 
46 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
47 
48 namespace Foam
49 {
50 
51 /*---------------------------------------------------------------------------*\
52  Class cyclicAMIFvPatch Declaration
53 \*---------------------------------------------------------------------------*/
54 
55 class cyclicAMIFvPatch
56 :
57  public coupledFvPatch,
59 {
60  // Private data
61 
62  const cyclicAMIPolyPatch& cyclicAMIPolyPatch_;
63 
64 
65 protected:
66 
67  // Protected Member functions
68 
69  //- Make patch weighting factors
70  void makeWeights(scalarField&) const;
71 
72  //- Correct patch deltaCoeffs
73  virtual void makeDeltaCoeffs(scalarField&) const;
74 
75  //- Correct patch non-ortho deltaCoeffs
76  virtual void makeNonOrthoDeltaCoeffs(scalarField&) const;
77 
78  //- Correct patch non-ortho correction vectors
79  virtual void makeNonOrthoCorrVectors(vectorField&) const;
80 
81  //- Correct patches after moving points
82  virtual void movePoints();
83 
84 
85 public:
86 
87  //- Runtime type information
88  TypeName(cyclicAMIPolyPatch::typeName_());
89 
90 
91  // Constructors
92 
93  //- Construct from polyPatch
95  :
96  coupledFvPatch(patch, bm),
98  cyclicAMIPolyPatch_(refCast<const cyclicAMIPolyPatch>(patch))
99  {}
100 
101 
102  // Member functions
103 
104  // Access
105 
106 // // Implicit treatment functions
107 //
108 // //- Return number of new internal of this polyPatch faces
109 // virtual void newInternalProcFaces(label&, label&) const;
110 //
111 // //- Return nbrCells
112 // virtual const labelUList& nbrCells() const
113 // {
114 // return cyclicAMIPolyPatch_.neighbPatch().faceCells();
115 // }
116 //
117 // //- Return nbr patch ID
118 // virtual label neighbPolyPatchID() const
119 // {
120 // return neighbPatchID();
121 // }
122 //
123 // //- Return collocated faces map
124 // virtual refPtr<labelListList> mapCollocatedFaces() const
125 // {
126 // const labelListList& sourceFaces =
127 // cyclicAMIPolyPatch_.AMI().srcAddress();
128 // return refPtr<labelListList>
129 // (
130 // new labelListList(sourceFaces)
131 // );
132 // }
133 //
134 // //- Return implicit master
135 // virtual bool masterImplicit() const
136 // {
137 // return owner();
138 // }
139 
140 
141  //- Return local reference cast into the cyclic patch
142  const cyclicAMIPolyPatch& cyclicAMIPatch() const
143  {
144  return cyclicAMIPolyPatch_;
145  }
146 
147  //- Return neighbour
148  virtual label neighbPatchID() const
149  {
150  return cyclicAMIPolyPatch_.neighbPatchID();
151  }
152 
153  virtual bool owner() const
154  {
155  return cyclicAMIPolyPatch_.owner();
156  }
157 
158  //- Return processor number
159  virtual const cyclicAMIFvPatch& neighbPatch() const
160  {
161  return refCast<const cyclicAMIFvPatch>
162  (
163  this->boundaryMesh()[cyclicAMIPolyPatch_.neighbPatchID()]
164  );
165  }
166 
167  //- Return a reference to the AMI interpolator
168  virtual const AMIPatchToPatchInterpolation& AMI() const
169  {
170  return cyclicAMIPolyPatch_.AMI();
171  }
172 
173  //- Return true if applying the low weight correction
174  virtual bool applyLowWeightCorrection() const
175  {
176  return cyclicAMIPolyPatch_.applyLowWeightCorrection();
177  }
178 
179 
180  //- Are the cyclic planes parallel
181  virtual bool parallel() const
182  {
183  return cyclicAMIPolyPatch_.parallel();
184  }
185 
186  //- Return face transformation tensor
187  virtual const tensorField& forwardT() const
188  {
189  return cyclicAMIPolyPatch_.forwardT();
190  }
191 
192  //- Return neighbour-cell transformation tensor
193  virtual const tensorField& reverseT() const
194  {
195  return cyclicAMIPolyPatch_.reverseT();
196  }
197 
198  const cyclicAMIFvPatch& neighbFvPatch() const
199  {
200  return refCast<const cyclicAMIFvPatch>
201  (
202  this->boundaryMesh()[cyclicAMIPolyPatch_.neighbPatchID()]
203  );
204  }
205 
206 
207  //- Return true if this patch is coupled. This is equivalent
208  //- to the coupledPolyPatch::coupled() if parallel running or
209  //- both sides present, false otherwise
210  virtual bool coupled() const;
211 
212  //- Return delta (P to N) vectors across coupled patch
213  virtual tmp<vectorField> delta() const;
214 
215  template<class Type>
217  (
218  const Field<Type>& fld,
219  const UList<Type>& defaultValues = UList<Type>()
220  ) const
221  {
222  return cyclicAMIPolyPatch_.interpolate(fld, defaultValues);
223  }
224 
225  template<class Type>
227  (
228  const tmp<Field<Type>>& tFld,
229  const UList<Type>& defaultValues = UList<Type>()
230  ) const
231  {
232  return cyclicAMIPolyPatch_.interpolate(tFld, defaultValues);
233  }
234 
235 
236  // Interface transfer functions
237 
238  //- Return the values of the given internal data adjacent to
239  //- the interface as a field
241  (
242  const labelUList& internalData
243  ) const;
244 
245 
246  //- Return the values of the given internal data adjacent to
247  //- the interface as a field using a mapping faceCell
249  (
250  const labelUList& internalData,
251  const labelUList& faceCells
252  ) const;
253 
254 
255  //- Return neighbour field
257  (
258  const Pstream::commsTypes commsType,
259  const labelUList& internalData
260  ) const;
261 
262 };
263 
264 
265 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
266 
267 } // End namespace Foam
268 
269 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
270 
271 #endif
272 
273 // ************************************************************************* //
Foam::cyclicAMIFvPatch::neighbPatchID
virtual label neighbPatchID() const
Return neighbour.
Definition: cyclicAMIFvPatch.H:147
Foam::cyclicAMIFvPatch::interpolate
tmp< Field< Type > > interpolate(const Field< Type > &fld, const UList< Type > &defaultValues=UList< Type >()) const
Definition: cyclicAMIFvPatch.H:216
Foam::cyclicAMIFvPatch::internalFieldTransfer
virtual tmp< labelField > internalFieldTransfer(const Pstream::commsTypes commsType, const labelUList &internalData) const
Return neighbour field.
Definition: cyclicAMIFvPatch.C:236
Foam::cyclicAMIPolyPatch::owner
virtual bool owner() const
Does this side own the patch?
Definition: cyclicAMIPolyPatch.C:964
Foam::cyclicAMIFvPatch::AMI
virtual const AMIPatchToPatchInterpolation & AMI() const
Return a reference to the AMI interpolator.
Definition: cyclicAMIFvPatch.H:167
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::cyclicAMIFvPatch::delta
virtual tmp< vectorField > delta() const
Return delta (P to N) vectors across coupled patch.
Definition: cyclicAMIFvPatch.C:157
Foam::cyclicAMILduInterface
An abstract base class for cyclic AMI coupled interfaces.
Definition: cyclicAMILduInterface.H:52
Foam::coupledPolyPatch::forwardT
virtual const tensorField & forwardT() const
Return face transformation tensor.
Definition: coupledPolyPatch.H:301
Foam::cyclicAMIFvPatch::makeNonOrthoDeltaCoeffs
virtual void makeNonOrthoDeltaCoeffs(scalarField &) const
Correct patch non-ortho deltaCoeffs.
Definition: cyclicAMIFvPatch.C:143
Foam::coupledPolyPatch::reverseT
virtual const tensorField & reverseT() const
Return neighbour-cell transformation tensor.
Definition: coupledPolyPatch.H:307
Foam::cyclicAMIFvPatch::neighbPatch
virtual const cyclicAMIFvPatch & neighbPatch() const
Return processor number.
Definition: cyclicAMIFvPatch.H:158
Foam::cyclicAMIFvPatch::neighbFvPatch
const cyclicAMIFvPatch & neighbFvPatch() const
Definition: cyclicAMIFvPatch.H:197
Foam::cyclicAMIFvPatch::applyLowWeightCorrection
virtual bool applyLowWeightCorrection() const
Return true if applying the low weight correction.
Definition: cyclicAMIFvPatch.H:173
Foam::cyclicAMIFvPatch::makeDeltaCoeffs
virtual void makeDeltaCoeffs(scalarField &) const
Correct patch deltaCoeffs.
Definition: cyclicAMIFvPatch.C:137
Foam::cyclicAMIFvPatch::makeNonOrthoCorrVectors
virtual void makeNonOrthoCorrVectors(vectorField &) const
Correct patch non-ortho correction vectors.
Definition: cyclicAMIFvPatch.C:150
coupledFvPatch.H
Foam::fvPatch::boundaryMesh
const fvBoundaryMesh & boundaryMesh() const
Return boundaryMesh reference.
Definition: fvPatch.H:203
Foam::fvBoundaryMesh
Foam::fvBoundaryMesh.
Definition: fvBoundaryMesh.H:57
Foam::coupledPolyPatch::parallel
virtual bool parallel() const
Are the cyclic planes parallel.
Definition: coupledPolyPatch.H:295
Foam::cyclicAMIFvPatch::TypeName
TypeName(cyclicAMIPolyPatch::typeName_())
Runtime type information.
Foam::Field< scalar >
Foam::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:68
Foam::cyclicAMIFvPatch::reverseT
virtual const tensorField & reverseT() const
Return neighbour-cell transformation tensor.
Definition: cyclicAMIFvPatch.H:192
fld
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< ' ';}gmvFile<< nl;for(const word &name :lagrangianScalarNames){ IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Definition: gmvOutputLagrangian.H:23
Foam::cyclicAMIPolyPatch::interpolate
tmp< Field< Type > > interpolate(const Field< Type > &fld, const UList< Type > &defaultValues=UList< Type >()) const
Interpolate field.
fvBoundaryMesh.H
Foam::cyclicAMIFvPatch::coupled
virtual bool coupled() const
Definition: cyclicAMIFvPatch.C:86
Foam::cyclicAMIPolyPatch::neighbPatchID
virtual label neighbPatchID() const
Neighbour patch ID.
Definition: cyclicAMIPolyPatch.C:902
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::refCast
To & refCast(From &r)
Reference type cast template function.
Definition: typeInfo.H:131
cyclicAMIPolyPatch.H
Foam::UPstream::commsTypes
commsTypes
Types of communications.
Definition: UPstream.H:69
Foam::cyclicAMIPolyPatch::AMI
const AMIPatchToPatchInterpolation & AMI() const
Return a reference to the AMI interpolator.
Definition: cyclicAMIPolyPatch.C:977
Foam::cyclicAMIFvPatch::interfaceInternalField
virtual tmp< labelField > interfaceInternalField(const labelUList &internalData) const
Definition: cyclicAMIFvPatch.C:217
Foam::cyclicAMIFvPatch::owner
virtual bool owner() const
Definition: cyclicAMIFvPatch.H:152
Foam::cyclicAMIPolyPatch::applyLowWeightCorrection
bool applyLowWeightCorrection() const
Return true if applying the low weight correction.
Definition: cyclicAMIPolyPatch.C:995
Foam::cyclicAMIFvPatch::makeWeights
void makeWeights(scalarField &) const
Make patch weighting factors.
Definition: cyclicAMIFvPatch.C:94
Foam::AMIInterpolation
Interpolation class dealing with transfer of data between two primitive patches with an arbitrary mes...
Definition: AMIInterpolation.H:79
Foam::UList< Type >
Foam::cyclicAMIFvPatch::cyclicAMIFvPatch
cyclicAMIFvPatch(const polyPatch &patch, const fvBoundaryMesh &bm)
Construct from polyPatch.
Definition: cyclicAMIFvPatch.H:93
Foam::fvPatch::patch
const polyPatch & patch() const
Return the polyPatch.
Definition: fvPatch.H:161
Foam::cyclicAMIFvPatch::cyclicAMIPatch
const cyclicAMIPolyPatch & cyclicAMIPatch() const
Return local reference cast into the cyclic patch.
Definition: cyclicAMIFvPatch.H:141
Foam::cyclicAMIFvPatch
Cyclic patch for Arbitrary Mesh Interface (AMI)
Definition: cyclicAMIFvPatch.H:54
Foam::coupledFvPatch
An abstract base class for patches that couple regions of the computational domain e....
Definition: coupledFvPatch.H:53
cyclicAMILduInterface.H
Foam::cyclicAMIFvPatch::forwardT
virtual const tensorField & forwardT() const
Return face transformation tensor.
Definition: cyclicAMIFvPatch.H:186
Foam::cyclicAMIFvPatch::parallel
virtual bool parallel() const
Are the cyclic planes parallel.
Definition: cyclicAMIFvPatch.H:180
Foam::faceCells
Smooth ATC in cells next to a set of patches supplied by type.
Definition: faceCells.H:56
Foam::cyclicAMIFvPatch::movePoints
virtual void movePoints()
Correct patches after moving points.
Definition: cyclicAMIFvPatch.C:245
Foam::cyclicAMIPolyPatch
Cyclic patch for Arbitrary Mesh Interface (AMI)
Definition: cyclicAMIPolyPatch.H:68