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