openfoam/applications/utilities/preProcessing/boxTurb/boxTurb.C
Mark Olesen 4b60453cf1 use while (runTime.loop() { .. } where possible in solvers
- change system/controlDict to use functions {..} instead of functions (..);
  * This is internally more efficient
- fixed formatting of system/controlDict functions entry

- pedantic change: use 'return 0' instead of 'return(0)' in the applications,
  since return is a C/C++ keyword, not a function.
2009-02-18 08:57:10 +01:00

80 lines
2.3 KiB
C

/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 1991-2009 OpenCFD Ltd.
\\/ M anipulation |
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
OpenFOAM is free software; you can redistribute it and/or modify it
under the terms of the GNU General Public License as published by the
Free Software Foundation; either version 2 of the License, or (at your
option) any later version.
OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
for more details.
You should have received a copy of the GNU General Public License
along with OpenFOAM; if not, write to the Free Software Foundation,
Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
Application
boxTurb3D
Description
Makes a box of turbulence which conforms to a given energy
spectrum and is divergence free.
\*---------------------------------------------------------------------------*/
#include "fvCFD.H"
#include "graph.H"
#include "OFstream.H"
#include "Kmesh.H"
#include "turbGen.H"
#include "calcEk.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
int main(int argc, char *argv[])
{
argList::noParallel();
# include "setRootCase.H"
# include "createTime.H"
# include "createMesh.H"
# include "createFields.H"
# include "readBoxTurbDict.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
Kmesh K(mesh);
turbGen Ugen(K, Ea, k0);
U.internalField() = Ugen.U();
U.correctBoundaryConditions();
Info<< "k("
<< runTime.timeName()
<< ") = "
<< 3.0/2.0*average(magSqr(U)).value() << endl;
U.write();
calcEk(U, K).write(runTime.timePath()/"Ek", runTime.graphFormat());
Info<< "end" << endl;
return 0;
}
// ************************************************************************* //