All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
turbGen.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2018 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "fft.H"
27 #include "turbGen.H"
28 #include "Kmesh.H"
29 #include "primitiveFields.H"
30 #include "Ek.H"
31 #include "mathematicalConstants.H"
32 
33 
34 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
35 
36 Foam::turbGen::turbGen(const Kmesh& k, const scalar EA, const scalar K0)
37 :
38  K(k),
39  Ea(EA),
40  k0(K0),
41  RanGen(label(0))
42 {}
43 
44 
45 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
46 
48 {
49  vectorField s(K.size());
50  scalarField rndPhases(K.size());
51 
52  forAll(K, i)
53  {
54  s[i] = RanGen.sample01<vector>();
55  rndPhases[i] = RanGen.scalar01();
56  }
57 
58  s = K ^ s;
59  s = s/(mag(s) + 1.0e-20);
60 
61  s = Ek(Ea, k0, mag(K))*s;
62 
64  (
66  (
68  sin(constant::mathematical::twoPi*rndPhases)*s),
69  K.nn()
70  )
71  );
72 
73  return ReImSum(up);
74 }
75 
76 
77 // ************************************************************************* //
complexField ComplexField(const UList< scalar > &re, const UList< scalar > &im)
Definition: complexFields.C:42
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
vectorField U()
Generate and return a velocity field.
Definition: turbGen.C:47
static tmp< complexField > reverseTransform(const tmp< complexField > &field, const labelList &nn)
Definition: fft.C:203
label k
Boltzmann constant.
tmp< scalarField > Ek(const scalar Ea, const scalar k0, const scalarField &k)
Definition: Ek.H:10
CGAL::Exact_predicates_exact_constructions_kernel K
scalarField ReImSum(const UList< complex > &cf)
Definition: complexFields.C:84
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Type sample01()
Advance the state and return a sample of a given type from a.
Definition: RandomI.H:70
dimensionedScalar cos(const dimensionedScalar &ds)
Calculate the wavenumber vector field corresponding to the space vector field of a finite volume mesh...
Definition: Kmesh.H:49
const scalar twoPi(2 *pi)
dimensionedScalar sin(const dimensionedScalar &ds)
Specialisations of Field<T> for scalar, vector and tensor.
turbGen(const Kmesh &k, const scalar EA, const scalar K0)
Construct from components.
Definition: turbGen.C:36
dimensioned< scalar > mag(const dimensioned< Type > &)
scalar Ea(const scalar p, const scalar T) const
Definition: HtoEthermo.H:20
scalar scalar01()
Advance the state and return a scalar sample from a uniform.
Definition: RandomI.H:57