fixedIncidentRadiationFvPatchScalarField.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) 2016-2019,2021 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
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
30#include "fvPatchFieldMapper.H"
31#include "volFields.H"
32#include "constants.H"
33#include "radiationModel.H"
35
36using namespace Foam::constant;
37
38// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
39
40Foam::radiation::fixedIncidentRadiationFvPatchScalarField::
41fixedIncidentRadiationFvPatchScalarField
42(
43 const fvPatch& p,
45)
46:
47 fixedGradientFvPatchScalarField(p, iF),
49 (
50 patch(),
51 "undefined",
52 "undefined",
53 "undefined-K",
54 "undefined-alpha"
55 ),
56 qrIncident_(p.size(), Zero)
57{}
58
59
60Foam::radiation::fixedIncidentRadiationFvPatchScalarField::
61fixedIncidentRadiationFvPatchScalarField
62(
64 const fvPatch& p,
66 const fvPatchFieldMapper& mapper
67)
68:
69 fixedGradientFvPatchScalarField(psf, p, iF, mapper),
70 temperatureCoupledBase(patch(), psf),
71 qrIncident_(psf.qrIncident_)
72{}
73
74
75Foam::radiation::fixedIncidentRadiationFvPatchScalarField::
76fixedIncidentRadiationFvPatchScalarField
77(
78 const fvPatch& p,
80 const dictionary& dict
81)
82:
83 fixedGradientFvPatchScalarField(p, iF),
85 qrIncident_("qrIncident", dict, p.size())
86{
87 if (dict.found("value") && dict.found("gradient"))
88 {
90 gradient() = Field<scalar>("gradient", dict, p.size());
91 }
92 else
93 {
94 // Still reading so cannot yet evaluate. Make up a value.
95 fvPatchField<scalar>::operator=(patchInternalField());
96 gradient() = 0.0;
97 }
98}
99
100
101Foam::radiation::fixedIncidentRadiationFvPatchScalarField::
102fixedIncidentRadiationFvPatchScalarField
103(
106)
107:
108 fixedGradientFvPatchScalarField(psf, iF),
109 temperatureCoupledBase(patch(), psf),
110 qrIncident_(psf.qrIncident_)
111{}
112
113
114Foam::radiation::fixedIncidentRadiationFvPatchScalarField::
115fixedIncidentRadiationFvPatchScalarField
116(
118)
119:
120 fixedGradientFvPatchScalarField(ptf),
121 temperatureCoupledBase(patch(), ptf),
122 qrIncident_(ptf.qrIncident_)
123{}
124
125
126// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
127
129(
130 const fvPatchFieldMapper& m
131)
132{
133 fixedGradientFvPatchScalarField::autoMap(m);
135 qrIncident_.autoMap(m);
136}
137
138
140(
141 const fvPatchScalarField& psf,
142 const labelList& addr
143)
144{
145 fixedGradientFvPatchScalarField::rmap(psf, addr);
146
148 refCast<const fixedIncidentRadiationFvPatchScalarField>
149 (
150 psf
151 );
152
153 temperatureCoupledBase::rmap(thftpsf, addr);
154 qrIncident_.rmap(thftpsf.qrIncident_, addr);
155}
156
157
159{
160 if (updated())
161 {
162 return;
163 }
164
165 scalarField intFld(patchInternalField());
166
168 db().lookupObject<radiation::radiationModel>("radiationProperties");
169
170 scalarField emissivity
171 (
172 radiation.absorptionEmission().e()().boundaryField()[patch().index()]
173 );
174
175 gradient() =
176 emissivity
177 *(
178 qrIncident_
180 )/kappa(*this);
181
182 fixedGradientFvPatchScalarField::updateCoeffs();
183
184 if (debug)
185 {
186 scalar qr = gSum(kappa(*this)*gradient()*patch().magSf());
187 Info<< patch().boundaryMesh().mesh().name() << ':'
188 << patch().name() << ':'
189 << this->internalField().name() << " -> "
190 << " radiativeFlux:" << qr
191 << " walltemperature "
192 << " min:" << gMin(*this)
193 << " max:" << gMax(*this)
194 << " avg:" << gAverage(*this)
195 << endl;
196 }
197}
198
199
201(
202 Ostream& os
203) const
204{
205 fixedGradientFvPatchScalarField::write(os);
207 qrIncident_.writeEntry("qrIncident", os);
208 writeEntry("value", os);
209}
210
211
212// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
213
214namespace Foam
215{
216namespace radiation
217{
219 (
222 );
223}
224}
225
226// ************************************************************************* //
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...
Generic templated field type.
Definition: Field.H:82
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
bool found(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Search for an entry (const access) with the given keyword.
Definition: dictionaryI.H:87
const Type & value() const
Return const reference to value.
virtual bool write()
Write the output fields.
A FieldMapper for finite-volume patch fields.
virtual void operator=(const UList< Type > &)
Definition: fvPatchField.C:408
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:71
Boundary condition for thermal coupling for solid regions. Used to emulate a fixed incident radiative...
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.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Top level model for radiation modelling.
Common functions used in temperature coupled boundaries.
virtual void autoMap(const fvPatchFieldMapper &)=0
Map (and resize as needed) from self given a mapping object.
virtual void rmap(const fvPatchField< scalar > &, const labelList &)=0
Reverse map the given fvPatchField onto this fvPatchField.
volScalarField & p
autoPtr< radiation::radiationModel > radiation(radiation::radiationModel::New(T))
OBJstream os(runTime.globalPath()/outputName)
#define makePatchTypeField(PatchTypeField, typePatchTypeField)
Definition: fvPatchField.H:676
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m2/K4].
Different types of constants.
Namespace for OpenFOAM.
Type gSum(const FieldField< Field, Type > &f)
messageStream Info
Information stream (stdout output on master, null elsewhere)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
dimensionedScalar pow4(const dimensionedScalar &ds)
Type gAverage(const FieldField< Field, Type > &f)
Type gMin(const FieldField< Field, Type > &f)
Type gMax(const FieldField< Field, Type > &f)
dictionary dict