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 | Copyright (C) 2011-2015 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
28 #include "fvPatchFieldMapper.H"
29 #include "volFields.H"
30 #include "radiationModel.H"
31 
32 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
33 
36 (
37  const fvPatch& p,
39 )
40 :
41  fixedValueFvPatchScalarField(p, iF),
42  Qro_(),
43  solarLoad_(false)
44 {}
45 
46 
49 (
51  const fvPatch& p,
53  const fvPatchFieldMapper& mapper
54 )
55 :
56  fixedValueFvPatchScalarField(ptf, p, iF, mapper),
57  Qro_(ptf.Qro_, mapper),
58  solarLoad_(ptf.solarLoad_)
59 {}
60 
61 
64 (
65  const fvPatch& p,
67  const dictionary& dict
68 )
69 :
70  fixedValueFvPatchScalarField(p, iF),
71  Qro_("Qro", dict, p.size()),
72  solarLoad_(dict.lookupOrDefault<bool>("solarLoad", false))
73 {
74  if (dict.found("value"))
75  {
76  fvPatchScalarField::operator=
77  (
78  scalarField("value", dict, p.size())
79  );
80 
81  }
82  else
83  {
84  fvPatchScalarField::operator=(0.0);
85  }
86 }
87 
88 
91 (
93 )
94 :
95  fixedValueFvPatchScalarField(ptf),
96  Qro_(ptf.Qro_),
97  solarLoad_(ptf.solarLoad_)
98 {}
99 
100 
103 (
106 )
107 :
108  fixedValueFvPatchScalarField(ptf, iF),
109  Qro_(ptf.Qro_),
110  solarLoad_(ptf.solarLoad_)
111 {}
112 
113 
114 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
115 
118 {
119  if (this->updated())
120  {
121  return;
122  }
123 
124 
125  // Do nothing
126 
127  if (debug)
128  {
129  scalar Q = gSum((*this)*patch().magSf());
130 
131  Info<< patch().boundaryMesh().mesh().name() << ':'
132  << patch().name() << ':'
133  << this->dimensionedInternalField().name() << " <- "
134  << " heat transfer rate:" << Q
135  << " wall radiative heat flux "
136  << " min:" << gMin(*this)
137  << " max:" << gMax(*this)
138  << " avg:" << gAverage(*this)
139  << endl;
140  }
141 }
142 
143 
146 {
147  tmp<scalarField> tQrt(new scalarField(Qro_));
148 
149  if (solarLoad_)
150  {
151  const radiationModel& radiation =
152  db().lookupObject<radiationModel>("radiationProperties");
153 
154  tQrt() += patch().lookupPatchField<volScalarField,scalar>
155  (
156  radiation.externalRadHeatFieldName_
157  );
158  }
159 
160  return tQrt;
161 }
162 
163 
165 write
166 (
167  Ostream& os
168 ) const
169 {
171  Qro_.writeEntry("Qro", os);
172  os.writeKeyword("solarLoad") << solarLoad_ << token::END_STATEMENT << nl;
173 }
174 
175 
176 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
177 
178 namespace Foam
179 {
180 namespace radiation
181 {
183  (
186  );
187 }
188 }
189 
190 
191 // ************************************************************************* //
Foam::fvPatchField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:65
volFields.H
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:48
Foam::token::END_STATEMENT
@ END_STATEMENT
Definition: token.H:99
p
p
Definition: pEqn.H:62
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::gAverage
Type gAverage(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:571
dimensionedInternalField
rDeltaT dimensionedInternalField()
Foam::radiation::makePatchTypeField
makePatchTypeField(fvPatchScalarField, boundaryRadiationPropertiesFvPatchField)
greyDiffusiveViewFactorFixedValueFvPatchScalarField.H
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
fvPatchFieldMapper.H
Foam::gSum
Type gSum(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:564
Foam::Q
This function object calculates and outputs the second invariant of the velocity gradient tensor [1/s...
Definition: Q.H:119
Foam::radiation::greyDiffusiveViewFactorFixedValueFvPatchScalarField::Qro_
scalarField Qro_
External radiative heat flux.
Definition: greyDiffusiveViewFactorFixedValueFvPatchScalarField.H:93
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
Foam::radiation::greyDiffusiveViewFactorFixedValueFvPatchScalarField::Qro
tmp< scalarField > Qro() const
Return external + solar load radiative heat flux.
Definition: greyDiffusiveViewFactorFixedValueFvPatchScalarField.C:145
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
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 any number of values (e....
Definition: dictionary.H:137
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::radiation::greyDiffusiveViewFactorFixedValueFvPatchScalarField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: greyDiffusiveViewFactorFixedValueFvPatchScalarField.C:117
Foam::radiation::greyDiffusiveViewFactorFixedValueFvPatchScalarField::solarLoad_
bool solarLoad_
Activate solar load.
Definition: greyDiffusiveViewFactorFixedValueFvPatchScalarField.H:96
Foam::radiation::greyDiffusiveViewFactorFixedValueFvPatchScalarField
This boundary condition provides a grey-diffuse condition for radiative heat flux,...
Definition: greyDiffusiveViewFactorFixedValueFvPatchScalarField.H:86
Foam::radiation::radiationModel
Top level model for radiation modelling.
Definition: radiationModel.H:73
scalarField
volScalarField scalarField(fieldObject, mesh)
Foam::Ostream::writeKeyword
Ostream & writeKeyword(const keyType &)
Write the keyword followed by an appropriate indentation.
Definition: Ostream.C:59
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:45
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:53
Foam::gMin
Type gMin(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:563
Foam::radiation::greyDiffusiveViewFactorFixedValueFvPatchScalarField::greyDiffusiveViewFactorFixedValueFvPatchScalarField
greyDiffusiveViewFactorFixedValueFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Definition: greyDiffusiveViewFactorFixedValueFvPatchScalarField.C:36
Foam::radiation::greyDiffusiveViewFactorFixedValueFvPatchScalarField::write
virtual void write(Ostream &) const
Write.
Definition: greyDiffusiveViewFactorFixedValueFvPatchScalarField.C:166
radiationModel.H
write
Tcoeff write()
Foam::gMax
Type gMax(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:562
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:51