greyDiffusiveViewFactorFixedValueFvPatchScalarField.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-2017 OpenFOAM Foundation
9  Copyright (C) 2016 OpenCFD Ltd
10 -------------------------------------------------------------------------------
11 License
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"
34 #include "viewFactor.H"
35 
36 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
37 
40 (
41  const fvPatch& p,
43 )
44 :
45  fixedValueFvPatchScalarField(p, iF),
46  qro_()
47 {}
48 
49 
52 (
54  const fvPatch& p,
56  const fvPatchFieldMapper& mapper
57 )
58 :
59  fixedValueFvPatchScalarField(ptf, p, iF, mapper),
60  qro_(ptf.qro_, mapper)
61 {}
62 
63 
66 (
67  const fvPatch& p,
69  const dictionary& dict
70 )
71 :
72  fixedValueFvPatchScalarField(p, iF, dict, IOobjectOption::NO_READ),
73  qro_("qro", dict, p.size())
74 {
75  if (!this->readValueEntry(dict))
76  {
78  }
79 }
80 
81 
84 (
86 )
87 :
88  fixedValueFvPatchScalarField(ptf),
89  qro_(ptf.qro_)
90 {}
91 
92 
95 (
98 )
99 :
100  fixedValueFvPatchScalarField(ptf, iF),
101  qro_(ptf.qro_)
102 {}
103 
104 
105 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
106 
108 autoMap
109 (
110  const fvPatchFieldMapper& m
111 )
112 {
113  fixedValueFvPatchScalarField::autoMap(m);
114  qro_.autoMap(m);
115 }
116 
117 
119 (
120  const fvPatchScalarField& ptf,
121  const labelList& addr
122 )
123 {
124  fixedValueFvPatchScalarField::rmap(ptf, addr);
125 
127  refCast<const greyDiffusiveViewFactorFixedValueFvPatchScalarField>(ptf);
128 
129  qro_.rmap(mrptf.qro_, addr);
130 }
131 
132 
134 updateCoeffs()
135 {
136  if (this->updated())
137  {
138  return;
139  }
140 
141  if (debug)
142  {
143  scalar Q = gSum((*this)*patch().magSf());
144 
145  Info<< patch().boundaryMesh().mesh().name() << ':'
146  << patch().name() << ':'
147  << this->internalField().name() << " <- "
148  << " heat transfer rate:" << Q
149  << " wall radiative heat flux "
150  << " min:" << gMin(*this)
151  << " max:" << gMax(*this)
152  << " avg:" << gAverage(*this)
153  << endl;
154  }
155 }
156 
157 
160 {
161  tmp<scalarField> tqrt(new scalarField(qro_));
162 
163  const viewFactor& radiation =
164  db().lookupObject<viewFactor>("radiationProperties");
165 
166  if (radiation.useSolarLoad())
167  {
168  tqrt.ref() += patch().lookupPatchField<volScalarField>
169  (
170  radiation.primaryFluxName_ + "_" + name(bandI)
171  );
172 
173  word qSecName = radiation.relfectedFluxName_ + "_" + name(bandI);
174 
175  if (this->db().foundObject<volScalarField>(qSecName))
176  {
177  const volScalarField& qSec =
178  this->db().lookupObject<volScalarField>(qSecName);
179 
180  tqrt.ref() += qSec.boundaryField()[patch().index()];
181  }
182  }
184  return tqrt;
185 }
186 
187 
189 write
190 (
191  Ostream& os
192 ) const
193 {
195  qro_.writeEntry("qro", os);
196 }
197 
198 
199 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
200 
201 namespace Foam
202 {
203 namespace radiation
204 {
206  (
209  );
210 }
211 }
212 
213 
214 // ************************************************************************* //
makePatchTypeField(fvPatchScalarField, greyDiffusiveRadiationMixedFvPatchScalarField)
dictionary dict
virtual void write(Ostream &) const
Write.
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)
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition: tmpI.H:235
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
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
Macros for easy insertion into run-time selection tables.
Type gSum(const FieldField< Field, Type > &f)
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
greyDiffusiveViewFactorFixedValueFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
A class for handling words, derived from Foam::string.
Definition: word.H:63
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
A FieldMapper for finite-volume patch fields.
tmp< scalarField > qro(label bandI=0) const
Return external + solar load radiative heat flux.
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)
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
This boundary condition provides a grey-diffuse condition for radiative heat flux, qr, for use with the view factor model.
Type gAverage(const FieldField< Field, Type > &f)
virtual void operator=(const UList< scalar > &)
Definition: fvPatchField.C:391
A simple container of IOobject preferences. Can also be used for general handling of read/no-read/rea...
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: areaFieldsFwd.H:42
const std::string patch
OpenFOAM patch number as a std::string.
messageStream Info
Information stream (stdout output on master, null elsewhere)
volScalarField & p
A class for managing temporary objects.
Definition: HashPtrTable.H:50
virtual void rmap(const fvPatchScalarField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
View factor radiation model. The system solved is: C q = b where: Cij = deltaij/Ej - (1/Ej - 1)Fij q ...
Definition: viewFactor.H:73
autoPtr< radiation::radiationModel > radiation(radiation::radiationModel::New(T))
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127