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  Copyright (C) 2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
31 #include "volFields.H"
32 #include "surfaceFields.H"
33 #include "SRFModel.H"
35 
36 
37 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
38 
41 (
42  const fvPatch& p,
44 )
45 :
46  inletOutletFvPatchVectorField(p, iF),
47  relative_(false),
48  UInf_(Zero)
49 {}
50 
51 
54 (
56  const fvPatch& p,
58  const fvPatchFieldMapper& mapper
59 )
60 :
61  inletOutletFvPatchVectorField(ptf, p, iF, mapper),
62  relative_(ptf.relative_),
63  UInf_(ptf.UInf_)
64 {}
65 
66 
69 (
70  const fvPatch& p,
72  const dictionary& dict
73 )
74 :
75  inletOutletFvPatchVectorField(p, iF),
76  relative_(dict.getOrDefault("relative", false)),
77  UInf_(dict.get<vector>("UInf"))
78 {
79  this->phiName_ = dict.getOrDefault<word>("phi", "phi");
80 
82 }
83 
84 
87 (
89 )
90 :
91  inletOutletFvPatchVectorField(srfvpvf),
92  relative_(srfvpvf.relative_),
93  UInf_(srfvpvf.UInf_)
94 {}
95 
96 
99 (
102 )
103 :
104  inletOutletFvPatchVectorField(srfvpvf, iF),
105  relative_(srfvpvf.relative_),
106  UInf_(srfvpvf.UInf_)
107 {}
108 
109 
110 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
111 
113 {
114  if (updated())
115  {
116  return;
117  }
118 
119  // Get reference to the SRF model
120  const SRF::SRFModel& srf =
121  db().lookupObject<SRF::SRFModel>("SRFProperties");
122 
123 
124  word ddtScheme
125  (
126  this->internalField().mesh()
127  .ddtScheme(this->internalField().name())
128  );
129 
130  if (ddtScheme == fv::steadyStateDdtScheme<scalar>::typeName)
131  {
132  // If not relative to the SRF include the effect of the SRF
133  if (!relative_)
134  {
135  refValue() = UInf_ - srf.velocity(patch().Cf());
136  }
137  // If already relative to the SRF simply supply the inlet value
138  // as a fixed value
139  else
140  {
141  refValue() = UInf_;
142  }
143  }
144  else
145  {
146  scalar time = this->db().time().value();
147  scalar theta = time*mag(srf.omega().value());
148 
149  refValue() =
150  cos(theta)*UInf_ + sin(theta)*(srf.axis() ^ UInf_)
151  - srf.velocity(patch().Cf());
152  }
153 
154  // Set the inlet-outlet choice based on the direction of the freestream
155  valueFraction() = 1.0 - pos0(refValue() & patch().Sf());
156 
158 }
159 
160 
162 {
164  os.writeEntry("relative", relative_);
165  os.writeEntry("UInf", UInf_);
166  os.writeEntry("phi", this->phiName_);
167  writeEntry("value", os);
168 }
169 
170 
171 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
172 
173 namespace Foam
174 {
176  (
178  SRFFreestreamVelocityFvPatchVectorField
179  );
180 }
181 
182 
183 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
Foam::surfaceFields.
dictionary dict
fvPatchField< vector > fvPatchVectorField
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:120
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:68
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:312
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:333
Macros for easy insertion into run-time selection tables.
dynamicFvMesh & mesh
dimensionedScalar cos(const dimensionedScalar &ds)
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:52
A class for handling words, derived from Foam::string.
Definition: word.H:63
A FieldMapper for finite-volume patch fields.
SRFFreestreamVelocityFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:55
dimensionedScalar pos0(const dimensionedScalar &ds)
dimensionedScalar sin(const dimensionedScalar &ds)
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
OBJstream os(runTime.globalPath()/outputName)
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:270
Freestream velocity condition to be used in conjunction with the single rotating frame (SRF) model (s...
virtual void operator=(const UList< Type > &)
Definition: fvPatchField.C:352
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: areaFieldsFwd.H:42
const std::string patch
OpenFOAM patch number as a std::string.
Field< vector > vectorField
Specialisation of Field<T> for vector.
volScalarField & p
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:157