movingWallSlipVelocityFvPatchVectorField.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) 2021-2024 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  directionMixedFvPatchVectorField(p, iF)
42 {
43  const vectorField n(p.nf());
44 
45  refValue() = Zero;
46  refGrad() = Zero;
47  valueFraction() = sqr(n);
48 }
49 
50 
53 (
54  const fvPatch& p,
56  const dictionary& dict
57 )
58 :
59  directionMixedFvPatchVectorField(p, iF)
60 {
62  (
63  vectorField("value", iF.dimensions(), dict, p.size())
64  );
65 
66  const vectorField n(p.nf());
67 
68  refValue() = sqr(n) & *this;
69  refGrad() = Zero;
70  valueFraction() = sqr(n);
71 }
72 
73 
76 (
78  const fvPatch& p,
80  const fieldMapper& mapper
81 )
82 :
83  directionMixedFvPatchVectorField(ptf, p, iF, mapper)
84 {}
85 
86 
89 (
92 )
93 :
94  directionMixedFvPatchVectorField(mwvpvf, iF)
95 {}
96 
97 
98 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
99 
101 {
102  if (updated())
103  {
104  return;
105  }
106 
107  const fvMesh& mesh = internalField().mesh();
108 
109  if (mesh.moving())
110  {
111  const fvPatch& p = patch();
112 
113  const volVectorField& U =
114  static_cast<const volVectorField&>(internalField());
115 
116  const scalarField phip(fvc::meshPhi(U, p.index()));
117 
118  const vectorField n(p.nf());
119  const scalarField& magSf = p.magSf();
120 
121  tmp<scalarField> Un = phip/(magSf + vSmall);
122 
123  refValue() = n*Un;
124  refGrad() = Zero;
125  valueFraction() = sqr(n);
126  }
127 
128  directionMixedFvPatchVectorField::updateCoeffs();
130 }
131 
132 
134 {
136  writeEntry(os, "value", *this);
137 }
138 
139 
140 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
141 
142 namespace Foam
143 {
145  (
148  );
149 }
150 
151 // ************************************************************************* //
label n
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...
const dimensionSet & dimensions() const
Return dimensions.
Generic GeometricField class.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
Abstract base class for field mapping.
Definition: fieldMapper.H:48
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:99
const polyMesh & mesh() const
Return reference to polyMesh.
Definition: fvMesh.H:441
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:88
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:229
friend Ostream & operator(Ostream &, const fvPatchField< Type > &)
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
This boundary condition provides a slip velocity condition for cases with moving walls.
movingWallSlipVelocityFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
bool moving() const
Is mesh moving.
Definition: polyMesh.H:476
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.
static const zero Zero
Definition: zero.H:97
dimensionedSymmTensor sqr(const dimensionedVector &dv)
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
Field< vector > vectorField
Specialisation of Field<T> for vector.
void evaluate(GeometricField< Type, PatchField, GeoMesh > &result, const Function1< Type > &func, const GeometricField< Type, PatchField, GeoMesh > &x)
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
dictionary dict
volScalarField & p
Foam::surfaceFields.