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 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),
49  (
50  patch(),
51  "undefined",
52  "undefined",
53  "undefined-K",
54  "undefined-alpha"
55  ),
56  qrIncident_(p.size(), Zero)
57 {}
58 
59 
62 (
64  const fvPatch& p,
66  const fvPatchFieldMapper& mapper
67 )
68 :
69  fixedGradientFvPatchScalarField(psf, p, iF, mapper),
71  qrIncident_(psf.qrIncident_)
72 {}
73 
74 
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 
103 (
106 )
107 :
108  fixedGradientFvPatchScalarField(psf, iF),
110  qrIncident_(psf.qrIncident_)
111 {}
112 
113 
116 (
118 )
119 :
120  fixedGradientFvPatchScalarField(ptf),
122  qrIncident_(ptf.qrIncident_)
123 {}
124 
125 
126 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
127 
129 (
130  const fvPatchFieldMapper& m
131 )
132 {
133  fixedGradientFvPatchScalarField::autoMap(m);
134  qrIncident_.autoMap(m);
135 }
136 
137 
139 (
140  const fvPatchScalarField& psf,
141  const labelList& addr
142 )
143 {
144  fixedGradientFvPatchScalarField::rmap(psf, addr);
145 
147  refCast<const fixedIncidentRadiationFvPatchScalarField>
148  (
149  psf
150  );
151 
152  qrIncident_.rmap(thftpsf.qrIncident_, addr);
153 }
154 
155 
157 {
158  if (updated())
159  {
160  return;
161  }
162 
163  scalarField intFld(patchInternalField());
164 
166  db().lookupObject<radiation::radiationModel>("radiationProperties");
167 
168  scalarField emissivity
169  (
170  radiation.absorptionEmission().e()().boundaryField()[patch().index()]
171  );
172 
173  gradient() =
174  emissivity
175  *(
176  qrIncident_
178  )/kappa(*this);
179 
180  fixedGradientFvPatchScalarField::updateCoeffs();
181 
182  if (debug)
183  {
184  scalar qr = gSum(kappa(*this)*gradient()*patch().magSf());
185  Info<< patch().boundaryMesh().mesh().name() << ':'
186  << patch().name() << ':'
187  << this->internalField().name() << " -> "
188  << " radiativeFlux:" << qr
189  << " walltemperature "
190  << " min:" << gMin(*this)
191  << " max:" << gMax(*this)
192  << " avg:" << gAverage(*this)
193  << endl;
194  }
195 }
196 
197 
199 (
200  Ostream& os
201 ) const
202 {
205  qrIncident_.writeEntry("qrIncident", os);
206  writeEntry("value", os);
207 }
208 
209 
210 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
211 
212 namespace Foam
213 {
214 namespace radiation
215 {
217  (
220  );
221 }
222 }
223 
224 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::fvPatchField< scalar >
volFields.H
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::radiation::fixedIncidentRadiationFvPatchScalarField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: fixedIncidentRadiationFvPatchScalarField.C:156
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::gAverage
Type gAverage(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:604
Foam::constant
Different types of constants.
Definition: atomicConstants.C:38
Foam::temperatureCoupledBase
Common functions used in temperature coupled boundaries.
Definition: temperatureCoupledBase.H:110
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:350
fvPatchFieldMapper.H
Foam::dimensioned::value
const Type & value() const
Return const reference to value.
Definition: dimensionedType.C:434
Foam::gSum
Type gSum(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:594
Foam::radiation::fixedIncidentRadiationFvPatchScalarField::write
virtual void write(Ostream &) const
Write.
Definition: fixedIncidentRadiationFvPatchScalarField.C:199
Foam::constant::electromagnetic::kappa
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
Foam::pow4
dimensionedScalar pow4(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:100
Foam::Field< scalar >
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:65
Foam::radiation::fixedIncidentRadiationFvPatchScalarField::fixedIncidentRadiationFvPatchScalarField
fixedIncidentRadiationFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Definition: fixedIncidentRadiationFvPatchScalarField.C:42
radiation
autoPtr< radiation::radiationModel > radiation(radiation::radiationModel::New(T))
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
Foam::radiation::fixedIncidentRadiationFvPatchScalarField::autoMap
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
Definition: fixedIncidentRadiationFvPatchScalarField.C:129
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::radiation::fixedIncidentRadiationFvPatchScalarField
Boundary condition for thermal coupling for solid regions. Used to emulate a fixed incident radiative...
Definition: fixedIncidentRadiationFvPatchScalarField.H:84
constants.H
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
absorptionEmissionModel.H
Foam::List< label >
Foam::radiation::radiationModel
Top level model for radiation modelling.
Definition: radiationModel.H:75
Foam::constant::physicoChemical::sigma
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m2/K4].
Foam::vtk::write
void write(vtk::formatter &fmt, const Type &val, const label n=1)
Component-wise write of a value (N times)
Definition: foamVtkOutputTemplates.C:35
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:47
Foam::radiation::fixedIncidentRadiationFvPatchScalarField::rmap
virtual void rmap(const fvPatchScalarField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
Definition: fixedIncidentRadiationFvPatchScalarField.C:139
Foam::temperatureCoupledBase::write
void write(Ostream &os) const
Write.
Definition: temperatureCoupledBase.C:426
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::makePatchTypeField
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Foam::gMin
Type gMin(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:593
fixedIncidentRadiationFvPatchScalarField.H
Foam::fvPatchField::operator=
virtual void operator=(const UList< Type > &)
Definition: fvPatchField.C:379
radiationModel.H
Foam::gMax
Type gMax(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:592
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54