solidRegionDiffNo.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 | Copyright (C) 2011-2017 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 "solidRegionDiffNo.H"
27 #include "surfaceInterpolate.H"
28 
29 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
30 
31 Foam::scalar Foam::solidRegionDiffNo
32 (
33  const fvMesh& mesh,
34  const Time& runTime,
35  const volScalarField& Cprho,
36  const volScalarField& kappa
37 )
38 {
39  surfaceScalarField kapparhoCpbyDelta
40  (
41  sqr(mesh.surfaceInterpolation::deltaCoeffs())
42  *fvc::interpolate(kappa)
43  /fvc::interpolate(Cprho)
44  );
45 
46  const scalar DiNum = max(kapparhoCpbyDelta).value()*runTime.deltaTValue();
47  const scalar meanDiNum =
48  average(kapparhoCpbyDelta).value()*runTime.deltaTValue();
49 
50  Info<< "Region: " << mesh.name() << " Diffusion Number mean: " << meanDiNum
51  << " max: " << DiNum << endl;
52 
53  return DiNum;
54 }
55 
56 // ************************************************************************* //
scalar DiNum
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
Calculates and outputs the mean and maximum Diffusion Numbers for the solid regions.
scalar solidRegionDiffNo(const fvMesh &mesh, const Time &runTime, const volScalarField &Cprho, const volScalarField &kappa)
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
messageStream Info
tmp< GeometricField< Type, fvPatchField, volMesh > > average(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Area-weighted average a surfaceField creating a volField.
Definition: fvcAverage.C:46
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField