compressibleCourantNo.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2018 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
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
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
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.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "compressibleCourantNo.H"
27 #include "fvc.H"
28 
29 Foam::scalar Foam::compressibleCourantNo
30 (
31  const fvMesh& mesh,
32  const Time& runTime,
33  const volScalarField& rho,
34  const surfaceScalarField& phi
35 )
36 {
37  scalarField sumPhi
38  (
39  fvc::surfaceSum(mag(phi))().primitiveField()
40  / rho.primitiveField()
41  );
42 
43  scalar CoNum = 0.5*gMax(sumPhi/mesh.V().field())*runTime.deltaTValue();
44 
45  scalar meanCoNum =
46  0.5*(gSum(sumPhi)/gSum(mesh.V().field()))*runTime.deltaTValue();
47 
48  Info<< "Region: " << mesh.name() << " Courant Number mean: " << meanCoNum
49  << " max: " << CoNum << endl;
50 
51  return CoNum;
52 }
53 
54 
55 // ************************************************************************* //
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
scalar compressibleCourantNo(const fvMesh &mesh, const Time &runTime, const volScalarField &rho, const surfaceScalarField &phi)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
Type gSum(const FieldField< Field, Type > &f)
scalar CoNum
volScalarField scalarField(fieldObject, mesh)
Type gMax(const FieldField< Field, Type > &f)
messageStream Info
dimensioned< scalar > mag(const dimensioned< Type > &)
tmp< GeometricField< Type, fvPatchField, volMesh > > surfaceSum(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField