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, false),
73  qro_("qro", dict, p.size())
74 {
75  if (dict.found("value"))
76  {
78  (
79  scalarField("value", dict, p.size())
80  );
81 
82  }
83  else
84  {
86  }
87 }
88 
89 
92 (
94 )
95 :
96  fixedValueFvPatchScalarField(ptf),
97  qro_(ptf.qro_)
98 {}
99 
100 
103 (
106 )
107 :
108  fixedValueFvPatchScalarField(ptf, iF),
109  qro_(ptf.qro_)
110 {}
111 
112 
113 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
114 
116 autoMap
117 (
118  const fvPatchFieldMapper& m
119 )
120 {
121  fixedValueFvPatchScalarField::autoMap(m);
122  qro_.autoMap(m);
123 }
124 
125 
127 (
128  const fvPatchScalarField& ptf,
129  const labelList& addr
130 )
131 {
132  fixedValueFvPatchScalarField::rmap(ptf, addr);
133 
135  refCast<const greyDiffusiveViewFactorFixedValueFvPatchScalarField>(ptf);
136 
137  qro_.rmap(mrptf.qro_, addr);
138 }
139 
140 
142 updateCoeffs()
143 {
144  if (this->updated())
145  {
146  return;
147  }
148 
149  if (debug)
150  {
151  scalar Q = gSum((*this)*patch().magSf());
152 
153  Info<< patch().boundaryMesh().mesh().name() << ':'
154  << patch().name() << ':'
155  << this->internalField().name() << " <- "
156  << " heat transfer rate:" << Q
157  << " wall radiative heat flux "
158  << " min:" << gMin(*this)
159  << " max:" << gMax(*this)
160  << " avg:" << gAverage(*this)
161  << endl;
162  }
163 }
164 
165 
168 {
169  tmp<scalarField> tqrt(new scalarField(qro_));
170 
171  const viewFactor& radiation =
172  db().lookupObject<viewFactor>("radiationProperties");
173 
174  if (radiation.useSolarLoad())
175  {
176  tqrt.ref() += patch().lookupPatchField<volScalarField, scalar>
177  (
178  radiation.primaryFluxName_ + "_" + name(bandI)
179  );
180 
181  word qSecName = radiation.relfectedFluxName_ + "_" + name(bandI);
182 
183  if (this->db().foundObject<volScalarField>(qSecName))
184  {
185  const volScalarField& qSec =
186  this->db().lookupObject<volScalarField>(qSecName);
187 
188  tqrt.ref() += qSec.boundaryField()[patch().index()];
189  }
190  }
191 
192  return tqrt;
193 }
194 
195 
197 write
198 (
199  Ostream& os
200 ) const
201 {
203  qro_.writeEntry("qro", os);
204 }
205 
206 
207 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
208 
209 namespace Foam
210 {
211 namespace radiation
212 {
214  (
217  );
218 }
219 }
220 
221 
222 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Foam::fvPatchField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: volSurfaceMapping.H:47
volFields.H
Foam::fvPatchField::operator=
virtual void operator=(const UList< Type > &)
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:46
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:63
Foam::radiation::makePatchTypeField
makePatchTypeField(fvPatchScalarField, greyDiffusiveRadiationMixedFvPatchScalarField)
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:57
Foam::gAverage
Type gAverage(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:597
Foam::radiation::greyDiffusiveViewFactorFixedValueFvPatchScalarField::autoMap
virtual void autoMap(const fvPatchFieldMapper &)
Definition: greyDiffusiveViewFactorFixedValueFvPatchScalarField.C:110
Foam::fvPatchField::operator
friend Ostream & operator(Ostream &, const fvPatchField< Type > &)
greyDiffusiveViewFactorFixedValueFvPatchScalarField.H
Foam::endl
Ostream & endl(Ostream &os)
Definition: Ostream.H:381
fvPatchFieldMapper.H
Foam::gSum
Type gSum(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:587
Foam::radiation::greyDiffusiveViewFactorFixedValueFvPatchScalarField::qro
tmp< scalarField > qro(label bandI=0) const
Definition: greyDiffusiveViewFactorFixedValueFvPatchScalarField.C:160
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:220
Foam::Info
messageStream Info
Foam::radiation::greyDiffusiveViewFactorFixedValueFvPatchScalarField::rmap
virtual void rmap(const fvPatchScalarField &, const labelList &)
Definition: greyDiffusiveViewFactorFixedValueFvPatchScalarField.C:120
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
Foam::radiation::viewFactor
View factor radiation model. The system solved is: C q = b where: Cij = deltaij/Ej - (1/Ej - 1)Fij q ...
Definition: viewFactor.H:69
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:53
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:119
os
OBJstream os(runTime.globalPath()/outputName)
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Definition: atmBoundaryLayer.C:26
Foam::radiation::greyDiffusiveViewFactorFixedValueFvPatchScalarField::updateCoeffs
virtual void updateCoeffs()
Definition: greyDiffusiveViewFactorFixedValueFvPatchScalarField.C:135
Foam::radiation::greyDiffusiveViewFactorFixedValueFvPatchScalarField
This boundary condition provides a grey-diffuse condition for radiative heat flux,...
Definition: greyDiffusiveViewFactorFixedValueFvPatchScalarField.H:84
Foam::foamVersion::patch
const std::string patch
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: BitOps.H:58
Foam::vtk::write
void write(vtk::formatter &fmt, const Type &val, const label n=1)
Definition: foamVtkOutputTemplates.C:29
Foam::name
word name(const expressions::valueTypeCode typeCode)
Definition: exprTraits.C:52
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:41
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:52
Foam::gMin
Type gMin(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:586
Foam::GeometricField
Generic GeometricField class.
Definition: areaFieldsFwd.H:49
Foam::radiation::greyDiffusiveViewFactorFixedValueFvPatchScalarField::greyDiffusiveViewFactorFixedValueFvPatchScalarField
greyDiffusiveViewFactorFixedValueFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Definition: greyDiffusiveViewFactorFixedValueFvPatchScalarField.C:33
Foam::radiation::greyDiffusiveViewFactorFixedValueFvPatchScalarField::write
virtual void write(Ostream &) const
Definition: greyDiffusiveViewFactorFixedValueFvPatchScalarField.C:191
radiationModel.H
Foam::gMax
Type gMax(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:585
Foam::GeometricField::boundaryField
const Boundary & boundaryField() const
Definition: GeometricFieldI.H:55
viewFactor.H
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:50