eddy.H
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) 2015 OpenFOAM Foundation
9  Copyright (C) 2016 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 Class
28  Foam::eddy
29 
30 Description
31  Class to describe eddies for the turbulentDFSEMInletFvPatchVectorField
32  boundary condition.
33 
34 SourceFiles
35  eddy.C
36  eddyI.H
37  eddyIO.C
38 
39 \*---------------------------------------------------------------------------*/
40 
41 #ifndef turbulentDFSEMInletFvPatchVectorField_eddy_H
42 #define turbulentDFSEMInletFvPatchVectorField_eddy_H
43 
44 #include "vector.H"
45 #include "point.H"
46 #include "tensor.H"
47 #include "Random.H"
48 #include "boundBox.H"
49 
50 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
51 
52 namespace Foam
53 {
54 
55 // Forward declarations
56 class eddy;
57 class Istream;
58 class Ostream;
59 
60 bool operator==(const eddy& a, const eddy& b);
61 bool operator!=(const eddy& a, const eddy& b);
62 Istream& operator>>(Istream& is, eddy& e);
63 Ostream& operator<<(Ostream& os, const eddy& e);
64 
65 
66 /*---------------------------------------------------------------------------*\
67  Class eddy Declaration
68 \*---------------------------------------------------------------------------*/
69 
70 class eddy
71 {
72  // Private data
73 
74  static label Gamma2Values[8];
75  static UList<label> Gamma2;
76 
77  //- Patch face index that spawned the eddy
78  label patchFaceI_;
79 
80  //- Reference position
81  point position0_;
82 
83  //- Distance from reference position in normal direction
84  scalar x_;
85 
86  //- Length scales in 3-D space
87  vector sigma_;
88 
89  //- Time-averaged intensity
90  vector alpha_;
91 
92  //- Coordinate system transformation from local to global axes
93  // X-direction aligned with max stress eigenvalue
94  tensor Rpg_;
95 
96  //- Model coefficient c1
97  scalar c1_;
98 
99  //- Index of streamwise direction
100  label dir1_;
101 
102 
103  // Private Member Functions
104 
105  //- Set the eddy scales: length, intensity
106  bool setScales
107  (
108  const scalar sigmaX,
109  const label gamma2,
110  const vector& e,
111  const vector& lambda,
112  vector& sigma,
113  vector& alpha
114  ) const;
115 
116  //- Return a number with zero mean and unit variance
117  inline scalar epsi(Random& rndGen) const;
118 
119 
120 public:
121 
122  // Constructors
123 
124  //- Construct null
125  eddy();
126 
127  //- Construct from Istream
128  eddy(Istream& is);
129 
130  //- Construct from components
131  eddy
132  (
133  const label patchFaceI, // patch face index
134  const point& position0, // reference position
135  const scalar x, // distance from reference position
136  const scalar sigmaX, // length scale
137  const symmTensor& R, // Stress tensor
138  Random& rndGen
139  );
140 
141  //- Construct copy
142  eddy(const eddy& e);
143 
144  static int debug;
145 
146 
147  // Public Member Functions
148 
149  // Access
150 
151  //- Return the patch face index that spawned the eddy
152  inline label patchFaceI() const;
153 
154  //- Return the reference position
155  inline const point& position0() const;
156 
157  //- Return the distance from the reference position
158  inline scalar x() const;
159 
160  //- Return the lLength scales in 3-D space
161  inline const vector& sigma() const;
162 
163  //- Return the time-averaged intensity
164  inline const vector& alpha() const;
165 
166  //- Return the coordinate system transformation from local
167  // principal to global axes
168  inline const tensor& Rpg() const;
169 
170  //- Return the model coefficient c1
171  inline scalar c1() const;
172 
173  //- Return the eddy position
174  inline point position(const vector& n) const;
175 
176  //- Return the index of the streamwise direction
177  inline label dir1() const;
178 
179  //- Return random vector of -1 and 1's
180  inline vector epsilon(Random& rndGen) const;
181 
182 
183  // Helper functions
184 
185  //- Volume
186  inline scalar volume() const;
187 
188  //- Move the eddy
189  inline void move(const scalar dx);
190 
191  //- Eddy bounds
192  inline boundBox bounds(const bool global = true) const;
193 
194 
195  // Evaluate
196 
197  //- Return the fluctuating velocity contribution at local point xp
198  vector uDash(const point& xp, const vector& n) const;
199 
200 
201  // Writing
202 
203  //- Write the eddy centre in OBJ format
204  void writeCentreOBJ(const vector& n, Ostream& os) const;
205 
206  //- Write the eddy surface in OBJ format
207  // Returns the number of points used to describe the eddy surface
208  label writeSurfaceOBJ
209  (
210  const label pointOffset,
211  const vector& n,
212  Ostream& os
213  ) const;
214 
215 
216  // Member Operators
217 
218  void operator=(const eddy& e);
219 
220 
221  // Friend Operators
222 
223  friend bool operator==(const eddy& a, const eddy& b)
224  {
225  return
226  a.patchFaceI_ == b.patchFaceI_
227  && a.position0_ == b.position0_
228  && a.x_ == b.x_
229  && a.sigma_ == b.sigma_
230  && a.alpha_ == b.alpha_
231  && a.Rpg_ == b.Rpg_
232  && a.c1_ == b.c1_
233  && a.dir1_ == b.dir1_;
234  }
235 
236  friend bool operator!=(const eddy& a, const eddy& b)
237  {
238  return !(a == b);
239  }
240 
241 
242  // IOstream Operators
243 
244  friend Istream& operator>>(Istream& is, eddy& e);
245  friend Ostream& operator<<(Ostream& os, const eddy& e);
246 };
247 
248 
249 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
250 
251 } // End namespace Foam
252 
253 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
254 
255 #include "eddyI.H"
256 
257 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
258 
259 #endif
260 
261 // ************************************************************************* //
Foam::Random
Random number generator.
Definition: Random.H:59
Foam::Tensor< scalar >
Foam::eddy::x
scalar x() const
Return the distance from the reference position.
Definition: eddyI.H:54
Foam::SymmTensor< scalar >
Foam::eddy::operator==
friend bool operator==(const eddy &a, const eddy &b)
Definition: eddy.H:222
Foam::eddy::debug
static int debug
Definition: eddy.H:143
Foam::eddy::operator<<
friend Ostream & operator<<(Ostream &os, const eddy &e)
Foam::eddy::writeCentreOBJ
void writeCentreOBJ(const vector &n, Ostream &os) const
Write the eddy centre in OBJ format.
Definition: eddy.C:254
Foam::eddy::epsilon
vector epsilon(Random &rndGen) const
Return random vector of -1 and 1's.
Definition: eddyI.H:90
point.H
Foam::eddy::bounds
boundBox bounds(const bool global=true) const
Eddy bounds.
Definition: eddyI.H:108
Foam::eddy::operator!=
friend bool operator!=(const eddy &a, const eddy &b)
Definition: eddy.H:235
Foam::eddy::volume
scalar volume() const
Volume.
Definition: eddyI.H:96
Foam::eddy::move
void move(const scalar dx)
Move the eddy.
Definition: eddyI.H:102
Foam::eddy
Class to describe eddies for the turbulentDFSEMInletFvPatchVectorField boundary condition.
Definition: eddy.H:69
Foam::operator>>
Istream & operator>>(Istream &, directionInfo &)
Definition: directionInfo.C:230
tensor.H
Foam::eddy::operator=
void operator=(const eddy &e)
Definition: eddyIO.C:51
Foam::eddy::writeSurfaceOBJ
label writeSurfaceOBJ(const label pointOffset, const vector &n, Ostream &os) const
Write the eddy surface in OBJ format.
Definition: eddy.C:265
Foam::eddy::Rpg
const tensor & Rpg() const
Return the coordinate system transformation from local.
Definition: eddyI.H:72
Foam::operator!=
bool operator!=(const eddy &a, const eddy &b)
Definition: eddy.H:235
Foam::operator<<
Ostream & operator<<(Ostream &, const boundaryPatch &p)
Write boundaryPatch as dictionary entries (without surrounding braces)
Definition: boundaryPatch.C:83
n
label n
Definition: TABSMDCalcMethod2.H:31
Foam::eddy::eddy
eddy()
Construct null.
Definition: eddy.C:115
R
#define R(A, B, C, D, E, F, K, M)
Foam::constant::physicoChemical::b
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
Foam::operator==
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
Foam::Istream
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:61
Foam::eddy::position
point position(const vector &n) const
Return the eddy position.
Definition: eddyI.H:78
lambda
dimensionedScalar lambda("lambda", dimTime/sqr(dimLength), laminarTransport)
Foam::eddy::dir1
label dir1() const
Return the index of the streamwise direction.
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::eddy::c1
scalar c1() const
Return the model coefficient c1.
Definition: eddyI.H:84
Random.H
boundBox.H
Foam::eddy::uDash
vector uDash(const point &xp, const vector &n) const
Return the fluctuating velocity contribution at local point xp.
Definition: eddy.C:229
eddyI.H
Foam::Vector< scalar >
Foam::UList< label >
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
vector.H
Foam::boundBox
A bounding box defined in terms of min/max extrema points.
Definition: boundBox.H:63
rndGen
Random rndGen
Definition: createFields.H:23
Foam::eddy::operator>>
friend Istream & operator>>(Istream &is, eddy &e)
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::eddy::alpha
const vector & alpha() const
Return the time-averaged intensity.
Definition: eddyI.H:66
Foam::eddy::sigma
const vector & sigma() const
Return the lLength scales in 3-D space.
Definition: eddyI.H:60
Foam::eddy::patchFaceI
label patchFaceI() const
Return the patch face index that spawned the eddy.
Definition: eddyI.H:42
Foam::eddy::position0
const point & position0() const
Return the reference position.
Definition: eddyI.H:48