Loading...
Searching...
No Matches
MarshakRadiationFixedTemperatureFvPatchScalarField.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 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 Trad_(p.size())
61 const fvPatchFieldMapper& mapper
63:
64 mixedFvPatchScalarField(ptf, p, iF, mapper),
65 Trad_(ptf.Trad_, mapper)
66{}
67
68
71(
72 const fvPatch& p,
74 const dictionary& dict
75)
76:
77 mixedFvPatchScalarField(p, iF),
78 Trad_("Trad", dict, p.size())
79{
80 // refValue updated on each call to updateCoeffs()
81 refValue() = 4.0*constant::physicoChemical::sigma.value()*pow4(Trad_);
91(
94:
95 mixedFvPatchScalarField(ptf),
96 Trad_(ptf.Trad_)
97{}
98
99
102(
105)
106:
107 mixedFvPatchScalarField(ptf, iF),
108 Trad_(ptf.Trad_)
109{}
110
111
112// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
113
116(
117 const fvPatchFieldMapper& m
119{
120 mixedFvPatchScalarField::autoMap(m);
121 Trad_.autoMap(m);
122}
123
124
126(
127 const fvPatchScalarField& ptf,
128 const labelList& addr
129)
130{
131 mixedFvPatchScalarField::rmap(ptf, addr);
132
135
136 Trad_.rmap(mrptf.Trad_, addr);
137}
138
139
142{
143 if (this->updated())
144 {
145 return;
146 }
147
148 // Since we're inside initEvaluate/evaluate there might be processor
149 // comms underway. Change the tag we use.
150 const int oldTag = UPstream::incrMsgType();
151
152 // Re-calc reference value
153 refValue() = 4.0*constant::physicoChemical::sigma.value()*pow4(Trad_);
154
155 // Diffusion coefficient - created by radiation model's ::updateCoeffs()
156 const scalarField& gamma =
157 patch().lookupPatchField<volScalarField>("gammaRad");
158
159 //const scalarField temissivity = emissivity();
160 const boundaryRadiationProperties& boundaryRadiation =
161 boundaryRadiationProperties::New(internalField().mesh());
162
163 const tmp<scalarField> temissivity
164 (
165 boundaryRadiation.emissivity(patch().index())
166 );
167
168 const scalarField& emissivity = temissivity();
169
170 const scalarField Ep(emissivity/(2.0*(scalar(2) - emissivity)));
171
172 // Set value fraction
173 valueFraction() = 1.0/(1.0 + gamma*patch().deltaCoeffs()/Ep);
175 UPstream::msgType(oldTag); // Restore tag
176
177 mixedFvPatchScalarField::updateCoeffs();
178}
179
180
182(
183 Ostream& os
184) const
185{
187 Trad_.writeEntry("Trad", os);
188}
189
191// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
192
193namespace Foam
194{
195namespace radiation
196{
198 (
201 );
202}
203}
204
205// ************************************************************************* //
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 autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
virtual void rmap(const fvPatchScalarField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
MarshakRadiationFixedTemperatureFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
A class for managing temporary objects.
Definition tmp.H:75
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.
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
List< label > labelList
A List of labels.
Definition List.H:62
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