1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 2011 OpenFOAM Foundation
7 -------------------------------------------------------------------------------
9 This file is part of OpenFOAM.
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.
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
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/>.
24 \*---------------------------------------------------------------------------*/
26 #include "MarshakRadiationFixedTMixedFvPatchScalarField.H"
27 #include "addToRunTimeSelectionTable.H"
28 #include "fvPatchFieldMapper.H"
29 #include "volFields.H"
32 #include "radiationModel.H"
33 #include "physicoChemicalConstants.H"
35 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
37 Foam::MarshakRadiationFixedTMixedFvPatchScalarField::
38 MarshakRadiationFixedTMixedFvPatchScalarField
41 const DimensionedField<scalar, volMesh>& iF
44 mixedFvPatchScalarField(p, iF),
45 radiationCoupledBase(p, "undefined", scalarField::null()),
50 valueFraction() = 0.0;
54 Foam::MarshakRadiationFixedTMixedFvPatchScalarField::
55 MarshakRadiationFixedTMixedFvPatchScalarField
57 const MarshakRadiationFixedTMixedFvPatchScalarField& ptf,
59 const DimensionedField<scalar, volMesh>& iF,
60 const fvPatchFieldMapper& mapper
63 mixedFvPatchScalarField(ptf, p, iF, mapper),
67 ptf.emissivityMethod(),
70 Trad_(ptf.Trad_, mapper)
74 Foam::MarshakRadiationFixedTMixedFvPatchScalarField::
75 MarshakRadiationFixedTMixedFvPatchScalarField
78 const DimensionedField<scalar, volMesh>& iF,
79 const dictionary& dict
82 mixedFvPatchScalarField(p, iF),
83 radiationCoupledBase(p, dict),
84 Trad_("Trad", dict, p.size())
86 // refValue updated on each call to updateCoeffs()
87 refValue() = 4.0*constant::physicoChemical::sigma.value()*pow4(Trad_);
92 valueFraction() = 1.0;
94 fvPatchScalarField::operator=(refValue());
98 Foam::MarshakRadiationFixedTMixedFvPatchScalarField::
99 MarshakRadiationFixedTMixedFvPatchScalarField
101 const MarshakRadiationFixedTMixedFvPatchScalarField& ptf
104 mixedFvPatchScalarField(ptf),
108 ptf.emissivityMethod(),
115 Foam::MarshakRadiationFixedTMixedFvPatchScalarField::
116 MarshakRadiationFixedTMixedFvPatchScalarField
118 const MarshakRadiationFixedTMixedFvPatchScalarField& ptf,
119 const DimensionedField<scalar, volMesh>& iF
122 mixedFvPatchScalarField(ptf, iF),
126 ptf.emissivityMethod(),
133 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
135 void Foam::MarshakRadiationFixedTMixedFvPatchScalarField::autoMap
137 const fvPatchFieldMapper& m
140 mixedFvPatchScalarField::autoMap(m);
145 void Foam::MarshakRadiationFixedTMixedFvPatchScalarField::rmap
147 const fvPatchScalarField& ptf,
148 const labelList& addr
151 mixedFvPatchScalarField::rmap(ptf, addr);
153 const MarshakRadiationFixedTMixedFvPatchScalarField& mrptf =
154 refCast<const MarshakRadiationFixedTMixedFvPatchScalarField>(ptf);
156 Trad_.rmap(mrptf.Trad_, addr);
160 void Foam::MarshakRadiationFixedTMixedFvPatchScalarField::updateCoeffs()
167 // Since we're inside initEvaluate/evaluate there might be processor
168 // comms underway. Change the tag we use.
169 int oldTag = UPstream::msgType();
170 UPstream::msgType() = oldTag+1;
172 // Re-calc reference value
173 refValue() = 4.0*constant::physicoChemical::sigma.value()*pow4(Trad_);
175 // Diffusion coefficient - created by radiation model's ::updateCoeffs()
176 const scalarField& gamma =
177 patch().lookupPatchField<volScalarField, scalar>("gammaRad");
179 const scalarField temissivity = emissivity();
181 const scalarField Ep(temissivity/(2.0*(scalar(2.0) - temissivity)));
183 // Set value fraction
184 valueFraction() = 1.0/(1.0 + gamma*patch().deltaCoeffs()/Ep);
187 UPstream::msgType() = oldTag;
189 mixedFvPatchScalarField::updateCoeffs();
193 void Foam::MarshakRadiationFixedTMixedFvPatchScalarField::write
198 mixedFvPatchScalarField::write(os);
199 radiationCoupledBase::write(os);
200 Trad_.writeEntry("Trad", os);
204 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
211 MarshakRadiationFixedTMixedFvPatchScalarField
216 // ************************************************************************* //