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-2016 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 "fvc.H"
28 
29 Foam::scalar Foam::solidRegionDiffNo
30 (
31  const fvMesh& mesh,
32  const Time& runTime,
33  const volScalarField& Cprho,
34  const volScalarField& kappa
35 )
36 {
37  scalar DiNum = 0.0;
38  scalar meanDiNum = 0.0;
39 
40  //- Take care: can have fluid domains with 0 cells so do not test for
41  // zero internal faces.
42  surfaceScalarField kapparhoCpbyDelta
43  (
44  mesh.surfaceInterpolation::deltaCoeffs()
45  * fvc::interpolate(kappa)
46  / fvc::interpolate(Cprho)
47  );
48 
49  DiNum = gMax(kapparhoCpbyDelta.primitiveField())*runTime.deltaT().value();
50 
51  meanDiNum = (average(kapparhoCpbyDelta)).value()*runTime.deltaT().value();
52 
53  Info<< "Region: " << mesh.name() << " Diffusion Number mean: " << meanDiNum
54  << " max: " << DiNum << endl;
55 
56  return DiNum;
57 }
58 
59 // ************************************************************************* //
scalar DiNum
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
tmp< surfaceScalarField > interpolate(const RhoType &rho)
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)
Type gMax(const FieldField< Field, Type > &f)
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