createBoxTurb.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) 2018 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 Application
27  createBoxTurb
28 
29 Description
30  Create a box of isotropic turbulence based on a user-specified
31  energy spectrum.
32 
33  Based on the reference
34  \verbatim
35  Saad, T., Cline, D., Stoll, R., Sutherland, J.C.
36  "Scalable Tools for Generating Synthetic Isotropic Turbulence with
37  Arbitrary Spectra"
38  AIAA Journal, Vol. 55, No. 1 (2017), pp. 327-331.
39  \endverbatim
40 
41  The \c -createBlockMesh option creates a block mesh and exits, which
42  can then be decomposed and the utility run in parallel.
43 
44 \*---------------------------------------------------------------------------*/
45 
46 #include "fvCFD.H"
47 #include "block.H"
48 #include "mathematicalConstants.H"
49 
50 using namespace Foam::constant;
51 
52 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
53 
54 Foam::vector randomUnitVector(Random& rndGen)
55 {
56  // Sample point on a sphere
57  scalar t = rndGen.globalPosition<scalar>(-1, 1);
58  scalar phim = rndGen.globalSample01<scalar>()*mathematical::twoPi;
59  scalar thetam = Foam::acos(t);
60 
61  return vector
62  (
63  Foam::sin(thetam*Foam::cos(phim)),
64  Foam::sin(thetam*Foam::sin(phim)),
65  Foam::cos(thetam)
66  );
67 }
68 
69 
70 int main(int argc, char *argv[])
71 {
73  (
74  "Create a box of isotropic turbulence based on a user-specified"
75  " energy spectrum."
76  );
77 
79  (
80  "createBlockMesh",
81  "create the block mesh and exit"
82  );
83 
84  #include "setRootCase.H"
85 
86  #include "createTime.H"
87  #include "createFields.H"
88 
89  if (args.found("createBlockMesh"))
90  {
91  // Create a box block mesh with cyclic patches
92  #include "createBlockMesh.H"
93  Info<< "\nEnd\n" << endl;
94  return 0;
95  }
96 
97  #include "createMesh.H"
98 
99  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
100 
101  // Minimum wave number
102  scalar kappa0 = mathematical::twoPi/cmptMin(L);
103 
104  // Maximum wave number
105  scalar kappaMax = mathematical::pi/cmptMin(delta);
106 
107  Info<< "Wave number min/max = " << kappa0 << ", " << kappaMax << endl;
108 
109  Info<< "Generating velocity field" << endl;
110 
112  (
113  IOobject
114  (
115  "U",
116  runTime.timeName(),
117  mesh,
120  ),
121  mesh,
123  );
124 
125  vectorField& Uc = U.primitiveFieldRef();
126  const scalar deltaKappa = (kappaMax - kappa0)/scalar(nModes - 1);
127  const vectorField& C(mesh.C());
128  for (label modei = 1; modei <= nModes; ++modei)
129  {
130  // Equidistant wave mode
131  scalar kappaM = kappa0 + deltaKappa*(modei-1);
132 
133  Info<< "Processing mode:" << modei << " kappaM:" << kappaM << endl;
134 
135  // Energy
136  scalar E = Ek->value(kappaM);
137 
138  // Wave amplitude
139  scalar qm = Foam::sqrt(E*deltaKappa);
140 
141  // Wave number unit vector
142  const vector kappaHatm(randomUnitVector(rndGen));
143 
144  vector kappaTildem(0.5*kappaM*cmptMultiply(kappaHatm, delta));
145  for (direction i = 0; i < 3; ++i)
146  {
147  kappaTildem[i] = 2/delta[i]*Foam::sin(kappaTildem[i]);
148  }
149 
150  // Intermediate unit vector zeta
151  const vector zetaHatm(randomUnitVector(rndGen));
152 
153  // Unit vector sigma
154  vector sigmaHatm(zetaHatm^kappaTildem);
155  sigmaHatm /= mag(kappaTildem);
156 
157  // Phase angle
158  scalar psim = 0.5*rndGen.position(-mathematical::pi, mathematical::pi);
159 
160  // Add the velocity contribution per mode
161  Uc += 2*qm*cos(kappaM*(kappaHatm & C) + psim)*sigmaHatm;
162  }
163 
164  U.write();
165 
166  {
167  Info<< "Generating kinetic energy field" << endl;
168  volScalarField k("k", 0.5*magSqr(U));
169  k.write();
170  Info<< "min/max/average k = "
171  << gMin(k) << ", " << gMax(k) << ", " << gAverage(k)
172  << endl;
173  }
174 
175  {
176  Info<< "Generating div(U) field" << endl;
178  divU.write();
179  Info<< "min/max/average div(U) = "
180  << gMin(divU) << ", " << gMax(divU) << ", " << gAverage(divU)
181  << endl;
182  }
183 
184  Info<< nl;
185  runTime.printExecutionTime(Info);
186 
187  Info<< "End\n" << endl;
188  return 0;
189 }
190 
191 
192 // ************************************************************************* //
Foam::IOobject::NO_WRITE
Definition: IOobject.H:195
runTime
engineTime & runTime
Definition: createEngineTime.H:13
L
const vector L(dict.get< vector >("L"))
mathematicalConstants.H
Foam::cmptMultiply
dimensioned< Type > cmptMultiply(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::dimVelocity
const dimensionSet dimVelocity
Foam::gAverage
Type gAverage(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:604
Foam::constant
Different types of constants.
Definition: atomicConstants.C:38
Foam::sin
dimensionedScalar sin(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:264
Foam::argList::addNote
static void addNote(const string &note)
Add extra notes for the usage information.
Definition: argList.C:412
Foam::fvc::div
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:49
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
C
volScalarField & C
Definition: readThermalProperties.H:102
Foam::vectorField
Field< vector > vectorField
Specialisation of Field<T> for vector.
Definition: primitiveFieldsFwd.H:54
Foam::magSqr
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
Foam::dimensionedVector
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
Definition: dimensionedVector.H:50
Foam::cmptMin
void cmptMin(FieldField< Field, typename FieldField< Field, Type >::cmptType > &cf, const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:302
createBlockMesh.H
nModes
const label nModes
Definition: createFields.H:22
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
delta
scalar delta
Definition: LISASMDCalcMethod2.H:8
divU
zeroField divU
Definition: alphaSuSp.H:3
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
Foam::constant::mathematical::twoPi
constexpr scalar twoPi(2 *M_PI)
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::Ek
tmp< scalarField > Ek(const scalar Ea, const scalar k0, const scalarField &k)
Definition: Ek.H:10
Foam::vector
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:51
Foam::volVectorField
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:62
Foam::argList::addBoolOption
static void addBoolOption(const word &optName, const string &usage="", bool advanced=false)
Add a bool option to validOptions with usage information.
Definition: argList.C:324
U
U
Definition: pEqn.H:72
setRootCase.H
Foam::constant::mathematical::pi
constexpr scalar pi(M_PI)
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Foam::Vector< scalar >
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:144
Foam::acos
dimensionedScalar acos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:268
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
createMesh.H
Required Variables.
k
label k
Boltzmann constant.
Definition: LISASMDCalcMethod2.H:41
Foam::direction
uint8_t direction
Definition: direction.H:52
createTime.H
rndGen
Random rndGen
Definition: createFields.H:23
fvCFD.H
Foam::gMin
Type gMin(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:593
Foam::IOobject::NO_READ
Definition: IOobject.H:188
args
Foam::argList args(argc, argv)
block.H
Foam::gMax
Type gMax(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:592
Foam::argList::found
bool found(const word &optName) const
Return true if the named option is found.
Definition: argListI.H:178
Foam::cos
dimensionedScalar cos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:265