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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2015 OpenFOAM Foundation
9  Copyright (C) 2016 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  Trad_(p.size())
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  Trad_(ptf.Trad_, mapper)
66 {}
67 
68 
71 (
72  const fvPatch& p,
74  const dictionary& dict
75 )
76 :
77  mixedFvPatchScalarField(p, iF),
78  Trad_("Trad", dict, p.size())
79 {
80  // refValue updated on each call to updateCoeffs()
81  refValue() = 4.0*constant::physicoChemical::sigma.value()*pow4(Trad_);
82 
83  // zero gradient
84  refGrad() = 0.0;
85 
86  valueFraction() = 1.0;
87 
88  fvPatchScalarField::operator=(refValue());
89 }
90 
91 
94 (
96 )
97 :
98  mixedFvPatchScalarField(ptf),
99  Trad_(ptf.Trad_)
100 {}
101 
102 
105 (
108 )
109 :
110  mixedFvPatchScalarField(ptf, iF),
111  Trad_(ptf.Trad_)
112 {}
113 
114 
115 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
116 
118 autoMap
119 (
120  const fvPatchFieldMapper& m
121 )
122 {
123  mixedFvPatchScalarField::autoMap(m);
124  Trad_.autoMap(m);
125 }
126 
127 
129 (
130  const fvPatchScalarField& ptf,
131  const labelList& addr
132 )
133 {
134  mixedFvPatchScalarField::rmap(ptf, addr);
135 
137  refCast<const MarshakRadiationFixedTemperatureFvPatchScalarField>(ptf);
138 
139  Trad_.rmap(mrptf.Trad_, addr);
140 }
141 
142 
145 {
146  if (this->updated())
147  {
148  return;
149  }
150 
151  // Since we're inside initEvaluate/evaluate there might be processor
152  // comms underway. Change the tag we use.
153  int oldTag = UPstream::msgType();
154  UPstream::msgType() = oldTag+1;
155 
156  // Re-calc reference value
157  refValue() = 4.0*constant::physicoChemical::sigma.value()*pow4(Trad_);
158 
159  // Diffusion coefficient - created by radiation model's ::updateCoeffs()
160  const scalarField& gamma =
161  patch().lookupPatchField<volScalarField, scalar>("gammaRad");
162 
163  //const scalarField temissivity = emissivity();
164  const boundaryRadiationProperties& boundaryRadiation =
165  boundaryRadiationProperties::New(internalField().mesh());
166 
167  const tmp<scalarField> temissivity
168  (
169  boundaryRadiation.emissivity(patch().index())
170  );
171 
172  const scalarField& emissivity = temissivity();
173 
174  const scalarField Ep(emissivity/(2.0*(scalar(2) - emissivity)));
175 
176  // Set value fraction
177  valueFraction() = 1.0/(1.0 + gamma*patch().deltaCoeffs()/Ep);
178 
179  // Restore tag
180  UPstream::msgType() = oldTag;
181 
182  mixedFvPatchScalarField::updateCoeffs();
183 }
184 
185 
187 (
188  Ostream& os
189 ) const
190 {
192  Trad_.writeEntry("Trad", os);
193 }
194 
195 
196 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
197 
198 namespace Foam
199 {
200 namespace radiation
201 {
203  (
206  );
207 }
208 }
209 
210 // ************************************************************************* //
Foam::fvPatchField< scalar >
volFields.H
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::radiation::makePatchTypeField
makePatchTypeField(fvPatchScalarField, greyDiffusiveRadiationMixedFvPatchScalarField)
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::radiation::MarshakRadiationFixedTemperatureFvPatchScalarField::MarshakRadiationFixedTemperatureFvPatchScalarField
MarshakRadiationFixedTemperatureFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Definition: MarshakRadiationFixedTemperatureFvPatchScalarField.C:41
Foam::radiation::MarshakRadiationFixedTemperatureFvPatchScalarField::write
virtual void write(Ostream &) const
Write.
Definition: MarshakRadiationFixedTemperatureFvPatchScalarField.C:187
Foam::MeshObject< fvMesh, Foam::GeometricMeshObject, boundaryRadiationProperties >::New
static const boundaryRadiationProperties & New(const fvMesh &mesh, Args &&... args)
Get existing or create a new MeshObject.
Definition: MeshObject.C:48
fvPatchFieldMapper.H
Foam::dimensioned::value
const Type & value() const
Return const reference to value.
Definition: dimensionedType.C:434
Foam::pow4
dimensionedScalar pow4(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:100
Foam::radiation::MarshakRadiationFixedTemperatureFvPatchScalarField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: MarshakRadiationFixedTemperatureFvPatchScalarField.C:144
Foam::Field< scalar >
MarshakRadiationFixedTemperatureFvPatchScalarField.H
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:65
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
radiation
autoPtr< radiation::radiationModel > radiation(radiation::radiationModel::New(T))
Foam::radiation::boundaryRadiationProperties
Boundary radiation properties holder.
Definition: boundaryRadiationProperties.H:56
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
os
OBJstream os(runTime.globalPath()/outputName)
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::radiation::MarshakRadiationFixedTemperatureFvPatchScalarField
Definition: MarshakRadiationFixedTemperatureFvPatchScalarField.H:88
Foam::UPstream::msgType
static int & msgType() noexcept
Message tag of standard messages.
Definition: UPstream.H:540
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
physicoChemicalConstants.H
Foam::List< label >
gamma
const scalar gamma
Definition: EEqn.H:9
Foam::constant::physicoChemical::sigma
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m2/K4].
Foam::vtk::write
void write(vtk::formatter &fmt, const Type &val, const label n=1)
Component-wise write of a value (N times)
Definition: foamVtkOutputTemplates.C:36
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:47
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
sigma
dimensionedScalar sigma("sigma", dimMass/sqr(dimTime), transportProperties)
Foam::radiation::MarshakRadiationFixedTemperatureFvPatchScalarField::autoMap
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
Definition: MarshakRadiationFixedTemperatureFvPatchScalarField.C:119
Foam::radiation::MarshakRadiationFixedTemperatureFvPatchScalarField::rmap
virtual void rmap(const fvPatchScalarField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
Definition: MarshakRadiationFixedTemperatureFvPatchScalarField.C:129
radiationModel.H
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:54