movingWallVelocityFvPatchVectorField.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2023 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  const dictionary& dict
40 )
41 :
42  fixedValueFvPatchVectorField(p, iF, dict)
43 {}
44 
45 
48 (
50  const fvPatch& p,
52  const fvPatchFieldMapper& mapper
53 )
54 :
55  fixedValueFvPatchVectorField(ptf, p, iF, mapper)
56 {}
57 
58 
61 (
64 )
65 :
66  fixedValueFvPatchVectorField(mwvpvf, iF)
67 {}
68 
69 
70 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
71 
73 {
74  if (updated())
75  {
76  return;
77  }
78 
79  const fvMesh& mesh = patch().boundaryMesh().mesh();
80 
81  if (mesh.moving())
82  {
83  const fvPatch& p = patch();
84 
85  const volVectorField& U =
86  static_cast<const volVectorField&>(internalField());
87 
88  const vectorField n(p.nf());
89  tmp<scalarField> Un = fvc::meshPhi(U, p.index())/(p.magSf() + vSmall);
90 
91  const pointField& oldPoints = mesh.oldPoints();
92 
93  const polyPatch& pp = p.patch();
94 
95  vectorField oldFc(pp.size());
96 
97  forAll(oldFc, i)
98  {
99  oldFc[i] = pp[i].centre(oldPoints);
100  }
101 
102  const vectorField Up
103  (
104  (pp.faceCentres() - oldFc)/mesh.time().deltaTValue()
105  );
106 
107  vectorField::operator=(Up + n*(Un - (n & Up)));
108  }
109 
110  fixedValueFvPatchVectorField::updateCoeffs();
111 }
112 
113 
115 {
117  writeEntry(os, "value", *this);
118 }
119 
120 
121 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
122 
123 namespace Foam
124 {
126  (
129  );
130 }
131 
132 // ************************************************************************* //
label n
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
void operator=(const Field< vector > &)
Definition: Field.C:526
Generic GeometricField class.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
scalar deltaTValue() const
Return time step value.
Definition: TimeStateI.H:34
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:101
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:406
Foam::fvPatchFieldMapper.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:87
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:231
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
This boundary condition provides a no-slip velocity condition for ridged moving walls or flexible mov...
movingWallVelocityFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &, const dictionary &)
Construct from patch, internal field and dictionary.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
const polyMesh & mesh() const
Return the mesh reference.
virtual const pointField & oldPoints() const
Return old points for mesh motion.
Definition: polyMesh.C:1399
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:403
bool moving() const
Is mesh moving.
Definition: polyMesh.H:475
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:70
const vectorField::subField faceCentres() const
Return face centres.
Definition: polyPatch.C:276
A class for managing temporary objects.
Definition: tmp.H:55
Calculate the mesh motion flux and convert fluxes from absolute to relative and back.
U
Definition: pEqn.H:72
tmp< surfaceScalarField > meshPhi(const volVectorField &U)
Definition: fvcMeshPhi.C:34
Namespace for OpenFOAM.
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
dictionary dict
volScalarField & p
Foam::surfaceFields.