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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2017 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "standardRadiation.H"
29 #include "volFields.H"
32 
33 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 namespace regionModels
38 {
39 namespace surfaceFilmModels
40 {
41 
42 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
43 
44 defineTypeNameAndDebug(standardRadiation, 0);
45 
47 (
48  filmRadiationModel,
49  standardRadiation,
50  dictionary
51 );
52 
53 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
54 
55 standardRadiation::standardRadiation
56 (
58  const dictionary& dict
59 )
60 :
61  filmRadiationModel(typeName, film, dict),
62  qinPrimary_
63  (
64  IOobject
65  (
66  "qin", // same name as qin on primary region to enable mapping
67  film.time().timeName(),
68  film.regionMesh(),
71  ),
72  film.regionMesh(),
74  film.mappedPushedFieldPatchTypes<scalar>()
75  ),
76  qrNet_
77  (
78  IOobject
79  (
80  "qrNet",
81  film.time().timeName(),
82  film.regionMesh(),
85  ),
86  film.regionMesh(),
88  zeroGradientFvPatchScalarField::typeName
89  ),
90  beta_(coeffDict_.get<scalar>("beta")),
91  kappaBar_(coeffDict_.get<scalar>("kappaBar"))
92 {}
93 
94 
95 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
96 
98 {
99  // Transfer qr from primary region
100  qinPrimary_.correctBoundaryConditions();
101 }
102 
103 
105 {
107  (
108  new volScalarField
109  (
110  IOobject
111  (
112  typeName + ":Shs",
113  film().time().timeName(),
114  film().regionMesh(),
117  ),
118  film().regionMesh(),
120  )
121  );
122 
123  scalarField& Shs = tShs.ref();
124  const scalarField& qinP = qinPrimary_;
125  const scalarField& delta = filmModel_.delta();
126  const scalarField& alpha = filmModel_.alpha();
127 
128  Shs = beta_*qinP*alpha*(1.0 - exp(-kappaBar_*delta));
129 
130  // Update net qr on local region
131  qrNet_.primitiveFieldRef() = qinP - Shs;
132  qrNet_.correctBoundaryConditions();
133 
134  return tShs;
135 }
136 
137 
138 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
139 
140 } // End namespace surfaceFilmModels
141 } // End namespace regionModels
142 } // End namespace Foam
143 
144 // ************************************************************************* //
Foam::IOobject::NO_WRITE
Definition: IOobject.H:195
volFields.H
Foam::regionModels::surfaceFilmModels::standardRadiation::correct
virtual void correct()
Correct.
Definition: standardRadiation.C:97
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:212
Foam::regionModels::surfaceFilmModels::surfaceFilmRegionModel::delta
virtual const volScalarField & delta() const =0
Return the film thickness [m].
Foam::Time::timeName
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:780
Foam::regionModels::surfaceFilmModels::standardRadiation::Shs
virtual tmp< volScalarField > Shs()
Return the radiation sensible enthalpy source.
Definition: standardRadiation.C:104
Foam::regionModels::regionModel::time
const Time & time() const
Return the reference to the time database.
Definition: regionModelI.H:40
Foam::regionModels::surfaceFilmModels::filmSubModelBase::film
const surfaceFilmRegionModel & film() const
Return const access to the film surface film model.
Definition: filmSubModelBaseI.H:39
Foam::exp
dimensionedScalar exp(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:261
Foam::regionModels::surfaceFilmModels::addToRunTimeSelectionTable
addToRunTimeSelectionTable(surfaceFilmRegionModel, kinematicSingleLayer, mesh)
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:227
Foam::Field< scalar >
Foam::regionModels::surfaceFilmModels::surfaceFilmRegionModel
Base class for surface film models.
Definition: surfaceFilmRegionModel.H:55
Foam::pow3
dimensionedScalar pow3(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:89
Foam::regionModels::regionModel::regionMesh
const fvMesh & regionMesh() const
Return the region mesh database.
Definition: regionModelI.H:64
delta
scalar delta
Definition: LISASMDCalcMethod2.H:8
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:42
Foam::regionModels::surfaceFilmModels::surfaceFilmRegionModel::alpha
virtual const volScalarField & alpha() const =0
Return the film coverage, 1 = covered, 0 = uncovered / [].
timeName
word timeName
Definition: getTimeIndex.H:3
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::regionModels::singleLayerRegion::mappedPushedFieldPatchTypes
wordList mappedPushedFieldPatchTypes() const
Return boundary types for pushed mapped field patches.
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::dimMass
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:51
Foam::GeometricField::primitiveFieldRef
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field.
Definition: GeometricField.C:766
Foam::GeometricField::correctBoundaryConditions
void correctBoundaryConditions()
Correct boundary field.
Definition: GeometricField.C:940
Foam::regionModels::surfaceFilmModels::defineTypeNameAndDebug
defineTypeNameAndDebug(kinematicSingleLayer, 0)
Foam::regionModels::surfaceFilmModels::filmSubModelBase::filmModel_
surfaceFilmRegionModel & filmModel_
Reference to the film surface film model.
Definition: filmSubModelBase.H:65
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::IOobject::NO_READ
Definition: IOobject.H:188
Foam::regionModels::surfaceFilmModels::filmRadiationModel
Base class for film radiation models.
Definition: filmRadiationModel.H:56
zeroGradientFvPatchFields.H
standardRadiation.H