fvPatch.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-2024 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::fvPatch
26 
27 Description
28  A finiteVolume patch using a polyPatch and a fvBoundaryMesh
29 
30 SourceFiles
31  fvPatch.C
32  fvPatchNew.C
33 
34 \*---------------------------------------------------------------------------*/
35 
36 #ifndef fvPatch_H
37 #define fvPatch_H
38 
39 #include "polyPatch.H"
40 #include "labelList.H"
41 #include "SubList.H"
42 #include "typeInfo.H"
43 #include "tmp.H"
44 #include "tmpNrc.H"
45 #include "primitiveFields.H"
46 #include "SubField.H"
47 #include "fvPatchFieldsFwd.H"
48 #include "autoPtr.H"
49 #include "runTimeSelectionTables.H"
50 
51 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
52 
53 namespace Foam
54 {
55 
56 class fvBoundaryMesh;
57 class objectRegistry;
58 class surfaceInterpolation;
59 
60 /*---------------------------------------------------------------------------*\
61  Class fvPatch Declaration
62 \*---------------------------------------------------------------------------*/
63 
64 class fvPatch
65 {
66  // Private Data
67 
68  //- Reference to the underlying polyPatch
69  const polyPatch& polyPatch_;
70 
71  //- Reference to boundary mesh
72  const fvBoundaryMesh& boundaryMesh_;
73 
74 
75 public:
76 
77  //- Runtime type information
78  TypeName(polyPatch::typeName_());
79 
80 
81  // Declare run-time constructor selection tables
82 
84  (
85  autoPtr,
86  fvPatch,
87  polyPatch,
88  (const polyPatch& patch, const fvBoundaryMesh& bm),
89  (patch, bm)
90  );
91 
92 
93  // Constructors
94 
95  //- Construct from polyPatch and fvBoundaryMesh
96  fvPatch(const polyPatch&, const fvBoundaryMesh&);
97 
98  //- Disallow default bitwise copy construction
99  fvPatch(const fvPatch&) = delete;
100 
101 
102  // Selectors
103 
104  //- Return a pointer to a new patch created on freestore from polyPatch
105  static autoPtr<fvPatch> New
106  (
107  const polyPatch&,
108  const fvBoundaryMesh&
109  );
110 
111 
112  //- Destructor
113  virtual ~fvPatch();
114 
115 
116  // Member Functions
117 
118  // Access
119 
120  //- Return the polyPatch
121  const polyPatch& patch() const
122  {
123  return polyPatch_;
124  }
125 
126  //- Return name
127  virtual const word& name() const
128  {
129  return polyPatch_.name();
130  }
131 
132  //- Return start label of this patch in the polyMesh face list
133  virtual label start() const
134  {
135  return polyPatch_.start();
136  }
137 
138  //- Return size
139  virtual label size() const
140  {
141  return polyPatch_.size();
142  }
143 
144  //- Return true if this patch is coupled
145  virtual bool coupled() const
146  {
147  return polyPatch_.coupled();
148  }
149 
150  //- Return true if the given type is a constraint type
151  static bool constraintType(const word& pt);
152 
153  //- Return a list of all the constraint patch types
154  static wordList constraintTypes();
155 
156  //- Return the index of this patch in the fvBoundaryMesh
157  label index() const
158  {
159  return polyPatch_.index();
160  }
161 
162  //- Return boundaryMesh reference
163  const fvBoundaryMesh& boundaryMesh() const
164  {
165  return boundaryMesh_;
166  }
167 
168  //- Return the local object registry
169  const objectRegistry& db() const;
170 
171  //- Slice list to patch
172  template<class T>
173  const typename List<T>::subList patchSlice(const List<T>& l) const
174  {
175  return typename List<T>::subList(l, size(), start());
176  }
177 
178  //- Return faceCells
179  virtual const labelUList& faceCells() const;
180 
181 
182  // Access functions for geometrical data
183 
184  //- Return face centres
185  const vectorField& Cf() const;
186 
187  //- Return neighbour cell centres
188  tmp<vectorField> Cn() const;
189 
190  //- Return face area vectors
191  const vectorField& Sf() const;
192 
193  //- Return face area magnitudes
194  const scalarField& magSf() const;
195 
196  //- Return face normals
197  tmp<vectorField> nf() const;
198 
199  //- Return cell-centre to face-centre vector
200  // except for coupled patches for which the cell-centre
201  // to coupled-cell-centre vector is returned
202  virtual tmp<vectorField> delta() const;
203 
204  //- Return the fraction of the poly-face that each fv-face in this
205  // patch covers. Will be equal to one, unless the mesh is
206  // non-conformal.
208 
209 
210  // Access functions for demand driven data
211 
212  //- Make patch weighting factors
213  virtual void makeWeights(scalarField&) const;
214 
215  //- Return patch weighting factors
216  const scalarField& weights() const;
217 
218  //- Return the face - cell distance coefficient
219  // except for coupled patches for which the cell-centre
220  // to coupled-cell-centre distance coefficient is returned
221  virtual const scalarField& deltaCoeffs() const;
222 
223 
224  // Evaluation functions
225 
226  //- Return given internal field next to patch as patch field
227  template<class Type>
229 
230  //- Return given internal field next to patch as patch field
231  template<class Type>
232  void patchInternalField(const UList<Type>&, Field<Type>&) const;
233 
234  //- Return the corresponding patchField of the named field
235  template<class GeometricField, class Type>
236  const typename GeometricField::Patch& patchField
237  (
238  const GeometricField&
239  ) const;
240 
241  //- Return the corresponding patchField reference of the named field
242  template<class GeometricField, class Type>
244  (
246  ) const;
247 
248  //- Lookup and return the patchField of the named field from the
249  // local objectRegistry.
250  template<class GeometricField, class Type>
252  (
253  const word& name
254  ) const;
255 
256 
257  // Member Operators
258 
259  //- Disallow default bitwise assignment
260  void operator=(const fvPatch&) = delete;
261 };
262 
263 
264 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
265 
266 } // End namespace Foam
267 
268 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
269 
270 #ifdef NoRepository
271  #include "fvPatchTemplates.C"
272 #endif
273 
274 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
275 
276 #endif
277 
278 // ************************************************************************* //
Generic GeometricField class.
PatchField< Type > Patch
Type of the patch field of which the Boundary is composed.
SubList< T > subList
Declare type of subList.
Definition: List.H:195
A List obtained as a section of another List.
Definition: SubList.H:56
label size() const
Return the number of elements in the UList.
Definition: UListI.H:311
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
Foam::fvBoundaryMesh.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
static wordList constraintTypes()
Return a list of all the constraint patch types.
Definition: fvPatch.C:67
virtual label size() const
Return size.
Definition: fvPatch.H:138
virtual ~fvPatch()
Destructor.
Definition: fvPatch.C:55
void operator=(const fvPatch &)=delete
Disallow default bitwise assignment.
virtual bool coupled() const
Return true if this patch is coupled.
Definition: fvPatch.H:144
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:156
virtual label start() const
Return start label of this patch in the polyMesh face list.
Definition: fvPatch.H:132
static autoPtr< fvPatch > New(const polyPatch &, const fvBoundaryMesh &)
Return a pointer to a new patch created on freestore from polyPatch.
Definition: fvPatchNew.C:32
const objectRegistry & db() const
Return the local object registry.
Definition: fvPatch.C:93
tmp< vectorField > Cn() const
Return neighbour cell centres.
Definition: fvPatch.C:111
const polyPatch & patch() const
Return the polyPatch.
Definition: fvPatch.H:120
virtual const word & name() const
Return name.
Definition: fvPatch.H:126
declareRunTimeSelectionTable(autoPtr, fvPatch, polyPatch,(const polyPatch &patch, const fvBoundaryMesh &bm),(patch, bm))
tmp< scalarField > polyFaceFraction() const
Return the fraction of the poly-face that each fv-face in this.
Definition: fvPatch.C:156
virtual tmp< vectorField > delta() const
Return cell-centre to face-centre vector.
Definition: fvPatch.C:148
static bool constraintType(const word &pt)
Return true if the given type is a constraint type.
Definition: fvPatch.C:61
const scalarField & magSf() const
Return face area magnitudes.
Definition: fvPatch.C:142
const scalarField & weights() const
Return patch weighting factors.
Definition: fvPatch.C:182
const List< T >::subList patchSlice(const List< T > &l) const
Slice list to patch.
Definition: fvPatch.H:172
fvPatch(const polyPatch &, const fvBoundaryMesh &)
Construct from polyPatch and fvBoundaryMesh.
Definition: fvPatch.C:46
tmp< vectorField > nf() const
Return face normals.
Definition: fvPatch.C:130
tmp< Field< Type > > patchInternalField(const UList< Type > &) const
Return given internal field next to patch as patch field.
const vectorField & Cf() const
Return face centres.
Definition: fvPatch.C:105
const fvBoundaryMesh & boundaryMesh() const
Return boundaryMesh reference.
Definition: fvPatch.H:162
virtual const scalarField & deltaCoeffs() const
Return the face - cell distance coefficient.
Definition: fvPatch.C:176
TypeName(polyPatch::typeName_())
Runtime type information.
const vectorField & Sf() const
Return face area vectors.
Definition: fvPatch.C:136
virtual const labelUList & faceCells() const
Return faceCells.
Definition: fvPatch.C:99
const GeometricField::Patch & lookupPatchField(const word &name) const
Lookup and return the patchField of the named field from the.
Registry of regIOobjects.
label index() const
Return the index of this patch in the boundaryMesh.
const word & name() const
Return name.
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:70
virtual bool coupled() const
Return true if this patch is geometrically coupled (i.e. faces and.
Definition: polyPatch.H:290
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:280
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
Namespace for OpenFOAM.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Specialisations of Field<T> for scalar, vector and tensor.
Macros to ease declaration of run-time selection tables.
Basic run-time type information using word as the type's name. Used to enhance the standard RTTI to c...