boxTurb.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-2016 OpenFOAM Foundation
9-------------------------------------------------------------------------------
10License
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
26Application
27 boxTurb
28
29Group
30 grpPreProcessingUtilities
31
32Description
33 Create a box of divergence-free turbulence conforming to a given
34 energy spectrum.
35
36\*---------------------------------------------------------------------------*/
37
38#include "fvCFD.H"
39#include "graph.H"
40#include "OFstream.H"
41#include "Kmesh.H"
42#include "turbGen.H"
43#include "calcEk.H"
44
45
46// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
47
48int main(int argc, char *argv[])
49{
50 argList::addNote
51 (
52 "Create a box of divergence-free turbulence conforming to a given"
53 " energy spectrum"
54 );
55
56 argList::noParallel();
57 #include "setRootCase.H"
58
59 #include "createTime.H"
60 #include "createNamedMesh.H"
61 #include "createFields.H"
62 #include "readBoxTurbDict.H"
63
64
65 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
66
67 Kmesh K(mesh);
68
69 turbGen Ugen(K, Ea, k0);
70
71 U.primitiveFieldRef() = Ugen.U();
72 U.correctBoundaryConditions();
73
74 Info<< "k("
75 << runTime.timeName()
76 << ") = "
77 << 3.0/2.0*average(magSqr(U)).value() << endl;
78
79 U.write();
80
81 calcEk(U, K).write
82 (
83 runTime.path()/"graphs"/runTime.timeName(),
84 "Ek",
85 runTime.graphFormat()
86 );
87
88 Info<< "End\n" << endl;
89
90 return 0;
91}
92
93
94// ************************************************************************* //
CGAL::Exact_predicates_exact_constructions_kernel K
scalar Ea(const scalar p, const scalar T) const
Definition: HtoEthermo.H:32
void write(Ostream &, const word &format) const
Write graph to stream in given format.
Definition: graph.C:266
U
Definition: pEqn.H:72
dynamicFvMesh & mesh
engineTime & runTime
Required Variables.
tmp< GeometricField< Type, faPatchField, areaMesh > > average(const GeometricField< Type, faePatchField, edgeMesh > &ssf)
Area-weighted average a edgeField creating a areaField.
Definition: facAverage.C:47
messageStream Info
Information stream (stdout output on master, null elsewhere)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
graph calcEk(const volVectorField &U, const Kmesh &K)
Definition: calcEk.C:27
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)