turbGen.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 OpenFOAM Foundation
9  Copyright (C) 2018 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 
29 #include "fft.H"
30 #include "turbGen.H"
31 #include "Kmesh.H"
32 #include "primitiveFields.H"
33 #include "Ek.H"
34 #include "mathematicalConstants.H"
35 
36 
37 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
38 
39 Foam::turbGen::turbGen(const Kmesh& k, const scalar EA, const scalar K0)
40 :
41  K(k),
42  Ea(EA),
43  k0(K0),
44  RanGen(label(0))
45 {}
46 
47 
48 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
49 
51 {
52  vectorField s(K.size());
53  scalarField rndPhases(K.size());
54 
55  forAll(K, i)
56  {
57  s[i] = RanGen.sample01<vector>();
58  rndPhases[i] = RanGen.sample01<scalar>();
59  }
60 
61  s = K ^ s;
62  s = s/(mag(s) + 1.0e-20);
63 
64  s = Ek(Ea, k0, mag(K))*s;
65 
66  label ntot = 1;
67  forAll(K.nn(), idim)
68  {
69  ntot *= K.nn()[idim];
70  }
71  const scalar recRootN = 1.0/sqrt(scalar(ntot));
72 
74  (
76  (
79  K.nn()
80  )*recRootN
81  );
82 
83  return ReImSum(up);
84 }
85 
86 
87 // ************************************************************************* //
Kmesh.H
mathematicalConstants.H
s
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.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))
Definition: gmvOutputSpray.H:25
Foam::sin
dimensionedScalar sin(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:264
Foam::ComplexField
complexField ComplexField(const UList< scalar > &re, const UList< scalar > &im)
Zip up two lists of values into a list of complex.
Definition: complexField.C:105
Foam::turbGen::turbGen
turbGen(const Kmesh &k, const scalar EA, const scalar K0)
Construct from components.
Definition: turbGen.C:39
primitiveFields.H
Specialisations of Field<T> for scalar, vector and tensor.
Foam::Kmesh
Calculate the wavenumber vector field corresponding to the space vector field of a finite volume mesh...
Definition: Kmesh.H:51
Ea
scalar Ea(const scalar p, const scalar T) const
Definition: HtoEthermo.H:32
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
K
CGAL::Exact_predicates_exact_constructions_kernel K
Definition: CGALTriangulation3DKernel.H:58
Foam::fft::reverseTransform
static tmp< complexField > reverseTransform(const tmp< complexField > &field, const UList< int > &nn)
Definition: fft.C:259
Foam::Field< vector >
Ek.H
Foam::constant::mathematical::twoPi
constexpr scalar twoPi(2 *M_PI)
Foam::Ek
tmp< scalarField > Ek(const scalar Ea, const scalar k0, const scalarField &k)
Definition: Ek.H:10
Foam::turbGen::U
vectorField U()
Generate and return a velocity field.
Definition: turbGen.C:50
Foam::ReImSum
scalarField ReImSum(const UList< complex > &cf)
Sum real and imag components.
Definition: complexField.C:146
Foam::Vector< scalar >
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:144
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
k
label k
Boltzmann constant.
Definition: LISASMDCalcMethod2.H:41
fft.H
turbGen.H
Foam::cos
dimensionedScalar cos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:265