Loading...
Searching...
No Matches
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-2022 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"
37// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38
39namespace Foam
40{
43 addToRunTimeSelectionTable(fvPatch, fvPatch, polyPatch);
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 }
60 return meshptr->boundary()[p.index()];
61}
62
63
64// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
65
66Foam::fvPatch::fvPatch(const polyPatch& p, const fvBoundaryMesh& bm)
67:
68 polyPatch_(p),
69 boundaryMesh_(bm)
70{}
71
72
73// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
76{} // fvBoundaryMesh was forward declared
77
78
79// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
80
81bool Foam::fvPatch::constraintType(const word& patchType)
82{
83 return
84 (
85 !patchType.empty()
88 );
89}
90
91
93{
94 const auto& cnstrTable = *polyPatchConstructorTablePtr_;
95
96 wordList cTypes(cnstrTable.size());
97
98 label i = 0;
99
100 forAllConstIters(cnstrTable, iter)
101 {
102 if (constraintType(iter.key()))
103 {
104 cTypes[i++] = iter.key();
105 }
106 }
108 cTypes.resize(i);
109
110 return cTypes;
111}
112
115{
116 return polyPatch_.faceCells();
117}
118
121{
122 return boundaryMesh().mesh().Cf().boundaryField()[index()];
123}
124
129}
130
133{
134 return boundaryMesh().mesh().Sf().boundaryField()[index()];
135}
136
139{
140 return boundaryMesh().mesh().magSf().boundaryField()[index()];
141}
142
147}
148
151{
152 return Sf()/magSf();
153}
154
155
158 // Use patch-normal delta for all non-coupled BCs
159 const vectorField nHat(nf());
160 return nHat*(nHat & (Cf() - Cn()));
161}
162
165{
166 w = 1.0;
168
169
172
173
176
177
180
181
184
185
187{}
188
191{
192 return boundaryMesh().mesh().deltaCoeffs().boundaryField()[index()];
193}
194
195
197{
198 return boundaryMesh().mesh().weights().boundaryField()[index()];
199}
200
201
202// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
void resize(const label len)
Adjust allocated size of list.
Definition ListI.H:153
bool found(const T &val, label pos=0) const
Same as contains().
Definition UList.H:983
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface....
const bMesh & mesh() const
A fvBoundaryMesh is a fvPatch list with a reference to the associated fvMesh, with additional search ...
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
const fvBoundaryMesh & boundary() const noexcept
Return reference to boundary mesh.
Definition fvMesh.H:395
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
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:85
virtual ~fvPatch()
Destructor.
Definition fvPatch.C:68
static const fvPatch & lookupPatch(const polyPatch &p)
Lookup the polyPatch index on corresponding fvMesh.
Definition fvPatch.C:42
virtual void makeWeights(scalarField &) const
Make patch weighting factors.
Definition fvPatch.C:157
virtual void makeNonOrthoDeltaCoeffs(scalarField &) const
Correct patch non-ortho deltaCoeffs.
Definition fvPatch.C:167
virtual void makeDeltaCoeffs(scalarField &) const
Correct patch deltaCoeffs.
Definition fvPatch.C:163
virtual void movePoints()
Correct patches after moving points.
Definition fvPatch.C:179
tmp< vectorField > Cn() const
Return neighbour cell centres.
Definition fvPatch.C:119
label index() const noexcept
The index of this patch in the boundary mesh.
Definition fvPatch.H:218
const fvBoundaryMesh & boundaryMesh() const noexcept
Return boundaryMesh reference.
Definition fvPatch.H:268
virtual tmp< vectorField > delta() const
Return cell-centre to face-centre vector except for coupled patches for which the cell-centre to coup...
Definition fvPatch.C:149
virtual void initMovePoints()
Initialise the patches for moving points.
Definition fvPatch.C:175
virtual void makeNonOrthoCorrVectors(vectorField &) const
Correct patch non-ortho correction vectors.
Definition fvPatch.C:171
const scalarField & magSf() const
Return face area magnitudes, like the fvMesh::magSf() method.
Definition fvPatch.C:131
const scalarField & weights() const
Return patch weighting factors.
Definition fvPatch.C:189
friend class fvBoundaryMesh
Definition fvPatch.H:140
tmp< vectorField > unitSf() const
Return face unit normals, like the fvMesh::unitSf() method. Same as nf().
Definition fvPatch.C:137
tmp< vectorField > nf() const
Return face unit normals, like the fvMesh::unitSf() method Same as unitSf().
Definition fvPatch.C:143
const vectorField & Cf() const
Return face centres.
Definition fvPatch.C:113
void patchInternalField(const UList< Type > &internalData, const labelUList &addressing, UList< Type > &pfld) const
Extract internal field next to patch using specified addressing.
const scalarField & deltaCoeffs() const
Return the face - cell distance coefficient except for coupled patches for which the cell-centre to c...
Definition fvPatch.C:183
const vectorField & Sf() const
Return face area vectors, like the fvMesh::Sf() method.
Definition fvPatch.C:125
static bool constraintType(const word &patchType)
Return true if the given type is a constraint type.
Definition fvPatch.C:74
virtual const labelUList & faceCells() const
Return faceCells.
Definition fvPatch.C:107
A patch is a list of labels that address the faces in the global face list.
Definition polyPatch.H:73
A class for managing temporary objects.
Definition tmp.H:75
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
volScalarField & p
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
Namespace for OpenFOAM.
List< word > wordList
List of word.
Definition fileName.H:60
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
Field< vector > vectorField
Specialisation of Field<T> for vector.
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
UList< label > labelUList
A UList of labels.
Definition UList.H:75
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
#define defineRunTimeSelectionTable(baseType, argNames)
Define run-time selection table.
#define forAllConstIters(container, iter)
Iterate across all elements of the container object with const access.
Definition stdFoam.H:235
Foam::surfaceFields.