Loading...
Searching...
No Matches
MarshakRadiationFvPatchScalarField.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-2015 OpenFOAM Foundation
9 Copyright (C) 2016,2020 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
31#include "fvPatchFieldMapper.H"
32#include "volFields.H"
33#include "radiationModel.H"
36
37// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
38
41(
42 const fvPatch& p,
44)
45:
46 mixedFvPatchScalarField(p, iF),
47 TName_("T")
49 refValue() = Zero;
50 refGrad() = Zero;
51 valueFraction() = 0.0;
52}
53
54
57(
59 const fvPatch& p,
61 const fvPatchFieldMapper& mapper
63:
64 mixedFvPatchScalarField(ptf, p, iF, mapper),
65 TName_(ptf.TName_)
66{}
67
68
71(
72 const fvPatch& p,
74 const dictionary& dict
75)
76:
77 mixedFvPatchScalarField(p, iF),
78 TName_(dict.getOrDefault<word>("T", "T"))
79{
80 if (!this->readValueEntry(dict))
81 {
83 }
84
85 refValue() = *this;
86 refGrad() = Zero; // zero gradient
87 valueFraction() = 1.0;
89}
90
91
94(
97:
98 mixedFvPatchScalarField(ptf),
99 TName_(ptf.TName_)
100{}
101
102
105(
108)
109:
110 mixedFvPatchScalarField(ptf, iF),
111 TName_(ptf.TName_)
112{}
113
114
115// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
116
118{
119 if (this->updated())
120 {
121 return;
122 }
123
124 // Since we're inside initEvaluate/evaluate there might be processor
125 // comms underway. Change the tag we use.
126 const int oldTag = UPstream::incrMsgType();
127
128 // Temperature field
129 const auto& Tp = patch().lookupPatchField<volScalarField>(TName_);
130
131 // Re-calc reference value
132 refValue() = 4.0*constant::physicoChemical::sigma.value()*pow4(Tp);
133
134 // Diffusion coefficient - created by radiation model's ::updateCoeffs()
135 const auto& gamma = patch().lookupPatchField<volScalarField>("gammaRad");
136
137 const boundaryRadiationProperties& boundaryRadiation =
138 boundaryRadiationProperties::New(internalField().mesh());
139
140 const tmp<scalarField> temissivity
141 (
142 boundaryRadiation.emissivity(patch().index(), 0, nullptr, &Tp)
143 );
144
145 const scalarField& emissivity = temissivity();
146
147 const scalarField Ep(emissivity/(2.0*(scalar(2) - emissivity)));
148
149 // Set value fraction
150 valueFraction() = 1.0/(1.0 + gamma*patch().deltaCoeffs()/Ep);
152 UPstream::msgType(oldTag); // Restore tag
153
154 mixedFvPatchScalarField::updateCoeffs();
155}
156
157
159(
160 Ostream& os
161) const
162{
164 os.writeEntryIfDifferent<word>("T", "T", TName_);
165}
166
168// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
169
170namespace Foam
171{
172namespace radiation
173{
175 (
178 );
179}
180}
181
182// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
static FOAM_NO_DANGLING_REFERENCE const boundaryRadiationProperties & New(const fvMesh &mesh, Args &&... args)
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
static int & msgType() noexcept
Message tag of standard messages.
Definition UPstream.H:1926
static int incrMsgType(int val=1) noexcept
Increment the message tag for standard messages.
Definition UPstream.H:1948
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
A FieldMapper for finite-volume patch fields.
virtual void operator=(const UList< scalar > &)
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition fvPatch.H:71
virtual void write(Ostream &) const
Write.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
MarshakRadiationFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
A class for managing temporary objects.
Definition tmp.H:75
A class for handling words, derived from Foam::string.
Definition word.H:66
volScalarField & p
const scalar gamma
Definition EEqn.H:9
dynamicFvMesh & mesh
OBJstream os(runTime.globalPath()/outputName)
#define makePatchTypeField(PatchTypeField, typePatchTypeField)
Define a concrete fvPatchField type and add to run-time tables Example, (fvPatchScalarField,...
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m2/K4].
Namespace for radiation modelling.
Namespace for OpenFOAM.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensionedScalar pow4(const dimensionedScalar &ds)
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
fvPatchField< scalar > fvPatchScalarField
dictionary dict