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-2021 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 
28 #include "volFields.H"
29 #include "surfaceFields.H"
30 
31 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
32 
33 void Foam::rotatingPressureInletOutletVelocityFvPatchVectorField::
34 calcTangentialVelocity()
35 {
36  const scalar t = this->db().time().userTimeValue();
37  vector om = omega_->value(t);
38 
39  vector axisHat = om/mag(om);
41  (
42  (-om) ^ (patch().Cf() - axisHat*(axisHat & patch().Cf()))
43  );
44 
45  const vectorField n(patch().nf());
46  refValue() = tangentialVelocity - n*(n & tangentialVelocity);
47 }
48 
49 
50 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
51 
54 (
55  const fvPatch& p,
57 )
58 :
60  omega_()
61 {}
62 
63 
66 (
67  const fvPatch& p,
69  const dictionary& dict
70 )
71 :
73  omega_(Function1<vector>::New("omega", dict))
74 {
75  calcTangentialVelocity();
76 }
77 
78 
81 (
83  const fvPatch& p,
85  const fvPatchFieldMapper& mapper
86 )
87 :
89  omega_(ptf.omega_, false)
90 {
91  calcTangentialVelocity();
92 }
93 
94 
97 (
100 )
101 :
103  omega_(rppvf.omega_, false)
104 {
105  calcTangentialVelocity();
106 }
107 
108 
109 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
110 
112 (
113  Ostream& os
114 ) const
115 {
117  writeEntry(os, "phi", phiName());
118  writeEntry(os, omega_());
119  writeEntry(os, "value", *this);
120 }
121 
122 
123 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
124 
125 namespace Foam
126 {
128  (
131  );
132 }
133 
134 // ************************************************************************* //
Foam::surfaceFields.
Run-time selectable general function of one variable.
Definition: Function1.H:52
This velocity inlet/outlet boundary condition is applied to patches in a rotating frame where the pre...
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
const autoPtr< Function1< vector > > & tangentialVelocity() const
Return the tangential velocity Function1.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:63
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
rotatingPressureInletOutletVelocityFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:243
Macros for easy insertion into run-time selection tables.
pressureInletOutletVelocityFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
Foam::fvPatchFieldMapper.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
dimensioned< scalar > mag(const dimensioned< Type > &)
label n
Field< vector > vectorField
Specialisation of Field<T> for vector.
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Namespace for OpenFOAM.