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;
57class Istream;
58class Ostream;
59
60bool operator==(const eddy& a, const eddy& b);
61bool operator!=(const eddy& a, const eddy& b);
62Istream& operator>>(Istream& is, eddy& e);
63Ostream& operator<<(Ostream& os, const eddy& e);
64
65
66/*---------------------------------------------------------------------------*\
67 Class eddy Declaration
68\*---------------------------------------------------------------------------*/
70class 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 //- Integral-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 (0,1,2)
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,
114 ) const;
115
116 //- Return a number with zero mean and unit variance
117 inline scalar epsi(Random& rndGen) const;
118
119
120public:
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, // integral-length scale
137 const symmTensor& R, // Reynolds stress tensor
139 );
140
141 //- Copy construct
142 eddy(const eddy& e);
143
144
145 // Public Data
146
147 //- Flag to activate debug statements
148 static int debug;
149
150
151 // Public Member Functions
152
153 // Access
154
155 //- Return the patch face index that spawned the eddy
156 inline label patchFaceI() const noexcept;
157
158 //- Return the reference position
159 inline const point& position0() const noexcept;
160
161 //- Return the distance from the reference position
162 inline scalar x() const noexcept;
163
164 //- Return the length scales in 3-D space
165 inline const vector& sigma() const noexcept;
166
167 //- Return the time-averaged intensity
168 inline const vector& alpha() const noexcept;
169
170 //- Return the coordinate system transformation from local
171 //- principal to global axes
172 inline const tensor& Rpg() const noexcept;
173
174 //- Return the model coefficient c1
175 inline scalar c1() const noexcept;
176
177 //- Return the eddy position
178 inline point position(const vector& n) const;
179
180 //- Return the index of the streamwise direction (0,1,2)
181 label dir1() const noexcept { return dir1_; }
182
183 //- Return random vector of -1 and 1's
184 inline vector epsilon(Random& rndGen) const;
185
186
187 // Helper functions
188
189 //- Volume
190 inline scalar volume() const;
191
192 //- Move the eddy
193 inline void move(const scalar dx);
194
195 //- Eddy bounds
196 inline boundBox bounds(const bool global = true) const;
197
198
199 // Evaluation
200
201 //- Return the fluctuating velocity contribution at local point xp
202 vector uPrime(const point& xp, const vector& n) const;
203
204
205 // Writing
206
207 //- Write the eddy centre in OBJ format
208 void writeCentreOBJ(const vector& n, Ostream& os) const;
209
210 //- Write the eddy surface in OBJ format
211 // Returns the number of points used to describe the eddy surface
212 label writeSurfaceOBJ
213 (
214 const label pointOffset,
215 const vector& n,
216 Ostream& os
217 ) const;
218
219
220 // Member Operators
221
222 void operator=(const eddy& e);
223
224
225 // Friend Operators
227 friend bool operator==(const eddy& a, const eddy& b)
228 {
229 return
230 a.patchFaceI_ == b.patchFaceI_
231 && a.position0_ == b.position0_
232 && a.x_ == b.x_
233 && a.sigma_ == b.sigma_
234 && a.alpha_ == b.alpha_
235 && a.Rpg_ == b.Rpg_
236 && a.c1_ == b.c1_
237 && a.dir1_ == b.dir1_;
238 }
240 friend bool operator!=(const eddy& a, const eddy& b)
241 {
242 return !(a == b);
243 }
244
245
246 // IOstream Operators
249 friend Ostream& operator<<(Ostream& os, const eddy& e);
250};
251
252
253// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
254
255} // End namespace Foam
256
257// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
258
259#include "eddyI.H"
260
261// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
262
263#endif
264
265// ************************************************************************* //
#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:64
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
Random number generator.
Definition: Random.H:60
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:94
A bounding box defined in terms of min/max extrema points.
Definition: boundBox.H:64
Class to describe eddies for the turbulentDFSEMInletFvPatchVectorField boundary condition.
Definition: eddy.H:70
friend bool operator==(const eddy &a, const eddy &b)
Definition: eddy.H:226
friend bool operator!=(const eddy &a, const eddy &b)
Definition: eddy.H:239
void move(const scalar dx)
Move the eddy.
Definition: eddyI.H:102
scalar c1() const noexcept
Return the model coefficient c1.
Definition: eddyI.H:78
const tensor & Rpg() const noexcept
Definition: eddyI.H:72
eddy()
Construct null.
Definition: eddy.C:117
const point & position0() const noexcept
Return the reference position.
Definition: eddyI.H:48
scalar volume() const
Volume.
Definition: eddyI.H:96
void writeCentreOBJ(const vector &n, Ostream &os) const
Write the eddy centre in OBJ format.
Definition: eddy.C:256
label dir1() const noexcept
Return the index of the streamwise direction (0,1,2)
Definition: eddy.H:180
void operator=(const eddy &e)
Definition: eddyIO.C:51
const vector & alpha() const noexcept
Return the time-averaged intensity.
Definition: eddyI.H:66
vector uPrime(const point &xp, const vector &n) const
Return the fluctuating velocity contribution at local point xp.
Definition: eddy.C:231
friend Istream & operator>>(Istream &is, eddy &e)
point position(const vector &n) const
Return the eddy position.
Definition: eddyI.H:84
label writeSurfaceOBJ(const label pointOffset, const vector &n, Ostream &os) const
Write the eddy surface in OBJ format.
Definition: eddy.C:267
static int debug
Flag to activate debug statements.
Definition: eddy.H:147
label patchFaceI() const noexcept
Return the patch face index that spawned the eddy.
Definition: eddyI.H:42
friend Ostream & operator<<(Ostream &os, const eddy &e)
scalar x() const noexcept
Return the distance from the reference position.
Definition: eddyI.H:54
const vector & sigma() const noexcept
Return the length scales in 3-D space.
Definition: eddyI.H:60
scalar epsilon
OBJstream os(runTime.globalPath()/outputName)
Namespace for OpenFOAM.
bool operator!=(const eddy &a, const eddy &b)
Definition: eddy.H:239
Ostream & operator<<(Ostream &, const boundaryPatch &p)
Write boundaryPatch as dictionary entries (without surrounding braces)
Definition: boundaryPatch.C:83
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
Istream & operator>>(Istream &, directionInfo &)
const direction noexcept
Definition: Scalar.H:223
volScalarField & b
Definition: createFields.H:27
volScalarField & e
Definition: createFields.H:11
dimensionedScalar lambda("lambda", dimTime/sqr(dimLength), laminarTransport)
Random rndGen
Definition: createFields.H:23