BUG: UListIO: byteSize overflowing on really big faceLists
[OpenFOAM-2.0.x.git] / applications / solvers / heatTransfer / chtMultiRegionFoam / solid / solidRegionDiffNo.C
blob3dbab74858368248c4abe6ec2d507f011e2ac0ff
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2011 OpenFOAM Foundation
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
9     This file is part of OpenFOAM.
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.
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.
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/>.
24 \*---------------------------------------------------------------------------*/
26 #include "solidRegionDiffNo.H"
27 #include "fvc.H"
29 Foam::scalar Foam::solidRegionDiffNo
31     const fvMesh& mesh,
32     const Time& runTime,
33     const volScalarField& Cprho,
34     const volScalarField& K
37     scalar DiNum = 0.0;
38     scalar meanDiNum = 0.0;
40     //- Take care: can have fluid domains with 0 cells so do not test for
41     //  zero internal faces.
42     surfaceScalarField KrhoCpbyDelta
43     (
44         mesh.surfaceInterpolation::deltaCoeffs()
45       * fvc::interpolate(K)
46       / fvc::interpolate(Cprho)
47     );
49     DiNum = gMax(KrhoCpbyDelta.internalField())*runTime.deltaT().value();
51     meanDiNum = (average(KrhoCpbyDelta)).value()*runTime.deltaT().value();
53     Info<< "Region: " << mesh.name() << " Diffusion Number mean: " << meanDiNum
54         << " max: " << DiNum << endl;
56     return DiNum;
60 Foam::scalar Foam::solidRegionDiffNo
62     const fvMesh& mesh,
63     const Time& runTime,
64     const volScalarField& Cprho,
65     const volSymmTensorField& Kdirectional
68     scalar DiNum = 0.0;
69     scalar meanDiNum = 0.0;
71     volScalarField K(mag(Kdirectional));
73     //- Take care: can have fluid domains with 0 cells so do not test for
74     //  zero internal faces.
75     surfaceScalarField KrhoCpbyDelta
76     (
77         mesh.surfaceInterpolation::deltaCoeffs()
78       * fvc::interpolate(K)
79       / fvc::interpolate(Cprho)
80     );
82     DiNum = gMax(KrhoCpbyDelta.internalField())*runTime.deltaT().value();
84     meanDiNum = (average(KrhoCpbyDelta)).value()*runTime.deltaT().value();
86     Info<< "Region: " << mesh.name() << " Diffusion Number mean: " << meanDiNum
87         << " max: " << DiNum << endl;
89     return DiNum;
93 // ************************************************************************* //