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 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
37\*---------------------------------------------------------------------------*/
38
39#ifndef Foam_fvPatch_H
40#define Foam_fvPatch_H
41
42#include "polyPatch.H"
43#include "labelList.H"
44#include "SubList.H"
45#include "SubField.H"
46#include "PtrList.H"
47#include "typeInfo.H"
48#include "autoPtr.H"
49#include "tmp.H"
50#include "primitiveFields.H"
51#include "fvPatchFieldsFwd.H"
53
54// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
55
56namespace Foam
57{
58
59// Forward Declarations
60class fvBoundaryMesh;
61class fvPatch;
62class surfaceInterpolation;
63
64//- Store lists of fvPatch as a PtrList
66
67/*---------------------------------------------------------------------------*\
68 Class fvPatch Declaration
69\*---------------------------------------------------------------------------*/
71class fvPatch
72{
73 // Private Data
74
75 //- Reference to the underlying polyPatch
76 const polyPatch& polyPatch_;
77
78 //- Reference to boundary mesh
79 const fvBoundaryMesh& boundaryMesh_;
80
81
82 // Private Member Functions
83
84 //- No copy construct
85 fvPatch(const fvPatch&) = delete;
86
87 //- No copy assignment
88 void operator=(const fvPatch&) = delete;
89
90
91public:
92
93 // Protected Member Functions
94
95 //- Make patch weighting factors
96 virtual void makeWeights(scalarField&) const;
97
98 //- Correct patch deltaCoeffs
99 virtual void makeDeltaCoeffs(scalarField&) const;
100
101 //- Correct patch non-ortho deltaCoeffs
102 virtual void makeNonOrthoDeltaCoeffs(scalarField&) const;
103
104 //- Correct patch non-ortho correction vectors
105 virtual void makeNonOrthoCorrVectors(vectorField&) const;
106
107 //- Initialise the patches for moving points
108 virtual void initMovePoints();
109
110 //- Correct patches after moving points
111 virtual void movePoints();
112
113
114public:
115
116 //- The boundary type associated with the patch
119 friend class fvBoundaryMesh;
120 friend class surfaceInterpolation;
121
122 //- Runtime type information
123 TypeName(polyPatch::typeName_());
124
125
126 // Declare run-time constructor selection tables
129 (
130 autoPtr,
131 fvPatch,
132 polyPatch,
133 (const polyPatch& patch, const fvBoundaryMesh& bm),
134 (patch, bm)
135 );
136
137
138 // Constructors
139
140 //- Construct from polyPatch and fvBoundaryMesh
141 fvPatch(const polyPatch&, const fvBoundaryMesh&);
142
143
144 // Selectors
145
146 //- Return a pointer to a new patch created on freestore from polyPatch
147 static autoPtr<fvPatch> New
148 (
149 const polyPatch&,
150 const fvBoundaryMesh&
151 );
152
153
154 //- Destructor
155 virtual ~fvPatch();
156
157
158 // Member Functions
159
160 //- Lookup the polyPatch index on corresponding fvMesh
161 // \note Fatal if the polyPatch is not associated with a fvMesh
162 static const fvPatch& lookupPatch(const polyPatch& p);
163
164
165 // Access
166
167 //- Return the polyPatch
168 const polyPatch& patch() const
169 {
170 return polyPatch_;
171 }
172
173 //- Return name
174 virtual const word& name() const
175 {
176 return polyPatch_.name();
177 }
178
179 //- Return start label of this patch in the polyMesh face list
180 virtual label start() const
181 {
182 return polyPatch_.start();
183 }
184
185 //- Return size
186 virtual label size() const
187 {
188 return polyPatch_.size();
189 }
190
191 //- Return true if this patch is coupled
192 virtual bool coupled() const
193 {
194 return polyPatch_.coupled();
195 }
196
197 //- Return true if the given type is a constraint type
198 static bool constraintType(const word& pt);
199
200 //- Return a list of all the constraint patch types
201 static wordList constraintTypes();
202
203 //- Return the index of this patch in the fvBoundaryMesh
204 label index() const
205 {
206 return polyPatch_.index();
207 }
208
209 //- Return boundaryMesh reference
210 const fvBoundaryMesh& boundaryMesh() const
211 {
212 return boundaryMesh_;
213 }
214
215 //- Slice list to patch
216 template<class T>
217 const typename List<T>::subList patchSlice(const List<T>& l) const
218 {
219 return typename List<T>::subList(l, size(), start());
220 }
221
222 //- Return faceCells
223 virtual const labelUList& faceCells() const;
224
225
226 // Access functions for geometrical data
227
228 //- Return face centres
229 const vectorField& Cf() const;
230
231 //- Return neighbour cell centres
232 tmp<vectorField> Cn() const;
233
234 //- Return face area vectors
235 const vectorField& Sf() const;
236
237 //- Return face area magnitudes
238 const scalarField& magSf() const;
239
240 //- Return face normals
241 tmp<vectorField> nf() const;
242
243 //- Return cell-centre to face-centre vector
244 //- except for coupled patches for which the cell-centre
245 //- to coupled-cell-centre vector is returned
246 virtual tmp<vectorField> delta() const;
247
248
249 // Access functions for demand driven data
250
251 //- Return patch weighting factors
252 const scalarField& weights() const;
253
254 //- Return the face - cell distance coefficient
255 //- except for coupled patches for which the cell-centre
256 //- to coupled-cell-centre distance coefficient is returned
257 const scalarField& deltaCoeffs() const;
258
259
260 // Evaluation functions
261
262 //- Return given internal field next to patch as patch field
263 template<class Type>
265
266 //- Return given internal field next to patch as patch field
267 //- using faceCells mapping
268 template<class Type>
270 (
271 const UList<Type>&,
272 const labelUList& faceCells
273 ) const;
274
275 //- Return given internal field next to patch as patch field
276 template<class Type>
277 void patchInternalField(const UList<Type>&, Field<Type>&) const;
278
279 //- Return the corresponding patchField of the named field
280 template<class GeometricField, class Type>
281 const typename GeometricField::Patch& patchField
282 (
283 const GeometricField&
284 ) const;
285
286 //- Lookup and return the patchField of the named field from the
287 //- local objectRegistry.
288 // N.B. The dummy pointer arguments are used if this function is
289 // instantiated within a templated function to avoid a bug in gcc.
290 // See inletOutletFvPatchField.C and outletInletFvPatchField.C
291 template<class GeometricField, class Type>
293 (
294 const word& name,
295 const GeometricField* = nullptr,
296 const Type* = nullptr
297 ) const;
298};
299
300
301// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
302
303} // End namespace Foam
304
305// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
306
307#ifdef NoRepository
308 #include "fvPatchTemplates.C"
309#endif
310
311// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
312
313#endif
314
315// ************************************************************************* //
Generic GeometricField class.
PatchField< Type > Patch
The patch field type for the GeometricBoundaryField.
SubList< T > subList
Declare type of subList.
Definition: List.H:111
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: PtrList.H:73
A List obtained as a section of another List.
Definition: SubList.H:70
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
Smooth ATC in cells next to a set of patches supplied by type.
Definition: faceCells.H:59
Foam::fvBoundaryMesh.
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:91
virtual label size() const
Return size.
Definition: fvPatch.H:185
virtual const word & name() const
Return name.
Definition: fvPatch.H:173
virtual ~fvPatch()
Destructor.
Definition: fvPatch.C:75
static const fvPatch & lookupPatch(const polyPatch &p)
Lookup the polyPatch index on corresponding fvMesh.
Definition: fvPatch.C:49
virtual bool coupled() const
Return true if this patch is coupled.
Definition: fvPatch.H:191
const GeometricField::Patch & patchField(const GeometricField &) const
Return the corresponding patchField of the named field.
virtual void makeWeights(scalarField &) const
Make patch weighting factors.
Definition: fvPatch.C:170
label index() const
Return the index of this patch in the fvBoundaryMesh.
Definition: fvPatch.H:203
virtual label start() const
Return start label of this patch in the polyMesh face list.
Definition: fvPatch.H:179
static autoPtr< fvPatch > New(const polyPatch &, const fvBoundaryMesh &)
Return a pointer to a new patch created on freestore from polyPatch.
Definition: fvPatchNew.C:35
virtual void makeNonOrthoDeltaCoeffs(scalarField &) const
Correct patch non-ortho deltaCoeffs.
Definition: fvPatch.C:180
virtual void makeDeltaCoeffs(scalarField &) const
Correct patch deltaCoeffs.
Definition: fvPatch.C:176
tmp< Field< Type > > patchInternalField(const UList< Type > &, const labelUList &faceCells) const
virtual void movePoints()
Correct patches after moving points.
Definition: fvPatch.C:192
tmp< vectorField > Cn() const
Return neighbour cell centres.
Definition: fvPatch.C:125
const GeometricField::Patch & lookupPatchField(const word &name, const GeometricField *=nullptr, const Type *=nullptr) const
declareRunTimeSelectionTable(autoPtr, fvPatch, polyPatch,(const polyPatch &patch, const fvBoundaryMesh &bm),(patch, bm))
tmp< Field< Type > > patchInternalField(const UList< Type > &) const
Return given internal field next to patch as patch field.
virtual tmp< vectorField > delta() const
Definition: fvPatch.C:162
virtual void initMovePoints()
Initialise the patches for moving points.
Definition: fvPatch.C:188
static bool constraintType(const word &pt)
Return true if the given type is a constraint type.
Definition: fvPatch.C:81
virtual void makeNonOrthoCorrVectors(vectorField &) const
Correct patch non-ortho correction vectors.
Definition: fvPatch.C:184
const scalarField & magSf() const
Return face area magnitudes.
Definition: fvPatch.C:156
const scalarField & weights() const
Return patch weighting factors.
Definition: fvPatch.C:202
const fvBoundaryMesh & boundaryMesh() const
Return boundaryMesh reference.
Definition: fvPatch.H:209
tmp< vectorField > nf() const
Return face normals.
Definition: fvPatch.C:144
fvBoundaryMesh BoundaryMesh
The boundary type associated with the patch.
Definition: fvPatch.H:116
const polyPatch & patch() const
Return the polyPatch.
Definition: fvPatch.H:167
const vectorField & Cf() const
Return face centres.
Definition: fvPatch.C:119
const scalarField & deltaCoeffs() const
Definition: fvPatch.C:196
const List< T >::subList patchSlice(const List< T > &l) const
Slice list to patch.
Definition: fvPatch.H:216
TypeName(polyPatch::typeName_())
Runtime type information.
const vectorField & Sf() const
Return face area vectors.
Definition: fvPatch.C:150
virtual const labelUList & faceCells() const
Return faceCells.
Definition: fvPatch.C:113
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:75
virtual bool coupled() const
Return true if this patch is geometrically coupled (i.e. faces and.
Definition: polyPatch.H:380
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:364
Cell to surface interpolation scheme. Included in fvMesh.
A class for managing temporary objects.
Definition: tmp.H:65
A class for handling words, derived from Foam::string.
Definition: word.H:68
volScalarField & p
Namespace for OpenFOAM.
PtrList< fvPatch > fvPatchList
Store lists of fvPatch as a PtrList.
Definition: fvPatch.H:64
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)
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition: typeInfo.H:73