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 | Copyright (C) 2011-2015 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 #include "radiationModel.H"
32 
33 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
34 
36 (
37  const fvPatch& p,
39 )
40 :
41  mixedFvPatchScalarField(p, iF),
42  radiationCoupledBase(p, "undefined", scalarField::null()),
43  TName_("T")
44 {
45  refValue() = 0.0;
46  refGrad() = 0.0;
47  valueFraction() = 0.0;
48 }
49 
50 
52 (
54  const fvPatch& p,
56  const fvPatchFieldMapper& mapper
57 )
58 :
59  mixedFvPatchScalarField(ptf, p, iF, mapper),
61  (
62  p,
63  ptf.emissivityMethod(),
64  ptf.emissivity_,
65  mapper
66  ),
67  TName_(ptf.TName_)
68 {}
69 
70 
72 (
73  const fvPatch& p,
75  const dictionary& dict
76 )
77 :
78  mixedFvPatchScalarField(p, iF),
79  radiationCoupledBase(p, dict),
80  TName_(dict.lookupOrDefault<word>("T", "T"))
81 {
82  if (dict.found("value"))
83  {
84  refValue() = scalarField("value", dict, p.size());
85  }
86  else
87  {
88  refValue() = 0.0;
89  }
90 
91  // zero gradient
92  refGrad() = 0.0;
93 
94  valueFraction() = 1.0;
95 
96  fvPatchScalarField::operator=(refValue());
97 }
98 
99 
101 (
103 )
104 :
105  mixedFvPatchScalarField(ptf),
107  (
108  ptf.patch(),
109  ptf.emissivityMethod(),
110  ptf.emissivity_
111  ),
112  TName_(ptf.TName_)
113 {}
114 
115 
117 (
120 )
121 :
122  mixedFvPatchScalarField(ptf, iF),
124  (
125  ptf.patch(),
126  ptf.emissivityMethod(),
127  ptf.emissivity_
128  ),
129  TName_(ptf.TName_)
130 {}
131 
132 
133 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
134 
136 (
137  const fvPatchFieldMapper& m
138 )
139 {
140  mixedFvPatchScalarField::autoMap(m);
141  radiationCoupledBase::autoMap(m);
142 }
143 
144 
146 (
147  const fvPatchScalarField& ptf,
148  const labelList& addr
149 )
150 {
151  mixedFvPatchScalarField::rmap(ptf, addr);
152  radiationCoupledBase::rmap(ptf, addr);
153 }
154 
155 
157 {
158  if (this->updated())
159  {
160  return;
161  }
162 
163  // Since we're inside initEvaluate/evaluate there might be processor
164  // comms underway. Change the tag we use.
165  int oldTag = UPstream::msgType();
166  UPstream::msgType() = oldTag+1;
167 
168  // Temperature field
169  const scalarField& Tp =
170  patch().lookupPatchField<volScalarField, scalar>(TName_);
171 
172  // Re-calc reference value
173  refValue() = 4.0*constant::physicoChemical::sigma.value()*pow4(Tp);
174 
175  // Diffusion coefficient - created by radiation model's ::updateCoeffs()
176  const scalarField& gamma =
177  patch().lookupPatchField<volScalarField, scalar>("gammaRad");
178 
179  const scalarField temissivity = emissivity();
180 
181  const scalarField Ep(temissivity/(2.0*(2.0 - temissivity)));
182 
183  // Set value fraction
184  valueFraction() = 1.0/(1.0 + gamma*patch().deltaCoeffs()/Ep);
185 
186  // Restore tag
187  UPstream::msgType() = oldTag;
188 
189  mixedFvPatchScalarField::updateCoeffs();
190 }
191 
192 
194 {
197  writeEntryIfDifferent<word>(os, "T", "T", TName_);
198 }
199 
200 
201 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
202 
203 namespace Foam
204 {
206  (
209  );
210 }
211 
212 
213 // ************************************************************************* //
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
virtual void rmap(const fvPatchScalarField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m2/K4].
scalarField emissivity() const
Calculate corresponding emissivity field.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:65
const Type & value() const
Return const reference to value.
static int & msgType()
Message tag of standard messages.
Definition: UPstream.H:464
scalarField emissivity_
Emissivity.
Macros for easy insertion into run-time selection tables.
void write(Ostream &) const
Write.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
makePatchTypeField(fvPatchVectorField, SRFFreestreamVelocityFvPatchVectorField)
A class for handling words, derived from string.
Definition: word.H:59
word emissivityMethod() const
Method to obtain emissivity.
Foam::fvPatchFieldMapper.
virtual label size() const
Return size.
Definition: fvPatch.H:161
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:306
Common functions to emissivity. It gets supplied from lookup into a dictionary or calculated by the s...
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
volScalarField scalarField(fieldObject, mesh)
MarshakRadiationFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
dimensionedScalar pow4(const dimensionedScalar &ds)
A &#39;mixed&#39; boundary condition that implements a Marshak condition for the incident radiation field (us...
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
runTime write()
Namespace for OpenFOAM.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.