Loading...
Searching...
No Matches
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-2025 OpenCFD Ltd
10-------------------------------------------------------------------------------
11License
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
32#include "volFields.H"
33#include "radiationModel.H"
34#include "viewFactor.H"
35
36// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
37
40(
41 const fvPatch& p,
56 const fvPatchFieldMapper& mapper
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, IOobjectOption::NO_READ),
73 qro_("qro", dict, p.size())
74{
84(
87:
88 fixedValueFvPatchScalarField(ptf),
89 qro_(ptf.qro_)
90{}
91
92
95(
98)
99:
100 fixedValueFvPatchScalarField(ptf, iF),
101 qro_(ptf.qro_)
102{}
103
104
105// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
106
109(
110 const fvPatchFieldMapper& m
112{
113 fixedValueFvPatchScalarField::autoMap(m);
114 qro_.autoMap(m);
115}
116
117
119(
120 const fvPatchScalarField& ptf,
121 const labelList& addr
122)
123{
124 fixedValueFvPatchScalarField::rmap(ptf, addr);
125
128
129 qro_.rmap(mrptf.qro_, addr);
130}
131
132
135{
136 if (this->updated())
137 {
138 return;
139 }
140
141 if (debug)
142 {
143 scalar Q = gWeightedSum(patch().magSf(), *this);
144
146
147 Info<< patch().boundaryMesh().mesh().name() << ':'
148 << patch().name() << ':'
149 << this->internalField().name() << " <- "
150 << " heat transfer rate:" << Q
151 << " wall radiative heat flux "
152 << " min:" << limits.min()
153 << " max:" << limits.max()
154 << " avg:" << gAverage(*this)
155 << endl;
156 }
157}
158
159
162{
163 auto tqrt = tmp<scalarField>::New(qro_);
164
165 const viewFactor& radiation =
166 db().lookupObject<viewFactor>("radiationProperties");
167
168 if (radiation.useSolarLoad())
169 {
170 tqrt.ref() += patch().lookupPatchField<volScalarField>
171 (
172 radiation.primaryFluxName_ + "_" + name(bandI)
173 );
174
175 if
176 (
177 const auto* qSec
178 = patch().cfindPatchField<volScalarField>
179 (
180 radiation.relfectedFluxName_ + "_" + name(bandI)
181 )
182 )
183 {
184 tqrt.ref() += *qSec;
186 }
187
188 return tqrt;
189}
190
191
193write
194(
195 Ostream& os
196) const
197{
199 qro_.writeEntry("qro", os);
200}
201
203// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
204
205namespace Foam
206{
207namespace radiation
208{
210 (
213 );
214}
215}
216
217
218// ************************************************************************* //
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...
A simple container of IOobject preferences. Can also be used for general handling of read/no-read/rea...
A min/max value pair with additional methods. In addition to conveniently storing values,...
Definition MinMax.H:106
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
virtual void write(Ostream &) const
Write includes "value" entry.
A FieldMapper for finite-volume patch fields.
virtual void operator=(const UList< scalar > &)
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition fvPatch.H:71
This boundary condition provides a grey-diffuse condition for radiative heat flux,...
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
virtual void rmap(const fvPatchScalarField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
greyDiffusiveViewFactorFixedValueFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
tmp< scalarField > qro(label bandI=0) const
Return external + solar load radiative heat flux.
View factor radiation model. The system solved is: C q = b where: Cij = deltaij/Ej - (1/Ej - 1)Fij q ...
Definition viewFactor.H:76
A class for managing temporary objects.
Definition tmp.H:75
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition tmp.H:215
volScalarField & p
auto limits
Definition setRDeltaT.H:186
OBJstream os(runTime.globalPath()/outputName)
#define makePatchTypeField(PatchTypeField, typePatchTypeField)
Define a concrete fvPatchField type and add to run-time tables Example, (fvPatchScalarField,...
auto & name
Namespace for handling debugging switches.
Definition debug.C:45
Namespace for radiation modelling.
Namespace for OpenFOAM.
Type gAverage(const FieldField< Field, Type > &f, const label comm)
The global arithmetic average of a FieldField.
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
List< label > labelList
A List of labels.
Definition List.H:62
GeometricField< scalar, fvPatchField, volMesh > volScalarField
messageStream Info
Information stream (stdout output on master, null elsewhere).
Type gWeightedSum(const UList< scalar > &weights, const UList< Type > &fld, const label comm)
The global weighted sum (integral) of a field, using the mag() of the weights.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
MinMax< Type > gMinMax(const FieldField< Field, Type > &f)
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
fvPatchField< scalar > fvPatchScalarField
dictionary dict