fvPatch.C
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-2022 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 \*---------------------------------------------------------------------------*/
25 
26 #include "fvPatch.H"
28 #include "fvBoundaryMesh.H"
29 #include "fvMesh.H"
30 #include "primitiveMesh.H"
31 #include "volFields.H"
32 #include "surfaceFields.H"
33 
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38  defineTypeNameAndDebug(fvPatch, 0);
39  defineRunTimeSelectionTable(fvPatch, polyPatch);
40  addToRunTimeSelectionTable(fvPatch, fvPatch, polyPatch);
41 }
42 
43 
44 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
45 
47 :
48  polyPatch_(p),
49  boundaryMesh_(bm)
50 {}
51 
52 
53 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
54 
56 {}
57 
58 
59 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
60 
62 {
64 }
65 
66 
68 {
69  wordList cTypes(polyPatchConstructorTablePtr_->size());
70 
71  label i = 0;
72 
73  for
74  (
75  polyPatchConstructorTable::iterator cstrIter =
76  polyPatchConstructorTablePtr_->begin();
77  cstrIter != polyPatchConstructorTablePtr_->end();
78  ++cstrIter
79  )
80  {
81  if (constraintType(cstrIter.key()))
82  {
83  cTypes[i++] = cstrIter.key();
84  }
85  }
86 
87  cTypes.setSize(i);
88 
89  return cTypes;
90 }
91 
92 
94 {
95  return boundaryMesh().mesh();
96 }
97 
98 
100 {
101  return polyPatch_.faceCells();
102 }
103 
104 
106 {
107  return boundaryMesh().mesh().Cf().boundaryField()[index()];
108 }
109 
110 
112 {
113  tmp<vectorField> tcc(new vectorField(size()));
114  vectorField& cc = tcc.ref();
115 
116  const labelUList& faceCells = this->faceCells();
117 
118  // get reference to global cell centres
119  const vectorField& gcc = boundaryMesh().mesh().cellCentres();
120 
121  forAll(faceCells, facei)
122  {
123  cc[facei] = gcc[faceCells[facei]];
124  }
125 
126  return tcc;
127 }
128 
129 
131 {
132  return Sf()/magSf();
133 }
134 
135 
137 {
138  return boundaryMesh().mesh().Sf().boundaryField()[index()];
139 }
140 
141 
143 {
144  return boundaryMesh().mesh().magSf().boundaryField()[index()];
145 }
146 
147 
149 {
150  // Use patch-normal delta for all non-coupled BCs
151  const vectorField nHat(nf());
152  return nHat*(nHat & (Cf() - Cn()));
153 }
154 
155 
157 {
158  w = 1.0;
159 }
160 
161 
163 {}
164 
165 
167 {}
168 
169 
171 {
172  return boundaryMesh().mesh().deltaCoeffs().boundaryField()[index()];
173 }
174 
175 
177 {
178  return boundaryMesh().mesh().weights().boundaryField()[index()];
179 }
180 
181 
182 // ************************************************************************* //
Foam::surfaceFields.
const fvBoundaryMesh & boundaryMesh() const
Return boundaryMesh reference.
Definition: fvPatch.H:181
virtual void initMovePoints()
Initialise the patches for moving points.
Definition: fvPatch.C:162
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
const surfaceVectorField & Sf() const
Return cell face area vectors.
const scalarField & magSf() const
Return face area magnitudes.
Definition: fvPatch.C:142
label index() const
Return the index of this patch in the fvBoundaryMesh.
Definition: fvPatch.H:175
const surfaceVectorField & Cf() const
Return face centres.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
tmp< vectorField > nf() const
Return face normals.
Definition: fvPatch.C:130
virtual ~fvPatch()
Destructor.
Definition: fvPatch.C:55
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
static bool constraintType(const word &pt)
Return true if the given type is a constraint type.
Definition: fvPatch.C:61
Macros for easy insertion into run-time selection tables.
const surfaceScalarField & weights() const
Return reference to linear difference weighting factors.
const labelUList & faceCells() const
Return face-cell addressing.
Definition: polyPatch.C:340
const fvMesh & mesh() const
Return the mesh reference.
A class for handling words, derived from string.
Definition: word.H:59
virtual const labelUList & faceCells() const
Return faceCells.
Definition: fvPatch.C:99
const scalarField & weights() const
Return patch weighting factors.
Definition: fvPatch.C:176
virtual void makeWeights(scalarField &) const
Make patch weighting factors.
Definition: fvPatch.C:156
const vectorField & Cf() const
Return face centres.
Definition: fvPatch.C:105
const vectorField & cellCentres() const
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: HashTable.H:60
virtual label size() const
Return size.
Definition: fvPatch.H:157
const surfaceScalarField & magSf() const
Return cell face area magnitudes.
const objectRegistry & db() const
Return the local object registry.
Definition: fvPatch.C:93
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
defineTypeNameAndDebug(combustionModel, 0)
const vectorField & Sf() const
Return face area vectors.
Definition: fvPatch.C:136
const scalarField & deltaCoeffs() const
Return the face - cell distance coeffient.
Definition: fvPatch.C:170
Foam::fvBoundaryMesh.
tmp< vectorField > Cn() const
Return neighbour cell centres.
Definition: fvPatch.C:111
Field< vector > vectorField
Specialisation of Field<T> for vector.
volScalarField & p
A class for managing temporary objects.
Definition: PtrList.H:53
Registry of regIOobjects.
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
virtual tmp< vectorField > delta() const
Return cell-centre to face-centre vector.
Definition: fvPatch.C:148
const surfaceScalarField & deltaCoeffs() const
Return reference to cell-centre difference coefficients.
virtual void movePoints()
Correct patches after moving points.
Definition: fvPatch.C:166
fvPatch(const polyPatch &, const fvBoundaryMesh &)
Construct from polyPatch and fvBoundaryMesh.
Definition: fvPatch.C:46
Namespace for OpenFOAM.
static wordList constraintTypes()
Return a list of all the constraint patch types.
Definition: fvPatch.C:67