Loading...
Searching...
No Matches
movingWallVelocityFvPatchVectorField.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) 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
30#include "volFields.H"
31#include "surfaceFields.H"
32#include "fvcMeshPhi.H"
34
35// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
36
39(
40 const fvPatch& p,
42)
43:
44 fixedValueFvPatchVectorField(p, iF)
45{}
46
47
50(
51 const fvPatch& p,
54)
55:
56 fixedValueFvPatchVectorField(p, iF, dict)
57{}
58
59
62(
64 const fvPatch& p,
66 const fvPatchFieldMapper& mapper
67)
68:
69 fixedValueFvPatchVectorField(ptf, p, iF, mapper)
70{}
71
72
75(
77)
78:
79 fixedValueFvPatchVectorField(mwvpvf)
80{}
81
82
85(
88)
89:
90 fixedValueFvPatchVectorField(mwvpvf, iF)
91{}
92
93
94// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
95
98{
99 const fvMesh& mesh = internalField().mesh();
100 const fvPatch& p = patch();
101 const polyPatch& pp = p.patch();
102 const pointField& oldPoints = mesh.oldPoints();
103
104 vectorField oldFc(pp.size());
105
106 forAll(oldFc, i)
107 {
108 oldFc[i] = pp[i].centre(oldPoints);
109 }
110
111 const scalar deltaT = mesh.time().deltaTValue();
112
113 const vectorField Up((pp.faceCentres() - oldFc)/deltaT);
114
115 const volVectorField& U =
116 static_cast<const volVectorField&>(internalField());
117
118 scalarField phip
119 (
120 p.patchField<surfaceScalarField>(fvc::meshPhi(U))
121 );
122
123 const vectorField n(p.nf());
124 const scalarField& magSf = p.magSf();
125 tmp<scalarField> Un = phip/(magSf + VSMALL);
126
127 return (Up + n*(Un - (n & Up)));
128}
129
130
132{
133 if (updated())
134 {
135 return;
136 }
137
138 const fvMesh& mesh = internalField().mesh();
139
140 if (mesh.moving())
141 {
142 const vectorField uwall(Uwall());
144 }
145
146 fixedValueFvPatchVectorField::updateCoeffs();
147}
148
149
151{
155
156
157// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
158
159namespace Foam
160{
162 (
164 movingWallVelocityFvPatchVectorField
165 );
166}
167
168// ************************************************************************* //
label n
Macros for easy insertion into run-time selection tables.
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
void operator=(const Field< vector > &)
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
A FieldMapper for finite-volume patch fields.
virtual void write(Ostream &) const
Write.
void writeValueEntry(Ostream &os) const
Write *this field as a "value" entry.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition fvPatch.H:71
This boundary condition provides a velocity condition for cases with moving walls.
movingWallVelocityFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
tmp< vectorField > Uwall() const
Return wall velocity field.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
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
U
Definition pEqn.H:72
volScalarField & p
dynamicFvMesh & mesh
OBJstream os(runTime.globalPath()/outputName)
#define makePatchTypeField(PatchTypeField, typePatchTypeField)
Define a concrete fvPatchField type and add to run-time tables Example, (fvPatchScalarField,...
Calculate the mesh motion flux and convert fluxes from absolute to relative and back.
tmp< surfaceScalarField > meshPhi(const volVectorField &U)
Definition fvcMeshPhi.C:29
Namespace for OpenFOAM.
GeometricField< vector, fvPatchField, volMesh > volVectorField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Field< vector > vectorField
Specialisation of Field<T> for vector.
vectorField pointField
pointField is a vectorField.
fvPatchField< vector > fvPatchVectorField
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
Foam::surfaceFields.