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-2023 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 protected:
76 
77  // Protected Member Functions
78 
79  //- Make patch weighting factors
80  virtual void makeWeights(scalarField&) const;
81 
82  //- Initialise the patches for moving points
83  virtual void initMovePoints();
84 
85  //- Correct patches after moving points
86  virtual void movePoints();
87 
88 
89 public:
90 
92 
93  friend class fvBoundaryMesh;
94  friend class surfaceInterpolation;
95 
96  //- Runtime type information
97  TypeName(polyPatch::typeName_());
98 
99 
100  // Declare run-time constructor selection tables
101 
103  (
104  autoPtr,
105  fvPatch,
106  polyPatch,
107  (const polyPatch& patch, const fvBoundaryMesh& bm),
108  (patch, bm)
109  );
110 
111 
112  // Constructors
113 
114  //- Construct from polyPatch and fvBoundaryMesh
115  fvPatch(const polyPatch&, const fvBoundaryMesh&);
116 
117  //- Disallow default bitwise copy construction
118  fvPatch(const fvPatch&);
119 
120 
121  // Selectors
122 
123  //- Return a pointer to a new patch created on freestore from polyPatch
124  static autoPtr<fvPatch> New
125  (
126  const polyPatch&,
127  const fvBoundaryMesh&
128  );
129 
130 
131  //- Destructor
132  virtual ~fvPatch();
133 
134 
135  // Member Functions
136 
137  // Access
138 
139  //- Return the polyPatch
140  const polyPatch& patch() const
141  {
142  return polyPatch_;
143  }
144 
145  //- Return name
146  virtual const word& name() const
147  {
148  return polyPatch_.name();
149  }
150 
151  //- Return start label of this patch in the polyMesh face list
152  virtual label start() const
153  {
154  return polyPatch_.start();
155  }
156 
157  //- Return size
158  virtual label size() const
159  {
160  return polyPatch_.size();
161  }
162 
163  //- Return true if this patch is coupled
164  virtual bool coupled() const
165  {
166  return polyPatch_.coupled();
167  }
168 
169  //- Return true if the given type is a constraint type
170  static bool constraintType(const word& pt);
171 
172  //- Return a list of all the constraint patch types
173  static wordList constraintTypes();
174 
175  //- Return the index of this patch in the fvBoundaryMesh
176  label index() const
177  {
178  return polyPatch_.index();
179  }
180 
181  //- Return boundaryMesh reference
182  const fvBoundaryMesh& boundaryMesh() const
183  {
184  return boundaryMesh_;
185  }
186 
187  //- Return the local object registry
188  const objectRegistry& db() const;
189 
190  //- Slice list to patch
191  template<class T>
192  const typename List<T>::subList patchSlice(const List<T>& l) const
193  {
194  return typename List<T>::subList(l, size(), start());
195  }
196 
197  //- Return faceCells
198  virtual const labelUList& faceCells() const;
199 
200 
201  // Access functions for geometrical data
202 
203  //- Return face centres
204  const vectorField& Cf() const;
205 
206  //- Return neighbour cell centres
207  tmp<vectorField> Cn() const;
208 
209  //- Return face area vectors
210  const vectorField& Sf() const;
211 
212  //- Return face area magnitudes
213  const scalarField& magSf() const;
214 
215  //- Return face normals
216  tmp<vectorField> nf() const;
217 
218  //- Return cell-centre to face-centre vector
219  // except for coupled patches for which the cell-centre
220  // to coupled-cell-centre vector is returned
221  virtual tmp<vectorField> delta() const;
222 
223 
224  // Access functions for demand driven data
225 
226  //- Return patch weighting factors
227  const scalarField& weights() const;
228 
229  //- Return the face - cell distance coefficient
230  // except for coupled patches for which the cell-centre
231  // to coupled-cell-centre distance coefficient is returned
232  virtual const scalarField& deltaCoeffs() const;
233 
234 
235  // Evaluation functions
236 
237  //- Return given internal field next to patch as patch field
238  template<class Type>
240 
241  //- Return given internal field next to patch as patch field
242  template<class Type>
243  void patchInternalField(const UList<Type>&, Field<Type>&) const;
244 
245  //- Return the corresponding patchField of the named field
246  template<class GeometricField, class Type>
247  const typename GeometricField::Patch& patchField
248  (
249  const GeometricField&
250  ) const;
251 
252  //- Return the corresponding patchField reference of the named field
253  template<class GeometricField, class Type>
255  (
257  ) const;
258 
259  //- Lookup and return the patchField of the named field from the
260  // local objectRegistry.
261  template<class GeometricField, class Type>
263  (
264  const word& name
265  ) const;
266 
267 
268  // Member Operators
269 
270  //- Disallow default bitwise assignment
271  void operator=(const fvPatch&);
272 };
273 
274 
275 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
276 
277 } // End namespace Foam
278 
279 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
280 
281 #ifdef NoRepository
282  #include "fvPatchTemplates.C"
283 #endif
284 
285 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
286 
287 #endif
288 
289 // ************************************************************************* //
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
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:157
virtual ~fvPatch()
Destructor.
Definition: fvPatch.C:55
virtual bool coupled() const
Return true if this patch is coupled.
Definition: fvPatch.H:163
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:156
label index() const
Return the index of this patch in the fvBoundaryMesh.
Definition: fvPatch.H:175
virtual label start() const
Return start label of this patch in the polyMesh face list.
Definition: fvPatch.H:151
static autoPtr< fvPatch > New(const polyPatch &, const fvBoundaryMesh &)
Return a pointer to a new patch created on freestore from polyPatch.
Definition: fvPatchNew.C:32
virtual void movePoints()
Correct patches after moving points.
Definition: fvPatch.C:166
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:139
virtual const word & name() const
Return name.
Definition: fvPatch.H:145
declareRunTimeSelectionTable(autoPtr, fvPatch, polyPatch,(const polyPatch &patch, const fvBoundaryMesh &bm),(patch, bm))
virtual tmp< vectorField > delta() const
Return cell-centre to face-centre vector.
Definition: fvPatch.C:148
virtual void initMovePoints()
Initialise the patches for moving points.
Definition: fvPatch.C:162
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:176
const List< T >::subList patchSlice(const List< T > &l) const
Slice list to patch.
Definition: fvPatch.H:191
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.
fvBoundaryMesh BoundaryMesh
Definition: fvPatch.H:90
const vectorField & Cf() const
Return face centres.
Definition: fvPatch.C:105
const fvBoundaryMesh & boundaryMesh() const
Return boundaryMesh reference.
Definition: fvPatch.H:181
virtual const scalarField & deltaCoeffs() const
Return the face - cell distance coefficient.
Definition: fvPatch.C:170
void operator=(const fvPatch &)
Disallow default bitwise assignment.
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
Cell to surface interpolation scheme. Included in fvMesh.
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.