BUG: UListIO: byteSize overflowing on really big faceLists
[OpenFOAM-2.0.x.git] / src / regionModels / surfaceFilmModels / submodels / thermo / filmRadiationModel / standardRadiation / standardRadiation.C
blobba1585713005857447375182f7aab00417f6bf5a
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 "standardRadiation.H"
27 #include "volFields.H"
28 #include "addToRunTimeSelectionTable.H"
29 #include "zeroGradientFvPatchFields.H"
31 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 namespace Foam
35 namespace regionModels
37 namespace surfaceFilmModels
40 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
42 defineTypeNameAndDebug(standardRadiation, 0);
44 addToRunTimeSelectionTable
46     filmRadiationModel,
47     standardRadiation,
48     dictionary
51 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
53 standardRadiation::standardRadiation
55     const surfaceFilmModel& owner,
56     const dictionary& dict
59     filmRadiationModel(typeName, owner, dict),
60     QrPrimary_
61     (
62         IOobject
63         (
64             "Qr", // same name as Qr on primary region to enable mapping
65             owner.time().timeName(),
66             owner.regionMesh(),
67             IOobject::NO_READ,
68             IOobject::NO_WRITE
69         ),
70         owner.regionMesh(),
71         dimensionedScalar("zero", dimMass/pow3(dimTime), 0.0),
72         owner.mappedPushedFieldPatchTypes<scalar>()
73     ),
74     QrNet_
75     (
76         IOobject
77         (
78             "QrNet",
79             owner.time().timeName(),
80             owner.regionMesh(),
81             IOobject::NO_READ,
82             IOobject::NO_WRITE
83         ),
84         owner.regionMesh(),
85         dimensionedScalar("zero", dimMass/pow3(dimTime), 0.0),
86         zeroGradientFvPatchScalarField::typeName
87     ),
88     delta_(owner.delta()),
89     deltaMin_(readScalar(coeffs_.lookup("deltaMin"))),
90     beta_(readScalar(coeffs_.lookup("beta"))),
91     kappaBar_(readScalar(coeffs_.lookup("kappaBar")))
95 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
97 standardRadiation::~standardRadiation()
101 // * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * * //
103 void standardRadiation::correct()
105     // Transfer Qr from primary region
106     QrPrimary_.correctBoundaryConditions();
110 tmp<volScalarField> standardRadiation::Shs()
112     tmp<volScalarField> tShs
113     (
114         new volScalarField
115         (
116             IOobject
117             (
118                 typeName + "::Shs",
119                 owner().time().timeName(),
120                 owner().regionMesh(),
121                 IOobject::NO_READ,
122                 IOobject::NO_WRITE
123             ),
124             owner().regionMesh(),
125             dimensionedScalar("zero", dimMass/pow3(dimTime), 0.0),
126             zeroGradientFvPatchScalarField::typeName
127         )
128     );
130     scalarField& Shs = tShs();
131     const scalarField& QrP = QrPrimary_.internalField();
132     const scalarField& delta = delta_.internalField();
134     Shs = beta_*(QrP*pos(delta - deltaMin_))*(1.0 - exp(-kappaBar_*delta));
136     // Update net Qr on local region
137     QrNet_.internalField() = QrP - Shs;
138     QrNet_.correctBoundaryConditions();
140     return tShs;
144 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
146 } // End namespace surfaceFilmModels
147 } // End namespace regionModels
148 } // End namespace Foam
150 // ************************************************************************* //