rotatingWallVelocityFvPatchVectorField.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 
28 #include "volFields.H"
29 #include "surfaceFields.H"
30 
31 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
32 
35 (
36  const fvPatch& p,
38  const dictionary& dict
39 )
40 :
41  fixedValueFvPatchField<vector>(p, iF, dict, false),
42  origin_(dict.lookup<vector>("origin")),
43  axis_(dict.lookup<vector>("axis")),
44  omega_(db().time(), dict)
45 {
46  if (dict.found("value"))
47  {
49  (
50  vectorField("value", iF.dimensions(), dict, p.size())
51  );
52  }
53  else
54  {
55  // Evaluate the wall velocity
56  updateCoeffs();
57  }
58 }
59 
60 
63 (
65  const fvPatch& p,
67  const fieldMapper& mapper
68 )
69 :
70  fixedValueFvPatchField<vector>(pvf, p, iF, mapper),
71  origin_(pvf.origin_),
72  axis_(pvf.axis_),
73  omega_(pvf.omega_)
74 {}
75 
76 
79 (
82 )
83 :
85  origin_(pvf.origin_),
86  axis_(pvf.axis_),
87  omega_(pvf.omega_)
88 {}
89 
90 
91 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
92 
94 {
95  if (updated())
96  {
97  return;
98  }
99 
100  const scalar omega = omega_.value(db().time().value());
101 
102  // Calculate the rotating wall velocity from the specification of the motion
103  const vectorField Up
104  (
105  (-omega)*((patch().Cf() - origin_) ^ (axis_/mag(axis_)))
106  );
107 
108  // Remove the component of Up normal to the wall
109  // just in case it is not exactly circular
110  const vectorField n(patch().nf());
111  vectorField::operator=(Up - n*(n & Up));
112 
113  fixedValueFvPatchVectorField::updateCoeffs();
114 }
115 
116 
118 {
120  writeEntry(os, "origin", origin_);
121  writeEntry(os, "axis", axis_);
122  writeEntry(os, omega_);
123  writeEntry(os, "value", *this);
124 }
125 
126 
127 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
128 
129 namespace Foam
130 {
132  (
135  );
136 }
137 
138 // ************************************************************************* //
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.
void operator=(const Field< vector > &)
Definition: Field.C:550
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
const Type & value() const
Return const reference to value.
Abstract base class for field mapping.
Definition: fieldMapper.H:48
This boundary condition supplies a fixed value constraint, and is the base class for a number of othe...
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
Convenience class to handle the input of constant rotational speed. Reads an omega entry with default...
Definition: omega.H:54
Condition on velocity for a boundary consisting of a rotating solid of revolution,...
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
rotatingWallVelocityFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &, const dictionary &)
Construct from patch, internal field and dictionary.
Namespace for OpenFOAM.
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.