Loading...
Searching...
No Matches
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-2025 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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
27Class
28 Foam::fvPatch
29
30Description
31 A finiteVolume patch using a polyPatch and a fvBoundaryMesh
32
33SourceFiles
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"
55
56// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
57
58namespace Foam
59{
60
61// Forward Declarations
62class fvBoundaryMesh;
63class fvPatch;
65
66//- Store lists of fvPatch as a PtrList
68
69/*---------------------------------------------------------------------------*\
70 Class fvPatch Declaration
71\*---------------------------------------------------------------------------*/
72
73class 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
93public:
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
116public:
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 //- The offset of this patch within the boundary face list
194 label offset() const noexcept
195 {
196 return polyPatch_.offset();
197 }
198
199 //- Patch size is the number of faces, but can be overloaded
200 virtual label size() const
201 {
202 return polyPatch_.size();
203 }
204
205 //- Return true if this patch is coupled
206 virtual bool coupled() const
207 {
208 return polyPatch_.coupled();
209 }
211 //- Return true if the given type is a constraint type
212 static bool constraintType(const word& patchType);
213
214 //- Return a list of all the constraint patch types
215 static wordList constraintTypes();
216
217 //- Return boundaryMesh reference
219 {
220 return boundaryMesh_;
221 }
222
223 //- This patch slice from the complete list, which has size
224 //- mesh::nFaces(), using the virtual patch size.
225 template<class T>
226 const typename List<T>::subList
227 patchSlice(const UList<T>& values) const
228 {
229 return typename List<T>::subList(values, size(), start());
230 }
231
232 //- This patch slice from the list of boundary values, which has size
233 //- mesh::nBoundaryFaces(), using the virtual patch size.
234 template<class T>
235 const typename List<T>::subList
236 boundarySlice(const UList<T>& values) const
237 {
238 return typename List<T>::subList(values, size(), offset());
239 }
240
241 //- Return faceCells
242 virtual const labelUList& faceCells() const;
243
244
245 // Access functions for geometrical data
246
247 //- Return face centres
248 const vectorField& Cf() const;
249
250 //- Return neighbour cell centres
251 tmp<vectorField> Cn() const;
252
253 //- Return face area vectors, like the fvMesh::Sf() method
254 const vectorField& Sf() const;
255
256 //- Return face area magnitudes, like the fvMesh::magSf() method
257 const scalarField& magSf() const;
258
259 //- Return face unit normals, like the fvMesh::unitSf() method.
260 //- Same as nf().
261 tmp<vectorField> unitSf() const;
262
263 //- Return face unit normals, like the fvMesh::unitSf() method
264 //- Same as unitSf().
265 tmp<vectorField> nf() const;
266
267 //- Return cell-centre to face-centre vector
268 //- except for coupled patches for which the cell-centre
269 //- to coupled-cell-centre vector is returned
270 virtual tmp<vectorField> delta() const;
271
272
273 // Access functions for demand driven data
274
275 //- Return patch weighting factors
276 const scalarField& weights() const;
277
278 //- Return the face - cell distance coefficient
279 //- except for coupled patches for which the cell-centre
280 //- to coupled-cell-centre distance coefficient is returned
281 const scalarField& deltaCoeffs() const;
282
283
284 // Evaluation Functions
285
286 //- Extract internal field next to patch using specified addressing
287 // \param internalData The internal field to extract from
288 // \param addressing Addressing from patch into internal field
289 // \param [out] pfld The extracted patch field.
290 // Should normally be sized according to the patch size(),
291 // which can be smaller than the addressing size
292 template<class Type>
293 inline void patchInternalField
294 (
295 const UList<Type>& internalData,
296 const labelUList& addressing,
297 UList<Type>& pfld
298 ) const;
299
300 //- Extract internal field next to patch as patch field
301 //- using faceCells() mapping.
302 // \param internalData The internal field to extract from
303 // \param [out] pfld The extracted patch field.
304 // Should normally be sized according to the patch size(),
305 // which can be smaller than the faceCells() size
306 template<class Type>
308 (
309 const UList<Type>& internalData,
310 UList<Type>& pfld
311 ) const;
312
313 //- Return given internal field next to patch as patch field
314 //- using faceCells() mapping.
315 // \param internalData The internal field to extract from
316 template<class Type>
317 [[nodiscard]] tmp<Field<Type>>
319 (
320 const UList<Type>& internalData
321 ) const;
322
323
324 // Lookup
325
326 //- Return the patch field of the GeometricField
327 //- corresponding to this patch.
328 template<class GeometricField, class AnyType = bool>
329 const typename GeometricField::Patch&
330 patchField(const GeometricField& gf) const;
331
332 //- Lookup the named field from the local registry and
333 //- return the patch field corresponding to this patch.
334 template<class GeometricField, class AnyType = bool>
335 const typename GeometricField::Patch&
336 lookupPatchField(const word& name) const;
337
338 //- Find the named field (if any) from the local registry and
339 //- return the patch field corresponding to this patch.
340 template<class GeometricField>
341 const typename GeometricField::Patch*
342 cfindPatchField(const word& name) const;
343};
344
345
346// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
347
348} // End namespace Foam
349
350// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
351
352#ifdef NoRepository
353 #include "fvPatchTemplates.C"
354#endif
355
356// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
357
358#endif
359
360// ************************************************************************* //
Generic GeometricField class.
PatchField< Type > Patch
The patch field type for the GeometricBoundaryField.
SubList< T > subList
Declare type of subList.
Definition List.H:129
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition PtrList.H:67
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition UList.H:89
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
A fvBoundaryMesh is a fvPatch list with a reference to the associated fvMesh, with additional search ...
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition fvPatch.H:71
static wordList constraintTypes()
Return a list of all the constraint patch types.
Definition fvPatch.C:85
virtual label size() const
Patch size is the number of faces, but can be overloaded.
Definition fvPatch.H:242
virtual const word & name() const
Return name.
Definition fvPatch.H:210
virtual ~fvPatch()
Destructor.
Definition fvPatch.C:68
static const fvPatch & lookupPatch(const polyPatch &p)
Lookup the polyPatch index on corresponding fvMesh.
Definition fvPatch.C:42
label offset() const noexcept
The offset of this patch within the boundary face list.
Definition fvPatch.H:234
const polyPatch & patch() const noexcept
Return the polyPatch.
Definition fvPatch.H:202
virtual bool coupled() const
Return true if this patch is coupled.
Definition fvPatch.H:250
virtual void makeWeights(scalarField &) const
Make patch weighting factors.
Definition fvPatch.C:157
static autoPtr< fvPatch > New(const polyPatch &, const fvBoundaryMesh &)
Return a pointer to a new patch created on freestore from polyPatch.
Definition fvPatchNew.C:28
virtual void makeNonOrthoDeltaCoeffs(scalarField &) const
Correct patch non-ortho deltaCoeffs.
Definition fvPatch.C:167
virtual void makeDeltaCoeffs(scalarField &) const
Correct patch deltaCoeffs.
Definition fvPatch.C:163
virtual void movePoints()
Correct patches after moving points.
Definition fvPatch.C:179
tmp< vectorField > Cn() const
Return neighbour cell centres.
Definition fvPatch.C:119
const GeometricField::Patch * cfindPatchField(const word &name) const
Find the named field (if any) from the local registry and return the patch field corresponding to thi...
const GeometricField::Patch & patchField(const GeometricField &gf) const
Return the patch field of the GeometricField corresponding to this patch.
label index() const noexcept
The index of this patch in the boundary mesh.
Definition fvPatch.H:218
declareRunTimeSelectionTable(autoPtr, fvPatch, polyPatch,(const polyPatch &patch, const fvBoundaryMesh &bm),(patch, bm))
const fvBoundaryMesh & boundaryMesh() const noexcept
Return boundaryMesh reference.
Definition fvPatch.H:268
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 initMovePoints()
Initialise the patches for moving points.
Definition fvPatch.C:175
virtual void makeNonOrthoCorrVectors(vectorField &) const
Correct patch non-ortho correction vectors.
Definition fvPatch.C:171
friend class surfaceInterpolation
Definition fvPatch.H:141
const scalarField & magSf() const
Return face area magnitudes, like the fvMesh::magSf() method.
Definition fvPatch.C:131
const scalarField & weights() const
Return patch weighting factors.
Definition fvPatch.C:189
friend class fvBoundaryMesh
Definition fvPatch.H:140
tmp< vectorField > unitSf() const
Return face unit normals, like the fvMesh::unitSf() method. Same as nf().
Definition fvPatch.C:137
label start() const noexcept
The patch start within the polyMesh face list.
Definition fvPatch.H:226
const List< T >::subList patchSlice(const UList< T > &values) const
This patch slice from the complete list, which has size mesh::nFaces(), using the virtual patch size.
Definition fvPatch.H:279
tmp< vectorField > nf() const
Return face unit normals, like the fvMesh::unitSf() method Same as unitSf().
Definition fvPatch.C:143
fvBoundaryMesh BoundaryMesh
The boundary type associated with the patch.
Definition fvPatch.H:138
const vectorField & Cf() const
Return face centres.
Definition fvPatch.C:113
void patchInternalField(const UList< Type > &internalData, const labelUList &addressing, UList< Type > &pfld) const
Extract internal field next to patch using specified addressing.
const List< T >::subList boundarySlice(const UList< T > &values) const
This patch slice from the list of boundary values, which has size mesh::nBoundaryFaces(),...
Definition fvPatch.H:290
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
TypeName(polyPatch::typeName_())
Runtime type information.
const vectorField & Sf() const
Return face area vectors, like the fvMesh::Sf() method.
Definition fvPatch.C:125
static bool constraintType(const word &patchType)
Return true if the given type is a constraint type.
Definition fvPatch.C:74
virtual const labelUList & faceCells() const
Return faceCells.
Definition fvPatch.C:107
tmp< Field< Type > > patchInternalField(const UList< Type > &internalData) const
Return given internal field next to patch as patch field using faceCells() mapping.
const GeometricField::Patch & lookupPatchField(const word &name) const
Lookup the named field from the local registry and return the patch field corresponding to this patch...
label index() const noexcept
The index of this patch in the boundaryMesh.
const word & name() const noexcept
The patch name.
A patch is a list of labels that address the faces in the global face list.
Definition polyPatch.H:73
virtual bool coupled() const
Return true if this patch is geometrically coupled (i.e. faces and.
Definition polyPatch.H:469
Cell to surface interpolation scheme. Included in fvMesh.
A class for managing temporary objects.
Definition tmp.H:75
A class for handling words, derived from Foam::string.
Definition word.H:66
volScalarField & p
Namespace for OpenFOAM.
List< word > wordList
List of word.
Definition fileName.H:60
PtrList< fvPatch > fvPatchList
Store lists of fvPatch as a PtrList.
Definition fvPatch.H:64
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Field< vector > vectorField
Specialisation of Field<T> for vector.
const direction noexcept
Definition scalarImpl.H:265
UList< label > labelUList
A UList of labels.
Definition UList.H:75
Specialisations of Field<T> for scalar, vector and tensor.
Macros to ease declaration of run-time selection tables.
#define declareRunTimeSelectionTable(ptrWrapper, baseType, argNames, argList, parList)
Declare a run-time selection (variables and adder classes).
Basic run-time type information using word as the type's name. Used to enhance the standard RTTI to c...
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition typeInfo.H:68