fvPatch.C
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-2016 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "fvPatch.H"
30 #include "fvBoundaryMesh.H"
31 #include "fvMesh.H"
32 #include "primitiveMesh.H"
33 #include "volFields.H"
34 #include "surfaceFields.H"
35 
36 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
40  defineTypeNameAndDebug(fvPatch, 0);
41  defineRunTimeSelectionTable(fvPatch, polyPatch);
42  addToRunTimeSelectionTable(fvPatch, fvPatch, polyPatch);
43 }
44 
45 
46 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
47 
48 Foam::fvPatch::fvPatch(const polyPatch& p, const fvBoundaryMesh& bm)
49 :
50  polyPatch_(p),
51  boundaryMesh_(bm)
52 {}
53 
54 
55 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
56 
58 {}
59 
60 
61 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
62 
64 {
66 }
67 
68 
70 {
71  wordList cTypes(polyPatchConstructorTablePtr_->size());
72 
73  label i = 0;
74 
75  forAllConstIters(*polyPatchConstructorTablePtr_, cstrIter)
76  {
77  if (constraintType(cstrIter.key()))
78  {
79  cTypes[i++] = cstrIter.key();
80  }
81  }
82 
83  cTypes.setSize(i);
84 
85  return cTypes;
86 }
87 
88 
90 {
91  return polyPatch_.faceCells();
92 }
93 
94 
96 {
97  return boundaryMesh().mesh().Cf().boundaryField()[index()];
98 }
99 
100 
102 {
103  tmp<vectorField> tcc(new vectorField(size()));
104  vectorField& cc = tcc.ref();
105 
106  const labelUList& faceCells = this->faceCells();
107 
108  // get reference to global cell centres
109  const vectorField& gcc = boundaryMesh().mesh().cellCentres();
110 
111  forAll(faceCells, facei)
112  {
113  cc[facei] = gcc[faceCells[facei]];
114  }
115 
116  return tcc;
117 }
118 
119 
121 {
122  return Sf()/magSf();
123 }
124 
125 
127 {
128  return boundaryMesh().mesh().Sf().boundaryField()[index()];
129 }
130 
131 
133 {
134  return boundaryMesh().mesh().magSf().boundaryField()[index()];
135 }
136 
137 
139 {
140  // Use patch-normal delta for all non-coupled BCs
141  const vectorField nHat(nf());
142  return nHat*(nHat & (Cf() - Cn()));
143 }
144 
145 
147 {
148  w = 1.0;
149 }
150 
151 
153 {}
154 
155 
157 {}
158 
159 
161 {
162  return boundaryMesh().mesh().deltaCoeffs().boundaryField()[index()];
163 }
164 
165 
167 {
168  return boundaryMesh().mesh().weights().boundaryField()[index()];
169 }
170 
171 
172 // ************************************************************************* //
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
Foam::fvPatchField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: volSurfaceMapping.H:50
volFields.H
Foam::boundaryMesh::mesh
const bMesh & mesh() const
Definition: boundaryMesh.H:205
Foam::fvPatch::Sf
const vectorField & Sf() const
Return face area vectors.
Definition: fvPatch.C:126
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:59
Foam::fvPatch::~fvPatch
virtual ~fvPatch()
Destructor.
Definition: fvPatch.C:57
Foam::defineRunTimeSelectionTable
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
Foam::fvPatch::delta
virtual tmp< vectorField > delta() const
Definition: fvPatch.C:138
primitiveMesh.H
surfaceFields.H
Foam::surfaceFields.
Foam::fvPatch::movePoints
virtual void movePoints()
Correct patches after moving points.
Definition: fvPatch.C:156
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
Foam::vectorField
Field< vector > vectorField
Specialisation of Field<T> for vector.
Definition: primitiveFieldsFwd.H:54
Foam::fvBoundaryMesh
Foam::fvBoundaryMesh.
Definition: fvBoundaryMesh.H:57
Foam::fvPatch::nf
tmp< vectorField > nf() const
Return face normals.
Definition: fvPatch.C:120
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:258
Foam::label
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:62
Foam::Field< vector >
Foam::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
Foam::fvPatch::initMovePoints
virtual void initMovePoints()
Initialise the patches for moving points.
Definition: fvPatch.C:152
Foam::fvPatch::constraintType
static bool constraintType(const word &pt)
Return true if the given type is a constraint type.
Definition: fvPatch.C:63
Foam::fvPatch::Cf
const vectorField & Cf() const
Return face centres.
Definition: fvPatch.C:95
fvBoundaryMesh.H
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
fvMesh.H
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
forAllConstIters
forAllConstIters(mixture.phases(), phase)
Definition: pEqn.H:28
Foam::fvPatch::weights
const scalarField & weights() const
Return patch weighting factors.
Definition: fvPatch.C:166
Foam::List< word >
Foam::fvPatch::deltaCoeffs
const scalarField & deltaCoeffs() const
Definition: fvPatch.C:160
Foam::UList< label >
Foam::fvPatch::faceCells
virtual const labelUList & faceCells() const
Return faceCells.
Definition: fvPatch.C:89
Foam::boundaryMesh
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface....
Definition: boundaryMesh.H:61
fvPatch.H
Foam::fvPatch::magSf
const scalarField & magSf() const
Return face area magnitudes.
Definition: fvPatch.C:132
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::faceCells
Smooth ATC in cells next to a set of patches supplied by type.
Definition: faceCells.H:56
Foam::fvPatch::makeWeights
virtual void makeWeights(scalarField &) const
Make patch weighting factors.
Definition: fvPatch.C:146
Foam::fvPatch::constraintTypes
static wordList constraintTypes()
Return a list of all the constraint patch types.
Definition: fvPatch.C:69
Foam::fvPatch::Cn
tmp< vectorField > Cn() const
Return neighbour cell centres.
Definition: fvPatch.C:101