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

compressibleCourantNo.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 \*---------------------------------------------------------------------------*/
00025 
00026 #include "compressibleCourantNo.H"
00027 #include <finiteVolume/fvc.H>
00028 
00029 Foam::scalar Foam::compressibleCourantNo
00030 (
00031     const fvMesh& mesh,
00032     const Time& runTime,
00033     const volScalarField& rho,
00034     const surfaceScalarField& phi
00035 )
00036 {
00037     scalar CoNum = 0.0;
00038     scalar meanCoNum = 0.0;
00039 
00040     //- Can have fluid domains with 0 cells so do not test.
00041     //if (mesh.nInternalFaces())
00042     {
00043         surfaceScalarField SfUfbyDelta =
00044             mesh.surfaceInterpolation::deltaCoeffs()
00045           * mag(phi)
00046           / fvc::interpolate(rho);
00047 
00048         CoNum = max(SfUfbyDelta/mesh.magSf())
00049             .value()*runTime.deltaT().value();
00050 
00051         meanCoNum = (sum(SfUfbyDelta)/sum(mesh.magSf()))
00052             .value()*runTime.deltaT().value();
00053     }
00054 
00055     Info<< "Region: " << mesh.name() << " Courant Number mean: " << meanCoNum
00056         << " max: " << CoNum << endl;
00057 
00058     return CoNum;
00059 }
00060 
00061 
00062 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines