Loading...
Searching...
No Matches
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-------------------------------------------------------------------------------
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"
35
36using namespace Foam::constant;
37
38// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
39
42(
43 const fvPatch& p,
45)
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)
62 fixedGradientFvPatchScalarField(psf, p, iF, mapper),
63 temperatureCoupledBase(patch(), psf),
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))
93)
95 fixedGradientFvPatchScalarField(psf, iF),
96 temperatureCoupledBase(patch(), psf),
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)
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
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 {
174
175 scalar qr = gSum(kappa(*this)*gradient()*patch().magSf());
176 Info<< patch().boundaryMesh().mesh().name() << ':'
177 << patch().name() << ':'
178 << this->internalField().name() << " -> "
179 << " radiativeFlux:" << qr
180 << " walltemperature "
181 << " min:" << limits.min()
182 << " max:" << limits.max()
183 << " avg:" << gAverage(*this)
184 << endl;
185 }
186}
187
188
190(
191 Ostream& os
192) const
193{
196 qrIncident_.writeEntry("qrIncident", os);
198}
199
201// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
202
203namespace Foam
204{
205namespace radiation
206{
208 (
211 );
212}
213}
214
215// ************************************************************************* //
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...
A min/max value pair with additional methods. In addition to conveniently storing values,...
Definition MinMax.H:106
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
virtual void write(Ostream &) const
Write.
A FieldMapper for finite-volume patch fields.
void writeValueEntry(Ostream &os) const
Write *this field as a "value" entry.
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.
fixedIncidentRadiationFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
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.
void write(Ostream &os) const
Write.
virtual void rmap(const fvPatchField< scalar > &, const labelList &)=0
Reverse map the given fvPatchField onto this fvPatchField.
virtual tmp< scalarField > kappa(const scalarField &Tp) const
Given patch temperature calculate corresponding K field.
temperatureCoupledBase(const fvPatch &patch, const KMethodType method=KMethodType::mtFluidThermo)
Default construct from patch, using fluidThermo (default) or specified method.
volScalarField & p
auto limits
Definition setRDeltaT.H:186
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].
Different types of constants.
Namespace for handling debugging switches.
Definition debug.C:45
Namespace for radiation modelling.
Namespace for OpenFOAM.
Type gAverage(const FieldField< Field, Type > &f, const label comm)
The global arithmetic average of a FieldField.
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
Type gSum(const FieldField< Field, Type > &f)
List< label > labelList
A List of labels.
Definition List.H:62
messageStream Info
Information stream (stdout output on master, null elsewhere).
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
dimensionedScalar pow4(const dimensionedScalar &ds)
MinMax< Type > gMinMax(const FieldField< Field, Type > &f)
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
fvPatchField< scalar > fvPatchScalarField
dictionary dict