FreeFOAM The Cross-Platform CFD Toolkit
Hosted by SourceForge:
Get FreeFOAM at SourceForge.net.
            Fast, secure and Free Open Source software downloads

dnsFoam.C

Go to the documentation of this file.
00001 /*---------------------------------------------------------------------------*\
00002   =========                 |
00003   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
00004    \\    /   O peration     |
00005     \\  /    A nd           | Copyright (C) 1991-2010 OpenCFD Ltd.
00006      \\/     M anipulation  |
00007 -------------------------------------------------------------------------------
00008 License
00009     This file is part of OpenFOAM.
00010 
00011     OpenFOAM is free software: you can redistribute it and/or modify it
00012     under the terms of the GNU General Public License as published by
00013     the Free Software Foundation, either version 3 of the License, or
00014     (at your option) any later version.
00015 
00016     OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
00017     ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
00018     FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
00019     for more details.
00020 
00021     You should have received a copy of the GNU General Public License
00022     along with OpenFOAM.  If not, see <http://www.gnu.org/licenses/>.
00023 
00024 Application
00025     dnsFoam
00026 
00027 Description
00028     Direct numerical simulation solver for boxes of isotropic turbulence
00029 
00030 Usage
00031     - dnsFoam [OPTION]
00032 
00033     @param -case <dir> \n
00034     Specify the case directory
00035 
00036     @param -parallel \n
00037     Run the case in parallel
00038 
00039     @param -help \n
00040     Display short usage message
00041 
00042     @param -doc \n
00043     Display Doxygen documentation page
00044 
00045     @param -srcDoc \n
00046     Display source code
00047 
00048 \*---------------------------------------------------------------------------*/
00049 
00050 #include <finiteVolume/fvCFD.H>
00051 #include <randomProcesses/Kmesh.H>
00052 #include <randomProcesses/UOprocess.H>
00053 #include <randomProcesses/fft.H>
00054 #include <randomProcesses/calcEk.H>
00055 #include <OpenFOAM/graph.H>
00056 
00057 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00058 
00059 int main(int argc, char *argv[])
00060 {
00061     #include <OpenFOAM/setRootCase.H>
00062 
00063     #include <OpenFOAM/createTime.H>
00064     #include <OpenFOAM/createMeshNoClear.H>
00065     #include "readTransportProperties.H"
00066     #include "createFields.H"
00067     #include "readTurbulenceProperties.H"
00068     #include <finiteVolume/initContinuityErrs.H>
00069 
00070     // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00071 
00072     Info<< nl << "Starting time loop" << endl;
00073 
00074     while (runTime.loop())
00075     {
00076         Info<< "Time = " << runTime.timeName() << nl << endl;
00077 
00078         #include <finiteVolume/readPISOControls.H>
00079 
00080         force.internalField() = ReImSum
00081         (
00082             fft::reverseTransform
00083             (
00084                 K/(mag(K) + 1.0e-6) ^ forceGen.newField(), K.nn()
00085             )
00086         );
00087 
00088         #include "globalProperties.H"
00089 
00090         fvVectorMatrix UEqn
00091         (
00092             fvm::ddt(U)
00093           + fvm::div(phi, U)
00094           - fvm::laplacian(nu, U)
00095          ==
00096             force
00097         );
00098 
00099         solve(UEqn == -fvc::grad(p));
00100 
00101 
00102         // --- PISO loop
00103 
00104         for (int corr=1; corr<=1; corr++)
00105         {
00106             volScalarField rUA = 1.0/UEqn.A();
00107 
00108             U = rUA*UEqn.H();
00109             phi = (fvc::interpolate(U) & mesh.Sf())
00110                 + fvc::ddtPhiCorr(rUA, U, phi);
00111 
00112             fvScalarMatrix pEqn
00113             (
00114                 fvm::laplacian(rUA, p) == fvc::div(phi)
00115             );
00116 
00117             pEqn.solve();
00118 
00119             phi -= pEqn.flux();
00120 
00121             #include <finiteVolume/continuityErrs.H>
00122 
00123             U -= rUA*fvc::grad(p);
00124             U.correctBoundaryConditions();
00125         }
00126 
00127         runTime.write();
00128 
00129         if (runTime.outputTime())
00130         {
00131             calcEk(U, K).write(runTime.timePath()/"Ek", runTime.graphFormat());
00132         }
00133 
00134         Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
00135             << "  ClockTime = " << runTime.elapsedClockTime() << " s"
00136             << nl << endl;
00137     }
00138 
00139     Info<< "End\n" << endl;
00140 
00141     return 0;
00142 }
00143 
00144 
00145 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines