fvPatch.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-2018 OpenFOAM Foundation
9  Copyright (C) 2020-2023 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::fvPatch
29 
30 Description
31  A finiteVolume patch using a polyPatch and a fvBoundaryMesh
32 
33 SourceFiles
34  fvPatch.C
35  fvPatchNew.C
36  fvPatchTemplates.C
37  fvPatchFvMeshTemplates.C
38 
39 \*---------------------------------------------------------------------------*/
40 
41 #ifndef Foam_fvPatch_H
42 #define Foam_fvPatch_H
43 
44 #include "polyPatch.H"
45 #include "labelList.H"
46 #include "SubList.H"
47 #include "SubField.H"
48 #include "PtrList.H"
49 #include "typeInfo.H"
50 #include "autoPtr.H"
51 #include "tmp.H"
52 #include "primitiveFields.H"
53 #include "fvPatchFieldsFwd.H"
54 #include "runTimeSelectionTables.H"
55 
56 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
57 
58 namespace Foam
59 {
60 
61 // Forward Declarations
62 class fvBoundaryMesh;
63 class fvPatch;
65 
66 //- Store lists of fvPatch as a PtrList
68 
69 /*---------------------------------------------------------------------------*\
70  Class fvPatch Declaration
71 \*---------------------------------------------------------------------------*/
72 
73 class fvPatch
74 {
75  // Private Data
76 
77  //- Reference to the underlying polyPatch
78  const polyPatch& polyPatch_;
79 
80  //- Reference to boundary mesh
81  const fvBoundaryMesh& boundaryMesh_;
82 
83 
84  // Private Member Functions
85 
86  //- No copy construct
87  fvPatch(const fvPatch&) = delete;
88 
89  //- No copy assignment
90  void operator=(const fvPatch&) = delete;
91 
92 
93 public:
94 
95  // Protected Member Functions
96 
97  //- Make patch weighting factors
98  virtual void makeWeights(scalarField&) const;
99 
100  //- Correct patch deltaCoeffs
101  virtual void makeDeltaCoeffs(scalarField&) const;
102 
103  //- Correct patch non-ortho deltaCoeffs
104  virtual void makeNonOrthoDeltaCoeffs(scalarField&) const;
105 
106  //- Correct patch non-ortho correction vectors
107  virtual void makeNonOrthoCorrVectors(vectorField&) const;
108 
109  //- Initialise the patches for moving points
110  virtual void initMovePoints();
111 
112  //- Correct patches after moving points
113  virtual void movePoints();
114 
115 
116 public:
117 
118  //- The boundary type associated with the patch
120 
121  friend class fvBoundaryMesh;
122  friend class surfaceInterpolation;
123 
124  //- Runtime type information
125  TypeName(polyPatch::typeName_());
126 
127 
128  // Declare run-time constructor selection tables
129 
131  (
132  autoPtr,
133  fvPatch,
134  polyPatch,
135  (const polyPatch& patch, const fvBoundaryMesh& bm),
136  (patch, bm)
137  );
139 
140  // Constructors
142  //- Construct from polyPatch and fvBoundaryMesh
143  fvPatch(const polyPatch&, const fvBoundaryMesh&);
144 
145 
146  // Selectors
147 
148  //- Return a pointer to a new patch created on freestore from polyPatch
149  static autoPtr<fvPatch> New
150  (
151  const polyPatch&,
152  const fvBoundaryMesh&
153  );
154 
155 
156  //- Destructor
157  virtual ~fvPatch();
158 
159 
160  // Member Functions
161 
162  //- Lookup the polyPatch index on corresponding fvMesh
163  // \note Fatal if the polyPatch is not associated with a fvMesh
164  static const fvPatch& lookupPatch(const polyPatch& p);
165 
166 
167  // Access
168 
169  //- Return the polyPatch
170  const polyPatch& patch() const noexcept
171  {
172  return polyPatch_;
173  }
174 
175  //- Return name
176  virtual const word& name() const
177  {
178  return polyPatch_.name();
179  }
180 
181  //- The index of this patch in the boundary mesh
182  label index() const noexcept
183  {
184  return polyPatch_.index();
185  }
186 
187  //- The patch start within the polyMesh face list
188  label start() const noexcept
189  {
190  return polyPatch_.start();
191  }
192 
193  //- Patch size is the number of faces, but can be overloaded
194  virtual label size() const
195  {
196  return polyPatch_.size();
197  }
198 
199  //- Return true if this patch is coupled
200  virtual bool coupled() const
201  {
202  return polyPatch_.coupled();
203  }
204 
205  //- Return true if the given type is a constraint type
206  static bool constraintType(const word& patchType);
207 
208  //- Return a list of all the constraint patch types
209  static wordList constraintTypes();
211  //- Return boundaryMesh reference
212  const fvBoundaryMesh& boundaryMesh() const noexcept
213  {
214  return boundaryMesh_;
215  }
216 
217  //- This patch slice from the complete list, which has size
218  //- mesh::nFaces(), using the virtual patch size.
219  template<class T>
220  const typename List<T>::subList
221  patchSlice(const List<T>& values) const
222  {
223  return typename List<T>::subList(values, size(), start());
224  }
225 
226  //- This patch slice from the list of boundary values, which has size
227  //- mesh::nBoundaryFaces(), using the virtual patch size.
228  template<class T>
229  const typename List<T>::subList
230  boundarySlice(const List<T>& values) const
231  {
232  return typename List<T>::subList
233  (
235  size(),
236  polyPatch_.offset()
237  );
238  }
239 
240  //- Return faceCells
241  virtual const labelUList& faceCells() const;
243 
244  // Access functions for geometrical data
245 
246  //- Return face centres
247  const vectorField& Cf() const;
248 
249  //- Return neighbour cell centres
250  tmp<vectorField> Cn() const;
251 
252  //- Return face area vectors, like the fvMesh::Sf() method
253  const vectorField& Sf() const;
254 
255  //- Return face area magnitudes, like the fvMesh::magSf() method
256  const scalarField& magSf() const;
257 
258  //- Return face unit normals, like the fvMesh::unitSf() method.
259  //- Same as nf().
261 
262  //- Return face unit normals, like the fvMesh::unitSf() method
263  //- Same as unitSf().
264  tmp<vectorField> nf() const;
265 
266  //- Return cell-centre to face-centre vector
267  //- except for coupled patches for which the cell-centre
268  //- to coupled-cell-centre vector is returned
269  virtual tmp<vectorField> delta() const;
270 
272  // Access functions for demand driven data
273 
274  //- Return patch weighting factors
275  const scalarField& weights() const;
276 
277  //- Return the face - cell distance coefficient
278  //- except for coupled patches for which the cell-centre
279  //- to coupled-cell-centre distance coefficient is returned
280  const scalarField& deltaCoeffs() const;
281 
283  // Evaluation Functions
284 
285  //- Extract internal field next to patch using specified addressing
286  // \param internalData The internal field to extract from
287  // \param addressing Addressing from patch into internal field
288  // \param [out] pfld The extracted patch field.
289  // It is always resized according to the patch size(),
290  // which can be smaller than the addressing size
291  template<class Type>
292  inline void patchInternalField
293  (
294  const UList<Type>& internalData,
295  const labelUList& addressing,
296  Field<Type>& pfld
297  ) const;
298 
299  //- Extract internal field next to patch as patch field
300  //- using faceCells() mapping.
301  // \param internalData The internal field to extract from
302  // \param [out] pfld The extracted patch field.
303  // It is always resized according to the patch size(),
304  // which can be smaller than the faceCells() size
305  template<class Type>
306  void patchInternalField
307  (
308  const UList<Type>& internalData,
309  Field<Type>& pfld
310  ) const;
311 
312  //- Return given internal field next to patch as patch field
313  //- using faceCells() mapping.
314  // \param internalData The internal field to extract from
315  template<class Type>
317  (
318  const UList<Type>& internalData
319  ) const;
320 
321  //- Return the patch field of the GeometricField
322  //- corresponding to this patch.
323  template<class GeometricField, class AnyType = bool>
324  const typename GeometricField::Patch& patchField
325  (
326  const GeometricField& gf
327  ) const;
328 
329  //- Lookup the named field from the local registry and
330  //- return the patch field corresponding to this patch.
331  // N.B. The dummy pointer arguments are used if this function is
332  // instantiated within a templated function to avoid a bug in gcc.
333  template<class GeometricField, class AnyType = bool>
335  (
336  const word& name,
337  const GeometricField* = nullptr,
338  const AnyType* = nullptr
339  ) const;
340 };
341 
342 
343 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
344 
345 } // End namespace Foam
346 
347 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
348 
349 #ifdef NoRepository
350  #include "fvPatchTemplates.C"
351 #endif
352 
353 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
354 
355 #endif
356 
357 // ************************************************************************* //
label start() const noexcept
The patch start within the polyMesh face list.
Definition: fvPatch.H:226
const GeometricField::Patch & patchField(const GeometricField &gf) const
Return the patch field of the GeometricField corresponding to this patch.
const List< T >::subList patchSlice(const List< T > &values) const
This patch slice from the complete list, which has size mesh::nFaces(), using the virtual patch size...
Definition: fvPatch.H:271
virtual void initMovePoints()
Initialise the patches for moving points.
Definition: fvPatch.C:175
const scalarField & magSf() const
Return face area magnitudes, like the fvMesh::magSf() method.
Definition: fvPatch.C:131
label start() const noexcept
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:441
static const fvPatch & lookupPatch(const polyPatch &p)
Lookup the polyPatch index on corresponding fvMesh.
Definition: fvPatch.C:42
tmp< vectorField > nf() const
Return face unit normals, like the fvMesh::unitSf() method Same as unitSf().
Definition: fvPatch.C:143
virtual ~fvPatch()
Destructor.
Definition: fvPatch.C:68
virtual void makeNonOrthoCorrVectors(vectorField &) const
Correct patch non-ortho correction vectors.
Definition: fvPatch.C:171
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:70
Generic GeometricField class.
Definition: areaFieldsFwd.H:50
SubList< T > subList
Declare type of subList.
Definition: List.H:144
Cell to surface interpolation scheme. Included in fvMesh.
const GeometricField::Patch & lookupPatchField(const word &name, const GeometricField *=nullptr, const AnyType *=nullptr) const
Lookup the named field from the local registry and return the patch field corresponding to this patch...
faPatchField<Type> abstract base class. This class gives a fat-interface to all derived classes cover...
Definition: areaFieldsFwd.H:56
declareRunTimeSelectionTable(autoPtr, fvPatch, polyPatch,(const polyPatch &patch, const fvBoundaryMesh &bm),(patch, bm))
virtual void makeNonOrthoDeltaCoeffs(scalarField &) const
Correct patch non-ortho deltaCoeffs.
Definition: fvPatch.C:167
tmp< vectorField > unitSf() const
Return face unit normals, like the fvMesh::unitSf() method. Same as nf().
Definition: fvPatch.C:137
List< T > values(const HashTable< T, Key, Hash > &tbl, const bool doSort=false)
List of values from HashTable, optionally sorted.
Definition: HashOps.H:164
void patchInternalField(const UList< Type > &internalData, const labelUList &addressing, Field< Type > &pfld) const
Extract internal field next to patch using specified addressing.
A List obtained as a section of another List.
Definition: SubList.H:50
virtual bool coupled() const
Return true if this patch is geometrically coupled (i.e. faces and.
Definition: polyPatch.H:464
fvBoundaryMesh BoundaryMesh
The boundary type associated with the patch.
Definition: fvPatch.H:138
A class for handling words, derived from Foam::string.
Definition: word.H:63
virtual const labelUList & faceCells() const
Return faceCells.
Definition: fvPatch.C:107
label offset() const noexcept
The offset where this patch starts in the boundary face list.
Definition: polyPatch.C:302
const scalarField & weights() const
Return patch weighting factors.
Definition: fvPatch.C:189
virtual void makeWeights(scalarField &) const
Make patch weighting factors.
Definition: fvPatch.C:157
const vectorField & Cf() const
Return face centres.
Definition: fvPatch.C:113
const List< T >::subList boundarySlice(const List< T > &values) const
This patch slice from the list of boundary values, which has size mesh::nBoundaryFaces(), using the virtual patch size.
Definition: fvPatch.H:282
const fvBoundaryMesh & boundaryMesh() const noexcept
Return boundaryMesh reference.
Definition: fvPatch.H:260
virtual label size() const
Patch size is the number of faces, but can be overloaded.
Definition: fvPatch.H:234
const word & name() const noexcept
The patch name.
const direction noexcept
Definition: Scalar.H:258
TypeName(polyPatch::typeName_())
Runtime type information.
static bool constraintType(const word &patchType)
Return true if the given type is a constraint type.
Definition: fvPatch.C:74
const vectorField & Sf() const
Return face area vectors, like the fvMesh::Sf() method.
Definition: fvPatch.C:125
Specialisations of Field<T> for scalar, vector and tensor.
const scalarField & deltaCoeffs() const
Return the face - cell distance coefficient except for coupled patches for which the cell-centre to c...
Definition: fvPatch.C:183
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers...
Definition: List.H:55
Foam::fvBoundaryMesh.
tmp< vectorField > Cn() const
Return neighbour cell centres.
Definition: fvPatch.C:119
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
virtual const word & name() const
Return name.
Definition: fvPatch.H:210
static autoPtr< fvPatch > New(const polyPatch &, const fvBoundaryMesh &)
Return a pointer to a new patch created on freestore from polyPatch.
Definition: fvPatchNew.C:28
Macros to ease declaration of run-time selection tables.
label index() const noexcept
The index of this patch in the boundaryMesh.
volScalarField & p
A class for managing temporary objects.
Definition: HashPtrTable.H:50
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:69
const polyPatch & patch() const noexcept
Return the polyPatch.
Definition: fvPatch.H:202
label index() const noexcept
The index of this patch in the boundary mesh.
Definition: fvPatch.H:218
virtual bool coupled() const
Return true if this patch is coupled.
Definition: fvPatch.H:242
virtual tmp< vectorField > delta() const
Return cell-centre to face-centre vector except for coupled patches for which the cell-centre to coup...
Definition: fvPatch.C:149
virtual void movePoints()
Correct patches after moving points.
Definition: fvPatch.C:179
Namespace for OpenFOAM.
static wordList constraintTypes()
Return a list of all the constraint patch types.
Definition: fvPatch.C:85
PtrList< fvPatch > fvPatchList
Store lists of fvPatch as a PtrList.
Definition: fvPatch.H:59
virtual void makeDeltaCoeffs(scalarField &) const
Correct patch deltaCoeffs.
Definition: fvPatch.C:163