00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026 #include "turbulentHeatFluxTemperatureFvPatchScalarField.H"
00027 #include <OpenFOAM/addToRunTimeSelectionTable.H>
00028 #include <finiteVolume/fvPatchFieldMapper.H>
00029 #include <finiteVolume/volFields.H>
00030 #include <compressibleRASModels/RASModel.H>
00031
00032
00033
00034 namespace Foam
00035 {
00036 namespace compressible
00037 {
00038
00039
00040
00041 template<>
00042 const char*
00043 NamedEnum<turbulentHeatFluxTemperatureFvPatchScalarField::heatSourceType, 2>::
00044 names[] =
00045 {
00046 "power",
00047 "flux"
00048 };
00049
00050 const
00051 NamedEnum<turbulentHeatFluxTemperatureFvPatchScalarField::heatSourceType, 2>
00052 turbulentHeatFluxTemperatureFvPatchScalarField::heatSourceTypeNames_;
00053
00054
00055
00056
00057 turbulentHeatFluxTemperatureFvPatchScalarField::
00058 turbulentHeatFluxTemperatureFvPatchScalarField
00059 (
00060 const fvPatch& p,
00061 const DimensionedField<scalar, volMesh>& iF
00062 )
00063 :
00064 fixedGradientFvPatchScalarField(p, iF),
00065 heatSource_(hsPower),
00066 q_(p.size(), 0.0)
00067 {}
00068
00069
00070 turbulentHeatFluxTemperatureFvPatchScalarField::
00071 turbulentHeatFluxTemperatureFvPatchScalarField
00072 (
00073 const turbulentHeatFluxTemperatureFvPatchScalarField& ptf,
00074 const fvPatch& p,
00075 const DimensionedField<scalar, volMesh>& iF,
00076 const fvPatchFieldMapper& mapper
00077 )
00078 :
00079 fixedGradientFvPatchScalarField(ptf, p, iF, mapper),
00080 heatSource_(ptf.heatSource_),
00081 q_(ptf.q_, mapper)
00082 {}
00083
00084
00085 turbulentHeatFluxTemperatureFvPatchScalarField::
00086 turbulentHeatFluxTemperatureFvPatchScalarField
00087 (
00088 const fvPatch& p,
00089 const DimensionedField<scalar, volMesh>& iF,
00090 const dictionary& dict
00091 )
00092 :
00093 fixedGradientFvPatchScalarField(p, iF),
00094 heatSource_(heatSourceTypeNames_.read(dict.lookup("heatSource"))),
00095 q_("q", dict, p.size())
00096 {
00097 fvPatchField<scalar>::operator=(patchInternalField());
00098 gradient() = 0.0;
00099 }
00100
00101
00102 turbulentHeatFluxTemperatureFvPatchScalarField::
00103 turbulentHeatFluxTemperatureFvPatchScalarField
00104 (
00105 const turbulentHeatFluxTemperatureFvPatchScalarField& thftpsf
00106 )
00107 :
00108 fixedGradientFvPatchScalarField(thftpsf),
00109 heatSource_(thftpsf.heatSource_),
00110 q_(thftpsf.q_)
00111 {}
00112
00113
00114 turbulentHeatFluxTemperatureFvPatchScalarField::
00115 turbulentHeatFluxTemperatureFvPatchScalarField
00116 (
00117 const turbulentHeatFluxTemperatureFvPatchScalarField& thftpsf,
00118 const DimensionedField<scalar, volMesh>& iF
00119 )
00120 :
00121 fixedGradientFvPatchScalarField(thftpsf, iF),
00122 heatSource_(thftpsf.heatSource_),
00123 q_(thftpsf.q_)
00124 {}
00125
00126
00127
00128
00129 void turbulentHeatFluxTemperatureFvPatchScalarField::autoMap
00130 (
00131 const fvPatchFieldMapper& m
00132 )
00133 {
00134 fixedGradientFvPatchScalarField::autoMap(m);
00135 q_.autoMap(m);
00136 }
00137
00138
00139 void turbulentHeatFluxTemperatureFvPatchScalarField::rmap
00140 (
00141 const fvPatchScalarField& ptf,
00142 const labelList& addr
00143 )
00144 {
00145 fixedGradientFvPatchScalarField::rmap(ptf, addr);
00146
00147 const turbulentHeatFluxTemperatureFvPatchScalarField& thftptf =
00148 refCast<const turbulentHeatFluxTemperatureFvPatchScalarField>
00149 (
00150 ptf
00151 );
00152
00153 q_.rmap(thftptf.q_, addr);
00154 }
00155
00156
00157 void turbulentHeatFluxTemperatureFvPatchScalarField::updateCoeffs()
00158 {
00159 if (updated())
00160 {
00161 return;
00162 }
00163
00164 const label patchI = patch().index();
00165
00166 const RASModel& rasModel = db().lookupObject<RASModel>("RASProperties");
00167
00168 const scalarField alphaEffp = rasModel.alphaEff()().boundaryField()[patchI];
00169
00170
00171 const scalarField& Tp = *this;
00172
00173 const scalarField Cpp = rasModel.thermo().Cp(Tp, patchI);
00174
00175 switch (heatSource_)
00176 {
00177 case hsPower:
00178 {
00179 const scalar Ap = gSum(patch().magSf());
00180 gradient() = q_/(Ap*Cpp*alphaEffp);
00181 break;
00182 }
00183 case hsFlux:
00184 {
00185 gradient() = q_/(Cpp*alphaEffp);
00186 break;
00187 }
00188 default:
00189 {
00190 FatalErrorIn
00191 (
00192 "turbulentHeatFluxTemperatureFvPatchScalarField"
00193 "("
00194 "const fvPatch&, "
00195 "const DimensionedField<scalar, volMesh>&, "
00196 "const dictionary&"
00197 ")"
00198 ) << "Unknown heat source type. Valid types are: "
00199 << heatSourceTypeNames_ << nl << exit(FatalError);
00200 }
00201 }
00202
00203 fixedGradientFvPatchScalarField::updateCoeffs();
00204 }
00205
00206
00207 void turbulentHeatFluxTemperatureFvPatchScalarField::write
00208 (
00209 Ostream& os
00210 ) const
00211 {
00212 fvPatchScalarField::write(os);
00213 os.writeKeyword("heatSource") << heatSourceTypeNames_[heatSource_]
00214 << token::END_STATEMENT << nl;
00215 q_.writeEntry("q", os);
00216 gradient().writeEntry("gradient", os);
00217 writeEntry("value", os);
00218 }
00219
00220
00221
00222
00223 makePatchTypeField
00224 (
00225 fvPatchScalarField,
00226 turbulentHeatFluxTemperatureFvPatchScalarField
00227 );
00228
00229
00230
00231
00232 }
00233 }
00234
00235
00236
00237