1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 1991-2009 OpenCFD Ltd.
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 the
13 Free Software Foundation; either version 2 of the License, or (at your
14 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, write to the Free Software Foundation,
23 Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 \*---------------------------------------------------------------------------*/
27 #include "turbulentHeatFluxTemperatureFvPatchScalarField.H"
28 #include "addToRunTimeSelectionTable.H"
29 #include "fvPatchFieldMapper.H"
30 #include "volFields.H"
33 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37 namespace compressible
40 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
44 NamedEnum<turbulentHeatFluxTemperatureFvPatchScalarField::heatSourceType, 2>::
52 NamedEnum<turbulentHeatFluxTemperatureFvPatchScalarField::heatSourceType, 2>
53 turbulentHeatFluxTemperatureFvPatchScalarField::heatSourceTypeNames_;
56 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
58 turbulentHeatFluxTemperatureFvPatchScalarField::
59 turbulentHeatFluxTemperatureFvPatchScalarField
62 const DimensionedField<scalar, volMesh>& iF
65 fixedGradientFvPatchScalarField(p, iF),
71 turbulentHeatFluxTemperatureFvPatchScalarField::
72 turbulentHeatFluxTemperatureFvPatchScalarField
74 const turbulentHeatFluxTemperatureFvPatchScalarField& ptf,
76 const DimensionedField<scalar, volMesh>& iF,
77 const fvPatchFieldMapper& mapper
80 fixedGradientFvPatchScalarField(ptf, p, iF, mapper),
81 heatSource_(ptf.heatSource_),
86 turbulentHeatFluxTemperatureFvPatchScalarField::
87 turbulentHeatFluxTemperatureFvPatchScalarField
90 const DimensionedField<scalar, volMesh>& iF,
91 const dictionary& dict
94 fixedGradientFvPatchScalarField(p, iF),
95 heatSource_(heatSourceTypeNames_.read(dict.lookup("heatSource"))),
96 q_("q", dict, p.size())
98 fvPatchField<scalar>::operator=(patchInternalField());
103 turbulentHeatFluxTemperatureFvPatchScalarField::
104 turbulentHeatFluxTemperatureFvPatchScalarField
106 const turbulentHeatFluxTemperatureFvPatchScalarField& thftpsf
109 fixedGradientFvPatchScalarField(thftpsf),
110 heatSource_(thftpsf.heatSource_),
115 turbulentHeatFluxTemperatureFvPatchScalarField::
116 turbulentHeatFluxTemperatureFvPatchScalarField
118 const turbulentHeatFluxTemperatureFvPatchScalarField& thftpsf,
119 const DimensionedField<scalar, volMesh>& iF
122 fixedGradientFvPatchScalarField(thftpsf, iF),
123 heatSource_(thftpsf.heatSource_),
128 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
130 void turbulentHeatFluxTemperatureFvPatchScalarField::autoMap
132 const fvPatchFieldMapper& m
135 fixedGradientFvPatchScalarField::autoMap(m);
140 void turbulentHeatFluxTemperatureFvPatchScalarField::rmap
142 const fvPatchScalarField& ptf,
143 const labelList& addr
146 fixedGradientFvPatchScalarField::rmap(ptf, addr);
148 const turbulentHeatFluxTemperatureFvPatchScalarField& thftptf =
149 refCast<const turbulentHeatFluxTemperatureFvPatchScalarField>
154 q_.rmap(thftptf.q_, addr);
158 void turbulentHeatFluxTemperatureFvPatchScalarField::updateCoeffs()
165 const label patchI = patch().index();
167 const RASModel& rasModel = db().lookupObject<RASModel>("RASProperties");
169 const scalarField alphaEffp = rasModel.alphaEff()().boundaryField()[patchI];
171 // const scalarField& Tp = thermo.T().boundaryField()[patchI];
172 const scalarField& Tp = *this;
174 const scalarField Cpp = rasModel.thermo().Cp(Tp, patchI);
180 const scalar Ap = gSum(patch().magSf());
181 gradient() = q_/(Ap*Cpp*alphaEffp);
186 gradient() = q_/(Cpp*alphaEffp);
193 "turbulentHeatFluxTemperatureFvPatchScalarField"
196 "const DimensionedField<scalar, volMesh>&, "
199 ) << "Unknown heat source type. Valid types are: "
200 << heatSourceTypeNames_ << nl << exit(FatalError);
204 fixedGradientFvPatchScalarField::updateCoeffs();
208 void turbulentHeatFluxTemperatureFvPatchScalarField::write
213 fvPatchScalarField::write(os);
214 os.writeKeyword("heatSource") << heatSourceTypeNames_[heatSource_]
215 << token::END_STATEMENT << nl;
216 q_.writeEntry("q", os);
217 gradient().writeEntry("gradient", os);
218 writeEntry("value", os);
222 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
227 turbulentHeatFluxTemperatureFvPatchScalarField
231 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
233 } // End namespace compressible
234 } // End namespace Foam
237 // ************************************************************************* //