rotatingPressureInletOutletVelocityFvPatchVectorField.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-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"
30 
31 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
32 
33 void Foam::rotatingPressureInletOutletVelocityFvPatchVectorField::
34 calcTangentialVelocity()
35 {
36  const scalar omega = omega_.value(db().time().value());
37 
39  (
40  (-omega)*((patch().Cf() - origin_) ^ (axis_/mag(axis_)))
41  );
42 
43  const vectorField n(patch().nf());
44 
45  refValue() = tangentialVelocity - n*(n & tangentialVelocity);
46 }
47 
48 
49 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
50 
53 (
54  const fvPatch& p,
56  const dictionary& dict
57 )
58 :
60  origin_(dict.lookup<vector>("origin", dimLength)),
61  axis_(dict.lookup<vector>("axis", dimless)),
62  omega_(db().time(), dict)
63 {
64  calcTangentialVelocity();
65 }
66 
67 
70 (
72  const fvPatch& p,
74  const fieldMapper& mapper
75 )
76 :
78  origin_(pvf.origin_),
79  axis_(pvf.axis_),
80  omega_(pvf.omega_)
81 {
82  calcTangentialVelocity();
83 }
84 
85 
88 (
91 )
92 :
94  origin_(pvf.origin_),
95  axis_(pvf.axis_),
96  omega_(pvf.omega_)
97 {
98  calcTangentialVelocity();
99 }
100 
101 
102 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
103 
105 (
106  Ostream& os
107 ) const
108 {
110  writeEntry(os, "phi", phiName());
111  writeEntry(os, "origin", origin_);
112  writeEntry(os, "axis", axis_);
113  writeEntry(os, omega_);
114  writeEntry(os, "value", *this);
115 }
116 
117 
118 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
119 
120 namespace Foam
121 {
123  (
126  );
127 }
128 
129 // ************************************************************************* //
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...
scalar value(const scalar t) const
Return value for time t.
Definition: omega1I.H:30
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
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
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
Velocity inlet/outlet boundary condition for patches where the pressure is specified in some manner,...
const autoPtr< Function1< vector > > & tangentialVelocity() const
Return the tangential velocity Function1.
This velocity inlet/outlet boundary condition is applied to patches in a rotating frame where the pre...
rotatingPressureInletOutletVelocityFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &, const dictionary &)
Construct from patch, internal field and dictionary.
const scalar omega
Namespace for OpenFOAM.
const dimensionSet dimless
const dimensionSet dimLength
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
dimensioned< scalar > mag(const dimensioned< Type > &)
Field< vector > vectorField
Specialisation of Field<T> for vector.
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
dictionary dict
volScalarField & p
Foam::surfaceFields.