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-2022 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 (
50  const fvPatch& p,
52  const dictionary& dict
53 )
54 :
55  fixedValueFvPatchScalarField(p, iF, dict, false),
56  radiationCoupledBase(p, dict),
57  qro_("qro", dict, p.size())
58 {
59  if (dict.found("value"))
60  {
61  fvPatchScalarField::operator=
62  (
63  scalarField("value", dict, p.size())
64  );
65 
66  }
67  else
68  {
69  fvPatchScalarField::operator=(0.0);
70  }
71 }
72 
73 
76 (
78  const fvPatch& p,
80  const fvPatchFieldMapper& mapper
81 )
82 :
83  fixedValueFvPatchScalarField(ptf, p, iF, mapper),
85  (
86  patch(),
87  ptf.emissivityMethod(),
88  ptf.emissivity_,
89  mapper
90  ),
91  qro_(mapper(ptf.qro_))
92 {}
93 
94 
97 (
100 )
101 :
102  fixedValueFvPatchScalarField(ptf, iF),
104  (
105  ptf.patch(),
106  ptf.emissivityMethod(),
107  ptf.emissivity_
108  ),
109  qro_(ptf.qro_)
110 {}
111 
112 
113 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
114 
116 (
117  const fvPatchFieldMapper& m
118 )
119 {
120  fixedValueFvPatchScalarField::autoMap(m);
121  radiationCoupledBase::autoMap(m);
122  m(qro_, qro_);
123 }
124 
125 
127 (
128  const fvPatchScalarField& ptf,
129  const labelList& addr
130 )
131 {
132  fixedValueFvPatchScalarField::rmap(ptf, addr);
133  radiationCoupledBase::rmap(ptf, addr);
135  refCast<const greyDiffusiveViewFactorFixedValueFvPatchScalarField>(ptf);
136 
137  qro_.rmap(mrptf.qro_, addr);
138 }
139 
140 
142 (
143  const fvPatchScalarField& ptf
144 )
145 {
146  fixedValueFvPatchScalarField::reset(ptf);
147  radiationCoupledBase::reset(ptf);
149  refCast<const greyDiffusiveViewFactorFixedValueFvPatchScalarField>(ptf);
150 
151  qro_.reset(mrptf.qro_);
152 }
153 
154 
156 {
157  if (this->updated())
158  {
159  return;
160  }
161 
162  if (debug)
163  {
164  scalar Q = gSum((*this)*patch().magSf());
165 
166  Info<< patch().boundaryMesh().mesh().name() << ':'
167  << patch().name() << ':'
168  << this->internalField().name() << " <- "
169  << " heat transfer rate:" << Q
170  << " wall radiative heat flux "
171  << " min:" << gMin(*this)
172  << " max:" << gMax(*this)
173  << " avg:" << gAverage(*this)
174  << endl;
175  }
176 
177  fixedValueFvPatchScalarField::updateCoeffs();
178 }
179 
180 
182 (
183  Ostream& os
184 ) const
185 {
188  writeEntry(os, "qro", qro_);
189 }
190 
191 
192 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
193 
194 namespace Foam
195 {
197  (
200  );
201 }
202 
203 
204 // ************************************************************************* //
This boundary condition provides a grey-diffuse condition for radiative heat flux, qr, for use with the view factor model.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:663
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
Type gMin(const FieldField< Field, Type > &f)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:63
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
scalarField emissivity_
Emissivity.
Macros for easy insertion into run-time selection tables.
Type gSum(const FieldField< Field, Type > &f)
greyDiffusiveViewFactorFixedValueFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
virtual void reset(const fvPatchScalarField &)
Reset the fvPatchField to the given fvPatchField.
void write(Ostream &) const
Write.
Foam::fvPatchFieldMapper.
void write(std::ostream &os, const bool binary, List< floatScalar > &fField)
Write floats ascii or binary.
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:157
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
volScalarField scalarField(fieldObject, mesh)
Type gMax(const FieldField< Field, Type > &f)
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
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
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
virtual void rmap(const fvPatchScalarField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
Namespace for OpenFOAM.