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-2023 OpenCFD Ltd.
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 
30 #include "fvPatchFieldMapper.H"
31 #include "volFields.H"
32 #include "constants.H"
33 #include "radiationModel.H"
35 
36 using namespace Foam::constant;
37 
38 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
39 
42 (
43  const fvPatch& p,
45 )
46 :
47  fixedGradientFvPatchScalarField(p, iF),
48  temperatureCoupledBase(patch()), // default method (fluidThermo)
49  qrIncident_(p.size(), Zero)
50 {}
51 
52 
55 (
57  const fvPatch& p,
59  const fvPatchFieldMapper& mapper
60 )
61 :
62  fixedGradientFvPatchScalarField(psf, p, iF, mapper),
64  qrIncident_(psf.qrIncident_)
65 {}
66 
67 
70 (
71  const fvPatch& p,
73  const dictionary& dict
74 )
75 :
76  fixedGradientFvPatchScalarField(p, iF), // Bypass dictionary constructor
78  qrIncident_("qrIncident", dict, p.size())
79 {
80  if (!this->readGradientEntry(dict) || !this->readValueEntry(dict))
81  {
82  extrapolateInternal();
83  gradient() = Zero;
84  }
85 }
86 
87 
90 (
93 )
94 :
95  fixedGradientFvPatchScalarField(psf, iF),
97  qrIncident_(psf.qrIncident_)
98 {}
99 
100 
103 (
105 )
106 :
107  fixedGradientFvPatchScalarField(ptf),
109  qrIncident_(ptf.qrIncident_)
110 {}
111 
112 
113 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
114 
116 (
117  const fvPatchFieldMapper& m
118 )
119 {
120  fixedGradientFvPatchScalarField::autoMap(m);
122  qrIncident_.autoMap(m);
123 }
124 
125 
127 (
128  const fvPatchScalarField& psf,
129  const labelList& addr
130 )
131 {
132  fixedGradientFvPatchScalarField::rmap(psf, addr);
133 
135  refCast<const fixedIncidentRadiationFvPatchScalarField>
136  (
137  psf
138  );
139 
140  temperatureCoupledBase::rmap(thftpsf, addr);
141  qrIncident_.rmap(thftpsf.qrIncident_, addr);
142 }
143 
144 
146 {
147  if (updated())
148  {
149  return;
150  }
151 
152  scalarField intFld(patchInternalField());
153 
155  db().lookupObject<radiation::radiationModel>("radiationProperties");
156 
157  scalarField emissivity
158  (
159  radiation.absorptionEmission().e()().boundaryField()[patch().index()]
160  );
161 
162  gradient() =
163  emissivity
164  *(
165  qrIncident_
166  - physicoChemical::sigma.value()*pow4(*this)
167  )/kappa(*this);
168 
169  fixedGradientFvPatchScalarField::updateCoeffs();
170 
171  if (debug)
172  {
173  scalar qr = gSum(kappa(*this)*gradient()*patch().magSf());
174  Info<< patch().boundaryMesh().mesh().name() << ':'
175  << patch().name() << ':'
176  << this->internalField().name() << " -> "
177  << " radiativeFlux:" << qr
178  << " walltemperature "
179  << " min:" << gMin(*this)
180  << " max:" << gMax(*this)
181  << " avg:" << gAverage(*this)
182  << endl;
183  }
184 }
185 
186 
188 (
189  Ostream& os
190 ) const
191 {
194  qrIncident_.writeEntry("qrIncident", os);
196 }
197 
198 
199 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
200 
201 namespace Foam
202 {
203 namespace radiation
204 {
206  (
208  fixedIncidentRadiationFvPatchScalarField
209  );
210 }
211 }
212 
213 // ************************************************************************* //
Different types of constants.
makePatchTypeField(fvPatchScalarField, greyDiffusiveRadiationMixedFvPatchScalarField)
dictionary dict
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
Type gMin(const FieldField< Field, Type > &f)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:70
virtual void autoMap(const fvPatchFieldMapper &)
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.
void writeValueEntry(Ostream &os) const
Write *this field as a "value" entry.
Definition: fvPatchField.H:375
virtual void autoMap(const fvPatchFieldMapper &)=0
Map (and resize as needed) from self given a mapping object.
void write(Ostream &os) const
Write.
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
dimensionedScalar sigma("sigma", dimMass/sqr(dimTime), transportProperties)
virtual void rmap(const fvPatchScalarField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
Macros for easy insertion into run-time selection tables.
fixedIncidentRadiationFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Boundary condition for thermal coupling for solid regions. Used to emulate a fixed incident radiative...
virtual void write(Ostream &) const
Write.
Type gSum(const FieldField< Field, Type > &f)
fvPatchField< scalar > fvPatchScalarField
A FieldMapper for finite-volume patch fields.
Top level model for radiation modelling.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
int debug
Static debugging option.
Type gMax(const FieldField< Field, Type > &f)
OBJstream os(runTime.globalPath()/outputName)
Common functions used in temperature coupled boundaries.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Type gAverage(const FieldField< Field, Type > &f)
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: areaFieldsFwd.H:42
dimensionedScalar pow4(const dimensionedScalar &ds)
const std::string patch
OpenFOAM patch number as a std::string.
messageStream Info
Information stream (stdout output on master, null elsewhere)
volScalarField & p
autoPtr< radiation::radiationModel > radiation(radiation::radiationModel::New(T))
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127