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-2023 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  const dictionary& dict
40 )
41 :
42  fixedValueFvPatchScalarField(p, iF, dict, false),
44  qro_("qro", dict, p.size())
45 {
46  if (dict.found("value"))
47  {
49  (
50  scalarField("value", dict, p.size())
51  );
52 
53  }
54  else
55  {
57  }
58 }
59 
60 
63 (
65  const fvPatch& p,
67  const fvPatchFieldMapper& mapper
68 )
69 :
70  fixedValueFvPatchScalarField(ptf, p, iF, mapper),
72  (
73  patch(),
74  ptf.emissivityMethod(),
75  ptf.emissivity_,
76  mapper
77  ),
78  qro_(mapper(ptf.qro_))
79 {}
80 
81 
84 (
87 )
88 :
89  fixedValueFvPatchScalarField(ptf, iF),
91  (
92  ptf.patch(),
93  ptf.emissivityMethod(),
94  ptf.emissivity_
95  ),
96  qro_(ptf.qro_)
97 {}
98 
99 
100 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
101 
103 (
104  const fvPatchScalarField& ptf,
105  const fvPatchFieldMapper& mapper
106 )
107 {
108  fixedValueFvPatchScalarField::map(ptf, mapper);
109  radiationCoupledBase::map(ptf, mapper);
111  refCast<const greyDiffusiveViewFactorFixedValueFvPatchScalarField>(ptf);
112 
113  mapper(qro_, mrptf.qro_);
114 }
115 
116 
118 (
119  const fvPatchScalarField& ptf
120 )
121 {
122  fixedValueFvPatchScalarField::reset(ptf);
125  refCast<const greyDiffusiveViewFactorFixedValueFvPatchScalarField>(ptf);
126 
127  qro_.reset(mrptf.qro_);
128 }
129 
130 
132 {
133  if (this->updated())
134  {
135  return;
136  }
137 
138  if (debug)
139  {
140  scalar Q = gSum((*this)*patch().magSf());
141 
142  Info<< patch().boundaryMesh().mesh().name() << ':'
143  << patch().name() << ':'
144  << this->internalField().name() << " <- "
145  << " heat transfer rate:" << Q
146  << " wall radiative heat flux "
147  << " min:" << gMin(*this)
148  << " max:" << gMax(*this)
149  << " avg:" << gAverage(*this)
150  << endl;
151  }
152 
153  fixedValueFvPatchScalarField::updateCoeffs();
154 }
155 
156 
158 (
159  Ostream& os
160 ) const
161 {
164  writeEntry(os, "qro", qro_);
165 }
166 
167 
168 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
169 
170 namespace Foam
171 {
173  (
176  );
177 }
178 
179 
180 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
Foam::fvPatchFieldMapper.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:87
friend Ostream & operator(Ostream &, const fvPatchField< Type > &)
virtual void operator=(const UList< Type > &)
Definition: fvPatchField.C:251
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
This boundary condition provides a grey-diffuse condition for radiative heat flux,...
virtual void reset(const fvPatchScalarField &)
Reset the fvPatchField to the given fvPatchField.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
virtual void map(const fvPatchScalarField &, const fvPatchFieldMapper &)
Map the given fvPatchField onto this fvPatchField.
greyDiffusiveViewFactorFixedValueFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &)
Construct from patch, internal field and dictionary.
Common functions to emissivity. It gets supplied from lookup into a dictionary or calculated by the s...
void write(Ostream &) const
Write.
virtual void reset(const fvPatchScalarField &)
Reset the fvPatchField to the given fvPatchField.
virtual void map(const fvPatchScalarField &, const fvPatchFieldMapper &)
Map the given fvPatchField onto this fvPatchField.
void write(std::ostream &os, const bool binary, List< floatScalar > &fField)
Write floats ascii or binary.
Namespace for OpenFOAM.
Type gSum(const FieldField< Field, Type > &f)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
messageStream Info
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
Type gAverage(const FieldField< Field, Type > &f)
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Type gMin(const FieldField< Field, Type > &f)
Type gMax(const FieldField< Field, Type > &f)
dictionary dict
volScalarField & p