MarshakRadiationFvPatchScalarField.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-2015 OpenFOAM Foundation
9  Copyright (C) 2016,2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
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 
31 #include "fvPatchFieldMapper.H"
32 #include "volFields.H"
33 #include "radiationModel.H"
36 
37 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
38 
41 (
42  const fvPatch& p,
44 )
45 :
46  mixedFvPatchScalarField(p, iF),
47  TName_("T")
48 {
49  refValue() = 0.0;
50  refGrad() = 0.0;
51  valueFraction() = 0.0;
52 }
53 
54 
57 (
59  const fvPatch& p,
61  const fvPatchFieldMapper& mapper
62 )
63 :
64  mixedFvPatchScalarField(ptf, p, iF, mapper),
65  TName_(ptf.TName_)
66 {}
67 
68 
71 (
72  const fvPatch& p,
74  const dictionary& dict
75 )
76 :
77  mixedFvPatchScalarField(p, iF),
78  TName_(dict.getOrDefault<word>("T", "T"))
79 {
80  if (dict.found("value"))
81  {
82  refValue() = scalarField("value", dict, p.size());
83  }
84  else
85  {
86  refValue() = 0.0;
87  }
88 
89  // zero gradient
90  refGrad() = 0.0;
91 
92  valueFraction() = 1.0;
93 
95 }
96 
97 
100 (
102 )
103 :
104  mixedFvPatchScalarField(ptf),
105  TName_(ptf.TName_)
106 {}
107 
108 
111 (
114 )
115 :
116  mixedFvPatchScalarField(ptf, iF),
117  TName_(ptf.TName_)
118 {}
119 
120 
121 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
122 
124 {
125  if (this->updated())
126  {
127  return;
128  }
129 
130  // Since we're inside initEvaluate/evaluate there might be processor
131  // comms underway. Change the tag we use.
132  int oldTag = UPstream::msgType();
133  UPstream::msgType() = oldTag+1;
134 
135  // Temperature field
136  const scalarField& Tp =
137  patch().lookupPatchField<volScalarField, scalar>(TName_);
138 
139  // Re-calc reference value
140  refValue() = 4.0*constant::physicoChemical::sigma.value()*pow4(Tp);
141 
142  // Diffusion coefficient - created by radiation model's ::updateCoeffs()
143  const scalarField& gamma =
144  patch().lookupPatchField<volScalarField, scalar>("gammaRad");
145 
146  const boundaryRadiationProperties& boundaryRadiation =
147  boundaryRadiationProperties::New(internalField().mesh());
148 
149  const tmp<scalarField> temissivity
150  (
151  boundaryRadiation.emissivity(patch().index())
152  );
153 
154  const scalarField& emissivity = temissivity();
155 
156  const scalarField Ep(emissivity/(2.0*(scalar(2) - emissivity)));
157 
158  // Set value fraction
159  valueFraction() = 1.0/(1.0 + gamma*patch().deltaCoeffs()/Ep);
160 
161  // Restore tag
162  UPstream::msgType() = oldTag;
163 
164  mixedFvPatchScalarField::updateCoeffs();
165 }
166 
167 
169 (
170  Ostream& os
171 ) const
172 {
174  os.writeEntryIfDifferent<word>("T", "T", TName_);
175 }
176 
177 
178 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
179 
180 namespace Foam
181 {
182 namespace radiation
183 {
185  (
188  );
189 }
190 }
191 
192 // ************************************************************************* //
makePatchTypeField(fvPatchScalarField, greyDiffusiveRadiationMixedFvPatchScalarField)
const Type & value() const noexcept
Return const reference to value.
dictionary dict
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:120
static const boundaryRadiationProperties & New(const fvMesh &mesh, Args &&... args)
Get existing or create a new MeshObject.
Definition: MeshObject.C:41
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:68
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m2/K4].
static int & msgType() noexcept
Message tag of standard messages.
Definition: UPstream.H:806
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:84
void write(vtk::formatter &fmt, const Type &val, const label n=1)
Component-wise write of a value (N times)
dynamicFvMesh & mesh
A class for handling words, derived from Foam::string.
Definition: word.H:63
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
A FieldMapper for finite-volume patch fields.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:55
OBJstream os(runTime.globalPath()/outputName)
MarshakRadiationFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
virtual void operator=(const UList< scalar > &)
Definition: fvPatchField.C:352
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: areaFieldsFwd.H:42
dimensionedScalar pow4(const dimensionedScalar &ds)
const std::string patch
OpenFOAM patch number as a std::string.
const scalar gamma
Definition: EEqn.H:9
volScalarField & p
A class for managing temporary objects.
Definition: HashPtrTable.H:50
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
autoPtr< radiation::radiationModel > radiation(radiation::radiationModel::New(T))
Namespace for OpenFOAM.