greyDiffusiveViewFactorFixedValueFvPatchScalarField.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2018 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 
31 
32 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
33 
36 (
37  const fvPatch& p,
39 )
40 :
41  fixedValueFvPatchScalarField(p, iF),
42  radiationCoupledBase(patch(), "undefined", scalarField::null()),
43  qro_(p.size(), 0.0)
44 {}
45 
46 
49 (
51  const fvPatch& p,
53  const fvPatchFieldMapper& mapper
54 )
55 :
56  fixedValueFvPatchScalarField(ptf, p, iF, mapper),
58  (
59  patch(),
60  ptf.emissivityMethod(),
61  ptf.emissivity_,
62  mapper
63  ),
64  qro_(ptf.qro_, mapper)
65 {}
66 
67 
70 (
71  const fvPatch& p,
73  const dictionary& dict
74 )
75 :
76  fixedValueFvPatchScalarField(p, iF, dict, false),
77  radiationCoupledBase(p, dict),
78  qro_("qro", dict, p.size())
79 {
80  if (dict.found("value"))
81  {
82  fvPatchScalarField::operator=
83  (
84  scalarField("value", dict, p.size())
85  );
86 
87  }
88  else
89  {
90  fvPatchScalarField::operator=(0.0);
91  }
92 }
93 
94 
97 (
99 )
100 :
101  fixedValueFvPatchScalarField(ptf),
103  (
104  ptf.patch(),
105  ptf.emissivityMethod(),
106  ptf.emissivity_
107  ),
108  qro_(ptf.qro_)
109 {}
110 
111 
114 (
117 )
118 :
119  fixedValueFvPatchScalarField(ptf, iF),
121  (
122  ptf.patch(),
123  ptf.emissivityMethod(),
124  ptf.emissivity_
125  ),
126  qro_(ptf.qro_)
127 {}
128 
129 
130 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
131 
133 autoMap
134 (
135  const fvPatchFieldMapper& m
136 )
137 {
138  fixedValueFvPatchScalarField::autoMap(m);
139  radiationCoupledBase::autoMap(m);
140  qro_.autoMap(m);
141 }
142 
143 
145 (
146  const fvPatchScalarField& ptf,
147  const labelList& addr
148 )
149 {
150  fixedValueFvPatchScalarField::rmap(ptf, addr);
151  radiationCoupledBase::rmap(ptf, addr);
153  refCast<const greyDiffusiveViewFactorFixedValueFvPatchScalarField>(ptf);
154 
155  qro_.rmap(mrptf.qro_, addr);
156 }
157 
160 {
161  if (this->updated())
162  {
163  return;
164  }
165 
166  if (debug)
167  {
168  scalar Q = gSum((*this)*patch().magSf());
169 
170  Info<< patch().boundaryMesh().mesh().name() << ':'
171  << patch().name() << ':'
172  << this->internalField().name() << " <- "
173  << " heat transfer rate:" << Q
174  << " wall radiative heat flux "
175  << " min:" << gMin(*this)
176  << " max:" << gMax(*this)
177  << " avg:" << gAverage(*this)
178  << endl;
179  }
180 
181  fixedValueFvPatchScalarField::updateCoeffs();
182 }
183 
184 
186 write
187 (
188  Ostream& os
189 ) const
190 {
193  qro_.writeEntry("qro", os);
194 }
195 
196 
197 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
198 
199 namespace Foam
200 {
201 namespace radiation
202 {
204  (
207  );
208 }
209 }
210 
211 
212 // ************************************************************************* //
makePatchTypeField(fvPatchScalarField, greyDiffusiveRadiationMixedFvPatchScalarField)
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:431
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
Type gMin(const FieldField< Field, Type > &f)
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
void writeEntry(const word &keyword, Ostream &os) const
Write the field as a dictionary entry.
Definition: Field.C:726
scalarField emissivity_
Emissivity.
Macros for easy insertion into run-time selection tables.
void write(Ostream &, const label, const dictionary &)
Write with dictionary lookup.
Type gSum(const FieldField< Field, Type > &f)
greyDiffusiveViewFactorFixedValueFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
void write(Ostream &) const
Write.
Foam::fvPatchFieldMapper.
Common functions to emissivity. It gets supplied from lookup into a dictionary or calculated by the s...
word emissivityMethod() const
Method to obtain emissivity.
virtual label size() const
Return size.
Definition: fvPatch.H:161
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
volScalarField scalarField(fieldObject, mesh)
Type gMax(const FieldField< Field, Type > &f)
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)
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
messageStream Info
virtual void rmap(const fvPatchScalarField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
autoPtr< radiation::radiationModel > radiation(radiation::radiationModel::New(T))
Namespace for OpenFOAM.