OpenFOAM
v2512
The open source CFD toolbox
Loading...
Searching...
No Matches
alphaCourantNo.H
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-2017 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
Global
27
alphaCourantNo
28
29
Description
30
Calculates and outputs the mean and maximum Courant Numbers.
31
32
\*---------------------------------------------------------------------------*/
33
34
scalar
maxAlphaCo
35
(
36
runTime
.controlDict().get<scalar>(
"maxAlphaCo"
)
37
);
38
39
scalar
alphaCoNum
= 0.0;
40
scalar
meanAlphaCoNum
= 0.0;
41
42
if
(
mesh
.nInternalFaces())
43
{
44
scalarField sumPhi
45
(
46
pos0(
alpha1
- 0.01)*pos0(0.99 -
alpha1
)
47
*fvc::surfaceSum(mag(
phi
))().primitiveField()
48
);
49
50
alphaCoNum
= 0.5*gMax(sumPhi/
mesh
.V().field())*
runTime
.deltaTValue();
51
52
meanAlphaCoNum
=
53
0.5*(gSum(sumPhi)/gSum(
mesh
.V().field()))*
runTime
.deltaTValue();
54
}
55
56
Info
<<
"Interface Courant Number mean: "
<<
meanAlphaCoNum
57
<<
" max: "
<<
alphaCoNum
<<
endl
;
58
59
// ************************************************************************* //
alpha1
const volScalarField & alpha1
Definition
setRegionFluidFields.H:6
phi
phi
Definition
correctPhiFaceMask.H:34
mesh
dynamicFvMesh & mesh
Definition
createDynamicFvMesh.H:6
runTime
engineTime & runTime
Definition
createEngineTime.H:13
meanAlphaCoNum
scalar meanAlphaCoNum
Definition
alphaCourantNo.H:42
alphaCoNum
scalar alphaCoNum
Definition
alphaCourantNo.H:41
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere).
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition
Ostream.H:519
maxAlphaCo
scalar maxAlphaCo
Definition
porousAlphaCourantNo.H:25
applications
solvers
multiphase
twoLiquidMixingFoam
alphaCourantNo.H
Generated by
1.16.1