standardRadiation.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-2013 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 "standardRadiation.H"
27 #include "volFields.H"
30 
31 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace regionModels
36 {
37 namespace surfaceFilmModels
38 {
39 
40 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
41 
43 
45 (
48  dictionary
49 );
50 
51 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
52 
54 (
55  surfaceFilmModel& owner,
56  const dictionary& dict
57 )
58 :
59  filmRadiationModel(typeName, owner, dict),
60  QinPrimary_
61  (
62  IOobject
63  (
64  "Qin", // same name as Qin on primary region to enable mapping
65  owner.time().timeName(),
66  owner.regionMesh(),
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(),
83  ),
84  owner.regionMesh(),
85  dimensionedScalar("zero", dimMass/pow3(dimTime), 0.0),
86  zeroGradientFvPatchScalarField::typeName
87  ),
88  beta_(readScalar(coeffDict_.lookup("beta"))),
89  kappaBar_(readScalar(coeffDict_.lookup("kappaBar")))
90 {}
91 
92 
93 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
94 
96 {}
97 
98 
99 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
100 
102 {
103  // Transfer Qr from primary region
105 }
106 
107 
109 {
111  (
112  new volScalarField
113  (
114  IOobject
115  (
116  typeName + ":Shs",
117  owner().time().timeName(),
118  owner().regionMesh(),
121  ),
122  owner().regionMesh(),
123  dimensionedScalar("zero", dimMass/pow3(dimTime), 0.0),
124  zeroGradientFvPatchScalarField::typeName
125  )
126  );
127 
128  scalarField& Shs = tShs();
129  const scalarField& QinP = QinPrimary_.internalField();
132 
133  Shs = beta_*QinP*alpha*(1.0 - exp(-kappaBar_*delta));
134 
135  // Update net Qr on local region
136  QrNet_.internalField() = QinP - Shs;
138 
139  return tShs;
140 }
141 
142 
143 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
144 
145 } // End namespace surfaceFilmModels
146 } // End namespace regionModels
147 } // End namespace Foam
148 
149 // ************************************************************************* //
Foam::regionModels::surfaceFilmModels::standardRadiation::beta_
scalar beta_
Beta coefficient.
Definition: standardRadiation.H:71
volFields.H
Foam::regionModels::surfaceFilmModels::standardRadiation::correct
virtual void correct()
Correct.
Definition: standardRadiation.C:101
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:216
Foam::regionModels::surfaceFilmModels::standardRadiation::Shs
virtual tmp< volScalarField > Shs()
Return the radiation sensible enthalpy source.
Definition: standardRadiation.C:108
Foam::regionModels::surfaceFilmModels::standardRadiation::QinPrimary_
volScalarField QinPrimary_
Radiative incident flux mapped from the primary region / [kg/s3].
Definition: standardRadiation.H:62
Foam::regionModels::regionModel::time
const Time & time() const
Return the reference to the time database.
Definition: regionModelI.H:37
Foam::regionModels::surfaceFilmModels::filmSubModelBase::owner
const surfaceFilmModel & owner() const
Return const access to the owner surface film model.
Definition: filmSubModelBaseI.H:37
Foam::regionModels::surfaceFilmModels::surfaceFilmModel::delta
virtual const volScalarField & delta() const =0
Return the film thickness [m].
Foam::exp
dimensionedScalar exp(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:252
Foam::IOobject::NO_WRITE
@ NO_WRITE
Definition: IOobject.H:118
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
Foam::regionModels::surfaceFilmModels::standardRadiation::QrNet_
volScalarField QrNet_
Remaining radiative flux after removing local contribution.
Definition: standardRadiation.H:65
filmRadiationModel
Base class for film radiation models.
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::pow3
dimensionedScalar pow3(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:87
Foam::IOobject::NO_READ
@ NO_READ
Definition: IOobject.H:111
Foam::regionModels::surfaceFilmModels::standardRadiation::~standardRadiation
virtual ~standardRadiation()
Destructor.
Definition: standardRadiation.C:95
delta
scalar delta
Definition: LISASMDCalcMethod2.H:8
standardRadiation
Standard radiation model.
Foam::GeometricField::internalField
InternalField & internalField()
Return internal field.
Definition: GeometricField.C:724
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:41
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:41
Foam::regionModels::singleLayerRegion::mappedPushedFieldPatchTypes
wordList mappedPushedFieldPatchTypes() const
Return boundary types for pushed mapped field patches.
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::dimMass
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
Foam::GeometricField::correctBoundaryConditions
void correctBoundaryConditions()
Correct boundary field.
Definition: GeometricField.C:885
Foam::regionModels::surfaceFilmModels::defineTypeNameAndDebug
defineTypeNameAndDebug(kinematicSingleLayer, 0)
Foam::readScalar
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if succesful.
Definition: doubleScalar.H:63
Foam::regionModels::surfaceFilmModels::standardRadiation::kappaBar_
scalar kappaBar_
Bar(kappa) coefficient.
Definition: standardRadiation.H:74
Foam::Time::timeName
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:741
Foam::regionModels::surfaceFilmModels::standardRadiation::standardRadiation
standardRadiation(const standardRadiation &)
Disallow default bitwise copy construct.
timeName
word timeName
Definition: getTimeIndex.H:3
Foam::regionModels::surfaceFilmModels::filmSubModelBase::owner_
surfaceFilmModel & owner_
Reference to the owner surface film model.
Definition: filmSubModelBase.H:63
Foam::regionModels::surfaceFilmModels::addToRunTimeSelectionTable
addToRunTimeSelectionTable(surfaceFilmModel, kinematicSingleLayer, mesh)
Foam::regionModels::surfaceFilmModels::surfaceFilmModel
Base class for surface film models.
Definition: surfaceFilmModel.H:61
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::regionModels::surfaceFilmModels::filmRadiationModel
Definition: filmRadiationModel.H:54
zeroGradientFvPatchFields.H
standardRadiation.H
Foam::regionModels::surfaceFilmModels::surfaceFilmModel::alpha
virtual const volScalarField & alpha() const =0
Return the film coverage, 1 = covered, 0 = uncovered / [].
Foam::regionModels::regionModel::regionMesh
const fvMesh & regionMesh() const
Return the region mesh database.
Definition: regionModelI.H:61