initial commit for version 1.5.x patch release
[OpenFOAM-1.5.x.git] / applications / solvers / heatTransfer / chtMultiRegionFoam / fluid / compressibleCourantNo.C
blob9e7aacd9b67cd5ffe71a84cf6d3beefa1763a81c
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2008 OpenCFD Ltd.
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
9     This file is part of OpenFOAM.
11     OpenFOAM is free software; you can redistribute it and/or modify it
12     under the terms of the GNU General Public License as published by the
13     Free Software Foundation; either version 2 of the License, or (at your
14     option) any later version.
16     OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17     ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18     FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
19     for more details.
21     You should have received a copy of the GNU General Public License
22     along with OpenFOAM; if not, write to the Free Software Foundation,
23     Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 Description
26     Calculates and outputs the mean and maximum Courant Numbers for the fluid
27     regions
29 \*---------------------------------------------------------------------------*/
31 scalar compressibleCourantNo
33     const fvMesh& mesh,
34     const Time& runTime,
35     const volScalarField& rho,
36     const surfaceScalarField& phi
39     scalar CoNum = 0.0;
40     scalar meanCoNum = 0.0;
42     if (mesh.nInternalFaces())
43     {
44         surfaceScalarField SfUfbyDelta =
45             mesh.surfaceInterpolation::deltaCoeffs()
46           * mag(phi)
47           / fvc::interpolate(rho);
49         CoNum = max(SfUfbyDelta/mesh.magSf())
50             .value()*runTime.deltaT().value();
52         meanCoNum = (sum(SfUfbyDelta)/sum(mesh.magSf()))
53             .value()*runTime.deltaT().value();
54     }
56     Info<< "Region: " << mesh.name() << " Courant Number mean: " << meanCoNum
57         << " max: " << CoNum << endl;
59     return CoNum;
63 // ************************************************************************* //