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 | Copyright (C) 2011-2016 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 
27 #include "volFields.H"
28 #include "surfaceFields.H"
29 #include "fvcMeshPhi.H"
31 
32 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
33 
36 (
37  const fvPatch& p,
39 )
40 :
41  fixedValueFvPatchVectorField(p, iF)
42 {}
43 
44 
47 (
48  const fvPatch& p,
50  const dictionary& dict
51 )
52 :
53  fixedValueFvPatchVectorField(p, iF)
54 {
55  fvPatchVectorField::operator=(vectorField("value", dict, p.size()));
56 }
57 
58 
61 (
63  const fvPatch& p,
65  const fvPatchFieldMapper& mapper
66 )
67 :
68  fixedValueFvPatchVectorField(ptf, p, iF, mapper)
69 {}
70 
71 
74 (
76 )
77 :
78  fixedValueFvPatchVectorField(mwvpvf)
79 {}
80 
81 
84 (
87 )
88 :
89  fixedValueFvPatchVectorField(mwvpvf, iF)
90 {}
91 
92 
93 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
94 
96 {
97  if (updated())
98  {
99  return;
100  }
101 
102  const fvMesh& mesh = internalField().mesh();
103 
104  if (mesh.moving())
105  {
106  const fvPatch& p = patch();
107  const polyPatch& pp = p.patch();
108  const pointField& oldPoints = mesh.oldPoints();
109 
110  vectorField oldFc(pp.size());
111 
112  forAll(oldFc, i)
113  {
114  oldFc[i] = pp[i].centre(oldPoints);
115  }
116 
117  const scalar deltaT = mesh.time().deltaTValue();
118 
119  const vectorField Up((pp.faceCentres() - oldFc)/deltaT);
120 
121  const volVectorField& U =
122  static_cast<const volVectorField&>(internalField());
123 
124  scalarField phip
125  (
127  );
128 
129  const vectorField n(p.nf());
130  const scalarField& magSf = p.magSf();
131  tmp<scalarField> Un = phip/(magSf + VSMALL);
132 
133 
134  vectorField::operator=(Up + n*(Un - (n & Up)));
135  }
136 
137  fixedValueFvPatchVectorField::updateCoeffs();
138 }
139 
140 
142 {
144  writeEntry("value", os);
145 }
146 
147 
148 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
149 
150 namespace Foam
151 {
153  (
156  );
157 }
158 
159 // ************************************************************************* //
Foam::surfaceFields.
const scalarField & magSf() const
Return face area magnitudes.
Definition: fvPatch.C:136
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
U
Definition: pEqn.H:83
This boundary condition provides a velocity condition for cases with moving walls.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:65
volVectorField vectorField(fieldObject, mesh)
const polyPatch & patch() const
Return the polyPatch.
Definition: fvPatch.H:143
Macros for easy insertion into run-time selection tables.
bool moving() const
Is mesh moving.
Definition: polyMesh.H:493
makePatchTypeField(fvPatchVectorField, SRFFreestreamVelocityFvPatchVectorField)
dynamicFvMesh & mesh
virtual const pointField & oldPoints() const
Return old points for mesh motion.
Definition: polyMesh.C:1029
Foam::fvPatchFieldMapper.
virtual label size() const
Return size.
Definition: fvPatch.H:161
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
void operator=(const Field< vector > &)
tmp< vectorField > nf() const
Return face normals.
Definition: fvPatch.C:124
scalar deltaTValue() const
Return time step value.
Definition: TimeStateI.H:41
Calculate the mesh motion flux and convert fluxes from absolute to relative and back.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
tmp< surfaceScalarField > meshPhi(const volVectorField &U)
Definition: fvcMeshPhi.C:34
label n
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
const vectorField::subField faceCentres() const
Return face centres.
Definition: polyPatch.C:284
movingWallVelocityFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
A class for managing temporary objects.
Definition: PtrList.H:54
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
const GeometricField::Patch & patchField(const GeometricField &) const
Return the corresponding patchField of the named field.
Namespace for OpenFOAM.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:243
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:363