MarshakRadiationFixedTemperatureFvPatchScalarField.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"
33 
34 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
35 
38 (
39  const fvPatch& p,
41 )
42 :
43  mixedFvPatchScalarField(p, iF),
44  Trad_(p.size())
45 {
46  refValue() = 0.0;
47  refGrad() = 0.0;
48  valueFraction() = 0.0;
49 }
50 
51 
54 (
56  const fvPatch& p,
58  const fvPatchFieldMapper& mapper
59 )
60 :
61  mixedFvPatchScalarField(ptf, p, iF, mapper),
62  Trad_(ptf.Trad_, mapper)
63 {}
64 
65 
68 (
69  const fvPatch& p,
71  const dictionary& dict
72 )
73 :
74  mixedFvPatchScalarField(p, iF),
75  Trad_("Trad", dict, p.size())
76 {
77  // refValue updated on each call to updateCoeffs()
78  refValue() = 4.0*constant::physicoChemical::sigma.value()*pow4(Trad_);
79 
80  // zero gradient
81  refGrad() = 0.0;
82 
83  valueFraction() = 1.0;
84 
85  fvPatchScalarField::operator=(refValue());
86 }
87 
88 
91 (
93 )
94 :
95  mixedFvPatchScalarField(ptf),
96  Trad_(ptf.Trad_)
97 {}
98 
99 
102 (
105 )
106 :
107  mixedFvPatchScalarField(ptf, iF),
108  Trad_(ptf.Trad_)
109 {}
110 
111 
112 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
113 
115 autoMap
116 (
117  const fvPatchFieldMapper& m
118 )
119 {
120  mixedFvPatchScalarField::autoMap(m);
121  Trad_.autoMap(m);
122 }
123 
124 
126 (
127  const fvPatchScalarField& ptf,
128  const labelList& addr
129 )
130 {
131  mixedFvPatchScalarField::rmap(ptf, addr);
132 
134  refCast<const MarshakRadiationFixedTemperatureFvPatchScalarField>(ptf);
135 
136  Trad_.rmap(mrptf.Trad_, addr);
137 }
138 
139 
142 {
143  if (this->updated())
144  {
145  return;
146  }
147 
148  // Since we're inside initEvaluate/evaluate there might be processor
149  // comms underway. Change the tag we use.
150  int oldTag = UPstream::msgType();
151  UPstream::msgType() = oldTag+1;
152 
153  // Re-calc reference value
154  refValue() = 4.0*constant::physicoChemical::sigma.value()*pow4(Trad_);
155 
156  // Diffusion coefficient - created by radiation model's ::updateCoeffs()
157  const scalarField& gamma =
158  patch().lookupPatchField<volScalarField, scalar>("gammaRad");
159 
160  //const scalarField temissivity = emissivity();
161  const boundaryRadiationProperties& boundaryRadiation =
163 
164  const tmp<scalarField> temissivity
165  (
166  boundaryRadiation.emissivity(patch().index())
167  );
168 
169  const scalarField& emissivity = temissivity();
170 
171  const scalarField Ep(emissivity/(2.0*(scalar(2.0) - emissivity)));
172 
173  // Set value fraction
174  valueFraction() = 1.0/(1.0 + gamma*patch().deltaCoeffs()/Ep);
175 
176  // Restore tag
177  UPstream::msgType() = oldTag;
178 
179  mixedFvPatchScalarField::updateCoeffs();
180 }
181 
182 
184 (
185  Ostream& os
186 ) const
187 {
189  Trad_.writeEntry("Trad", os);
190 }
191 
192 
193 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
194 
195 namespace Foam
196 {
197 namespace radiation
198 {
200  (
203  );
204 }
205 }
206 
207 // ************************************************************************* //
Foam::fvPatchField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:65
volFields.H
p
p
Definition: pEqn.H:62
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::radiation::MarshakRadiationFixedTemperatureFvPatchScalarField::MarshakRadiationFixedTemperatureFvPatchScalarField
MarshakRadiationFixedTemperatureFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Definition: MarshakRadiationFixedTemperatureFvPatchScalarField.C:38
dimensionedInternalField
rDeltaT dimensionedInternalField()
Foam::radiation::MarshakRadiationFixedTemperatureFvPatchScalarField::write
virtual void write(Ostream &) const
Write.
Definition: MarshakRadiationFixedTemperatureFvPatchScalarField.C:184
Foam::radiation::makePatchTypeField
makePatchTypeField(fvPatchScalarField, boundaryRadiationPropertiesFvPatchField)
fvPatchFieldMapper.H
Foam::dimensioned::value
const Type & value() const
Return const reference to value.
Definition: dimensionedType.C:261
Foam::MeshObject< fvMesh, Foam::GeometricMeshObject, boundaryRadiationProperties >::New
static const boundaryRadiationProperties & New(const fvMesh &mesh)
Definition: MeshObject.C:44
Foam::pow4
dimensionedScalar pow4(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:98
Foam::radiation::MarshakRadiationFixedTemperatureFvPatchScalarField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: MarshakRadiationFixedTemperatureFvPatchScalarField.C:141
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::radiation::MarshakRadiationFixedTemperatureFvPatchScalarField::Trad_
scalarField Trad_
Radiation temperature field.
Definition: MarshakRadiationFixedTemperatureFvPatchScalarField.H:96
MarshakRadiationFixedTemperatureFvPatchScalarField.H
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
radiation
autoPtr< radiation::radiationModel > radiation(radiation::radiationModel::New(T))
Foam::radiation::boundaryRadiationProperties
Definition: boundaryRadiationProperties.H:56
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::radiation::MarshakRadiationFixedTemperatureFvPatchScalarField
Definition: MarshakRadiationFixedTemperatureFvPatchScalarField.H:87
Foam::UPstream::msgType
static int & msgType()
Message tag of standard messages.
Definition: UPstream.H:452
physicoChemicalConstants.H
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
Foam::constant::physicoChemical::sigma
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m2/K4].
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:45
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:53
Foam::radiation::MarshakRadiationFixedTemperatureFvPatchScalarField::autoMap
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
Definition: MarshakRadiationFixedTemperatureFvPatchScalarField.C:116
Foam::radiation::MarshakRadiationFixedTemperatureFvPatchScalarField::rmap
virtual void rmap(const fvPatchScalarField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
Definition: MarshakRadiationFixedTemperatureFvPatchScalarField.C:126
radiationModel.H
write
Tcoeff write()
boundaryRadiationProperties.H
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:51