mirror of
https://develop.openfoam.com/Development/openfoam.git
synced 2025-12-28 03:37:59 +00:00
59 lines
2.1 KiB
C
59 lines
2.1 KiB
C
/*---------------------------------------------------------------------------*\
|
|
========= |
|
|
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
|
\\ / O peration |
|
|
\\ / A nd | www.openfoam.com
|
|
\\/ M anipulation |
|
|
-------------------------------------------------------------------------------
|
|
Copyright (C) 2011-2017 OpenFOAM Foundation
|
|
-------------------------------------------------------------------------------
|
|
License
|
|
This file is part of OpenFOAM.
|
|
|
|
OpenFOAM is free software: you can redistribute it and/or modify it
|
|
under the terms of the GNU General Public License as published by
|
|
the Free Software Foundation, either version 3 of the License, or
|
|
(at your option) any later version.
|
|
|
|
OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
|
|
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
|
|
FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
|
|
for more details.
|
|
|
|
You should have received a copy of the GNU General Public License
|
|
along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
|
|
|
|
\*---------------------------------------------------------------------------*/
|
|
|
|
#include "solidRegionDiffNo.H"
|
|
#include "surfaceInterpolate.H"
|
|
|
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
|
|
|
Foam::scalar Foam::solidRegionDiffNo
|
|
(
|
|
const fvMesh& mesh,
|
|
const Time& runTime,
|
|
const volScalarField& Cprho,
|
|
const volScalarField& kappa
|
|
)
|
|
{
|
|
surfaceScalarField kapparhoCpbyDelta
|
|
(
|
|
sqr(mesh.surfaceInterpolation::deltaCoeffs())
|
|
*fvc::interpolate(kappa)
|
|
/fvc::interpolate(Cprho)
|
|
);
|
|
|
|
const scalar DiNum = max(kapparhoCpbyDelta).value()*runTime.deltaTValue();
|
|
const scalar meanDiNum =
|
|
average(kapparhoCpbyDelta).value()*runTime.deltaTValue();
|
|
|
|
Info<< "Region: " << mesh.name() << " Diffusion Number mean: " << meanDiNum
|
|
<< " max: " << DiNum << endl;
|
|
|
|
return DiNum;
|
|
}
|
|
|
|
// ************************************************************************* //
|