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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2017 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
30 #include "volFields.H"
31 #include "surfaceFields.H"
32 #include "SRFModel.H"
33 #include "steadyStateDdtScheme.H"
34 
35 
36 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
37 
40 (
41  const fvPatch& p,
43 )
44 :
45  inletOutletFvPatchVectorField(p, iF),
46  relative_(false),
47  UInf_(Zero)
48 {}
49 
50 
53 (
55  const fvPatch& p,
57  const fvPatchFieldMapper& mapper
58 )
59 :
60  inletOutletFvPatchVectorField(ptf, p, iF, mapper),
61  relative_(ptf.relative_),
62  UInf_(ptf.UInf_)
63 {}
64 
65 
68 (
69  const fvPatch& p,
71  const dictionary& dict
72 )
73 :
74  inletOutletFvPatchVectorField(p, iF),
75  relative_(dict.lookupOrDefault("relative", false)),
76  UInf_(dict.lookup("UInf"))
77 {
78  this->phiName_ = dict.lookupOrDefault<word>("phi","phi");
79 
80  fvPatchVectorField::operator=(vectorField("value", dict, p.size()));
81 }
82 
83 
86 (
88 )
89 :
90  inletOutletFvPatchVectorField(srfvpvf),
91  relative_(srfvpvf.relative_),
92  UInf_(srfvpvf.UInf_)
93 {}
94 
95 
98 (
101 )
102 :
103  inletOutletFvPatchVectorField(srfvpvf, iF),
104  relative_(srfvpvf.relative_),
105  UInf_(srfvpvf.UInf_)
106 {}
107 
108 
109 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
110 
112 {
113  if (updated())
114  {
115  return;
116  }
117 
118  // Get reference to the SRF model
119  const SRF::SRFModel& srf =
120  db().lookupObject<SRF::SRFModel>("SRFProperties");
121 
122 
123  word ddtScheme
124  (
125  this->internalField().mesh()
126  .ddtScheme(this->internalField().name())
127  );
128 
130  {
131  // If not relative to the SRF include the effect of the SRF
132  if (!relative_)
133  {
134  refValue() = UInf_ - srf.velocity(patch().Cf());
135  }
136  // If already relative to the SRF simply supply the inlet value
137  // as a fixed value
138  else
139  {
140  refValue() = UInf_;
141  }
142  }
143  else
144  {
145  scalar time = this->db().time().value();
146  scalar theta = time*mag(srf.omega().value());
147 
148  refValue() =
149  cos(theta)*UInf_ + sin(theta)*(srf.axis() ^ UInf_)
150  - srf.velocity(patch().Cf());
151  }
152 
153  // Set the inlet-outlet choice based on the direction of the freestream
154  valueFraction() = 1.0 - pos0(refValue() & patch().Sf());
155 
157 }
158 
159 
161 {
163  os.writeEntry("relative", relative_);
164  os.writeEntry("UInf", UInf_);
165  os.writeEntry("phi", this->phiName_);
166  writeEntry("value", os);
167 }
168 
169 
170 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
171 
172 namespace Foam
173 {
175  (
178  );
179 }
180 
181 
182 // ************************************************************************* //
Foam::fvPatchField< vector >
volFields.H
Foam::fvPatchField< vector >::write
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:364
p
volScalarField & p
Definition: createFieldRefs.H:8
steadyStateDdtScheme.H
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::Zero
static constexpr const zero Zero
Global zero.
Definition: zero.H:128
Foam::sin
dimensionedScalar sin(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:264
surfaceFields.H
Foam::surfaceFields.
Foam::dimensioned::value
const Type & value() const
Return const reference to value.
Definition: dimensionedType.C:404
Foam::pos0
dimensionedScalar pos0(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:188
Foam::SRF::SRFModel::velocity
vectorField velocity(const vectorField &positions) const
Return velocity vector from positions.
Definition: SRFModel.C:168
Foam::SRF::SRFModel
Top level model for single rotating frame.
Definition: SRFModel.H:64
Foam::vectorField
Field< vector > vectorField
Specialisation of Field<T> for vector.
Definition: primitiveFieldsFwd.H:54
Foam::SRFFreestreamVelocityFvPatchVectorField
Freestream velocity condition to be used in conjunction with the single rotating frame (SRF) model (s...
Definition: SRFFreestreamVelocityFvPatchVectorField.H:128
Foam::SRF::SRFModel::axis
const vector & axis() const
Return the axis of rotation.
Definition: SRFModel.C:106
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:63
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
Foam::SRFFreestreamVelocityFvPatchVectorField::write
virtual void write(Ostream &) const
Write.
Definition: SRFFreestreamVelocityFvPatchVectorField.C:160
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::SRF::SRFModel::omega
const dimensionedVector & omega() const
Return the angular velocity field [rad/s].
Definition: SRFModel.C:112
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::fvPatchField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:313
SRFModel.H
Foam::fv::steadyStateDdtScheme
SteadyState implicit/explicit ddt which returns 0.
Definition: steadyStateDdtScheme.H:59
Foam::SRFFreestreamVelocityFvPatchVectorField::SRFFreestreamVelocityFvPatchVectorField
SRFFreestreamVelocityFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
Definition: SRFFreestreamVelocityFvPatchVectorField.C:40
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::SRFFreestreamVelocityFvPatchVectorField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: SRFFreestreamVelocityFvPatchVectorField.C:111
Foam::Ostream::writeEntry
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:219
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:47
SRFFreestreamVelocityFvPatchVectorField.H
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::makePatchTypeField
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54
Foam::cos
dimensionedScalar cos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:265