compressibleCourantNo.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 -------------------------------------------------------------------------------
10 License
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 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "compressibleCourantNo.H"
29 #include "fvc.H"
30 
31 Foam::scalar Foam::compressibleCourantNo
32 (
33  const fvMesh& mesh,
34  const Time& runTime,
35  const volScalarField& rho,
36  const surfaceScalarField& phi
37 )
38 {
39  scalarField sumPhi
40  (
41  fvc::surfaceSum(mag(phi))().primitiveField()
42  / rho.primitiveField()
43  );
44 
45  scalar CoNum = 0.5*gMax(sumPhi/mesh.V().field())*runTime.deltaTValue();
46 
47  scalar meanCoNum =
48  0.5*(gSum(sumPhi)/gSum(mesh.V().field()))*runTime.deltaTValue();
49 
50  Info<< "Region: " << mesh.name() << " Courant Number mean: " << meanCoNum
51  << " max: " << CoNum << endl;
52 
53  return CoNum;
54 }
55 
56 
57 // ************************************************************************* //
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:487
scalar compressibleCourantNo(const fvMesh &mesh, const Time &runTime, const volScalarField &rho, const surfaceScalarField &phi)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:84
dynamicFvMesh & mesh
Type gSum(const FieldField< Field, Type > &f)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
scalar CoNum
Type gMax(const FieldField< Field, Type > &f)
scalar meanCoNum
messageStream Info
Information stream (stdout output on master, null elsewhere)
tmp< GeometricField< Type, fvPatchField, volMesh > > surfaceSum(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField