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