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 Copyright (C) 2020-2021 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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\*---------------------------------------------------------------------------*/
28
29#include "fvPatch.H"
31#include "fvBoundaryMesh.H"
32#include "fvMesh.H"
33#include "primitiveMesh.H"
34#include "volFields.H"
35#include "surfaceFields.H"
36
37// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38
39namespace Foam
40{
44}
45
46
47// * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
48
50{
51 const fvMesh* meshptr = isA<fvMesh>(p.boundaryMesh().mesh());
52
53 if (!meshptr)
54 {
56 << "The polyPatch is not attached to a base fvMesh" << nl
57 << exit(FatalError);
58 }
59
60 return meshptr->boundary()[p.index()];
61}
62
63
64// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
65
67:
68 polyPatch_(p),
69 boundaryMesh_(bm)
70{}
71
72
73// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
74
76{} // fvBoundaryMesh was forward declared
77
78
79// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
80
82{
83 return
84 (
87 );
88}
89
90
92{
93 const auto& cnstrTable = *polyPatchConstructorTablePtr_;
94
95 wordList cTypes(cnstrTable.size());
96
97 label i = 0;
98
99 forAllConstIters(cnstrTable, iter)
100 {
101 if (constraintType(iter.key()))
102 {
103 cTypes[i++] = iter.key();
104 }
105 }
106
107 cTypes.setSize(i);
108
109 return cTypes;
110}
111
112
114{
115 return polyPatch_.faceCells();
116}
117
118
120{
121 return boundaryMesh().mesh().Cf().boundaryField()[index()];
122}
123
124
126{
127 auto tcc = tmp<vectorField>::New(size());
128 auto& cc = tcc.ref();
129
130 const labelUList& faceCells = this->faceCells();
131
132 // get reference to global cell centres
133 const vectorField& gcc = boundaryMesh().mesh().cellCentres();
134
135 forAll(faceCells, facei)
136 {
137 cc[facei] = gcc[faceCells[facei]];
138 }
139
140 return tcc;
141}
142
143
145{
146 return Sf()/magSf();
147}
148
149
151{
152 return boundaryMesh().mesh().Sf().boundaryField()[index()];
153}
154
155
157{
158 return boundaryMesh().mesh().magSf().boundaryField()[index()];
159}
160
161
163{
164 // Use patch-normal delta for all non-coupled BCs
165 const vectorField nHat(nf());
166 return nHat*(nHat & (Cf() - Cn()));
167}
168
169
171{
172 w = 1.0;
173}
174
175
177{}
178
179
181{}
182
183
185{}
186
187
189{}
190
191
193{}
194
195
197{
198 return boundaryMesh().mesh().deltaCoeffs().boundaryField()[index()];
199}
200
201
203{
204 return boundaryMesh().mesh().weights().boundaryField()[index()];
205}
206
207
208// ************************************************************************* //
bool found
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
const Mesh & mesh() const
Return mesh.
void setSize(const label n)
Alias for resize()
Definition: List.H:218
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface....
Definition: boundaryMesh.H:63
const bMesh & mesh() const
Definition: boundaryMesh.H:206
virtual const word & constraintType() const
Return the constraint type this pointPatch implements.
Smooth ATC in cells next to a set of patches supplied by type.
Definition: faceCells.H:59
Foam::fvBoundaryMesh.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:712
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:82
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:71
static wordList constraintTypes()
Return a list of all the constraint patch types.
Definition: fvPatch.C:91
virtual ~fvPatch()
Destructor.
Definition: fvPatch.C:75
static const fvPatch & lookupPatch(const polyPatch &p)
Lookup the polyPatch index on corresponding fvMesh.
Definition: fvPatch.C:49
virtual void makeWeights(scalarField &) const
Make patch weighting factors.
Definition: fvPatch.C:170
virtual void makeNonOrthoDeltaCoeffs(scalarField &) const
Correct patch non-ortho deltaCoeffs.
Definition: fvPatch.C:180
virtual void makeDeltaCoeffs(scalarField &) const
Correct patch deltaCoeffs.
Definition: fvPatch.C:176
virtual void movePoints()
Correct patches after moving points.
Definition: fvPatch.C:192
tmp< vectorField > Cn() const
Return neighbour cell centres.
Definition: fvPatch.C:125
virtual tmp< vectorField > delta() const
Definition: fvPatch.C:162
virtual void initMovePoints()
Initialise the patches for moving points.
Definition: fvPatch.C:188
virtual void makeNonOrthoCorrVectors(vectorField &) const
Correct patch non-ortho correction vectors.
Definition: fvPatch.C:184
const scalarField & magSf() const
Return face area magnitudes.
Definition: fvPatch.C:156
const scalarField & weights() const
Return patch weighting factors.
Definition: fvPatch.C:202
tmp< vectorField > nf() const
Return face normals.
Definition: fvPatch.C:144
const vectorField & Cf() const
Return face centres.
Definition: fvPatch.C:119
const scalarField & deltaCoeffs() const
Definition: fvPatch.C:196
const vectorField & Sf() const
Return face area vectors.
Definition: fvPatch.C:150
virtual const labelUList & faceCells() const
Return faceCells.
Definition: fvPatch.C:113
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:75
A class for managing temporary objects.
Definition: tmp.H:65
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
volScalarField & p
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Namespace for OpenFOAM.
error FatalError
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
#define defineRunTimeSelectionTable(baseType, argNames)
Define run-time selection table.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
#define forAllConstIters(container, iter)
Iterate across all elements of the container object with const access.
Definition: stdFoam.H:278
Foam::surfaceFields.