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-2018 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().timeOutputValue();
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 (
68  const fvPatch& p,
70  const fvPatchFieldMapper& mapper
71 )
72 :
74  omega_(ptf.omega_, false)
75 {
76  calcTangentialVelocity();
77 }
78 
79 
82 (
83  const fvPatch& p,
85  const dictionary& dict
86 )
87 :
89  omega_(Function1<vector>::New("omega", dict))
90 {
91  calcTangentialVelocity();
92 }
93 
94 
97 (
99 )
100 :
102  omega_(rppvf.omega_, false)
103 {
104  calcTangentialVelocity();
105 }
106 
107 
110 (
113 )
114 :
116  omega_(rppvf.omega_, false)
117 {
118  calcTangentialVelocity();
119 }
120 
121 
122 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
123 
125 (
126  Ostream& os
127 ) const
128 {
130  os.writeKeyword("phi") << phiName() << token::END_STATEMENT << nl;
131  omega_->writeData(os);
132  writeEntry("value", os);
133 }
134 
135 
136 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
137 
138 namespace Foam
139 {
141  (
144  );
145 }
146 
147 // ************************************************************************* //
Foam::surfaceFields.
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
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: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: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:362
Macros for easy insertion into run-time selection tables.
pressureInletOutletVelocityFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
const vectorField & tangentialVelocity() const
Return the tangential velocity.
Foam::fvPatchFieldMapper.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
static const char nl
Definition: Ostream.H:265
Ostream & writeKeyword(const keyType &)
Write the keyword followed by an appropriate indentation.
Definition: Ostream.C:54
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.