SRFFreestreamVelocityFvPatchVectorField.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------* \
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2011-2016 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 #include "SRFModel.H"
31 #include "steadyStateDdtScheme.H"
32 
33 
34 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
35 
38 (
39  const fvPatch& p,
41 )
42 :
43  inletOutletFvPatchVectorField(p, iF),
44  relative_(false),
45  UInf_(Zero)
46 {}
47 
48 
51 (
53  const fvPatch& p,
55  const fvPatchFieldMapper& mapper
56 )
57 :
58  inletOutletFvPatchVectorField(ptf, p, iF, mapper),
59  relative_(ptf.relative_),
60  UInf_(ptf.UInf_)
61 {}
62 
63 
66 (
67  const fvPatch& p,
69  const dictionary& dict
70 )
71 :
72  inletOutletFvPatchVectorField(p, iF),
73  relative_(dict.lookupOrDefault("relative", false)),
74  UInf_(dict.lookup("UInf"))
75 {
76  this->phiName_ = dict.lookupOrDefault<word>("phi","phi");
77 
78  fvPatchVectorField::operator=(vectorField("value", dict, p.size()));
79 }
80 
81 
84 (
86 )
87 :
88  inletOutletFvPatchVectorField(srfvpvf),
89  relative_(srfvpvf.relative_),
90  UInf_(srfvpvf.UInf_)
91 {}
92 
93 
96 (
99 )
100 :
101  inletOutletFvPatchVectorField(srfvpvf, iF),
102  relative_(srfvpvf.relative_),
103  UInf_(srfvpvf.UInf_)
104 {}
105 
106 
107 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
108 
110 {
111  if (updated())
112  {
113  return;
114  }
115 
116  // Get reference to the SRF model
117  const SRF::SRFModel& srf =
118  db().lookupObject<SRF::SRFModel>("SRFProperties");
119 
120 
121  word ddtScheme
122  (
123  this->internalField().mesh()
124  .ddtScheme(this->internalField().name())
125  );
126 
128  {
129  // If not relative to the SRF include the effect of the SRF
130  if (!relative_)
131  {
132  refValue() = UInf_ - srf.velocity(patch().Cf());
133  }
134  // If already relative to the SRF simply supply the inlet value
135  // as a fixed value
136  else
137  {
138  refValue() = UInf_;
139  }
140  }
141  else
142  {
143  scalar time = this->db().time().value();
144  scalar theta = time*mag(srf.omega().value());
145 
146  refValue() =
147  cos(theta)*UInf_ + sin(theta)*(srf.axis() ^ UInf_)
148  - srf.velocity(patch().Cf());
149  }
150 
151  // Set the inlet-outlet choice based on the direction of the freestream
152  valueFraction() = 1.0 - pos(refValue() & patch().Sf());
153 
155 }
156 
157 
159 {
161  os.writeKeyword("relative") << relative_ << token::END_STATEMENT << nl;
162  os.writeKeyword("UInf") << UInf_ << token::END_STATEMENT << nl;
163  os.writeKeyword("phi") << this->phiName_ << token::END_STATEMENT << nl;
164  writeEntry("value", os);
165 }
166 
167 
168 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
169 
170 namespace Foam
171 {
173  (
176  );
177 }
178 
179 
180 // ************************************************************************* //
Foam::surfaceFields.
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:65
volVectorField vectorField(fieldObject, mesh)
const Type & value() const
Return const reference to value.
const vector & axis() const
Return the axis of rotation.
Definition: SRFModel.C:109
Macros for easy insertion into run-time selection tables.
makePatchTypeField(fvPatchVectorField, SRFFreestreamVelocityFvPatchVectorField)
vectorField velocity(const vectorField &positions) const
Return velocity vector from positions.
Definition: SRFModel.C:171
dimensionedScalar pos(const dimensionedScalar &ds)
dynamicFvMesh & mesh
dimensionedScalar cos(const dimensionedScalar &ds)
A class for handling words, derived from string.
Definition: word.H:59
Foam::fvPatchFieldMapper.
SRFFreestreamVelocityFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
virtual label size() const
Return size.
Definition: fvPatch.H:161
static const zero Zero
Definition: zero.H:91
SteadyState implicit/explicit ddt which returns 0.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
dimensionedScalar sin(const dimensionedScalar &ds)
static const char nl
Definition: Ostream.H:262
Ostream & writeKeyword(const keyType &)
Write the keyword followed by an appropriate indentation.
Definition: Ostream.C:54
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:312
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
dimensioned< scalar > mag(const dimensioned< Type > &)
const dimensionedVector & omega() const
Return the angular velocity field [rad/s].
Definition: SRFModel.C:115
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
Top level model for single rotating frame.
Definition: SRFModel.H:62
Namespace for OpenFOAM.
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:363
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:451