Go to the documentation of this file.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 "fixedUnburntEnthalpyFvPatchScalarField.H"
00027 #include <OpenFOAM/addToRunTimeSelectionTable.H>
00028 #include <finiteVolume/fvPatchFieldMapper.H>
00029 #include <finiteVolume/volFields.H>
00030 #include <reactionThermophysicalModels/hhuCombustionThermo.H>
00031
00032
00033
00034 namespace Foam
00035 {
00036
00037
00038
00039 fixedUnburntEnthalpyFvPatchScalarField::fixedUnburntEnthalpyFvPatchScalarField
00040 (
00041 const fvPatch& p,
00042 const DimensionedField<scalar, volMesh>& iF
00043 )
00044 :
00045 fixedValueFvPatchScalarField(p, iF)
00046 {}
00047
00048
00049 fixedUnburntEnthalpyFvPatchScalarField::fixedUnburntEnthalpyFvPatchScalarField
00050 (
00051 const fixedUnburntEnthalpyFvPatchScalarField& ptf,
00052 const fvPatch& p,
00053 const DimensionedField<scalar, volMesh>& iF,
00054 const fvPatchFieldMapper& mapper
00055 )
00056 :
00057 fixedValueFvPatchScalarField(ptf, p, iF, mapper)
00058 {}
00059
00060
00061 fixedUnburntEnthalpyFvPatchScalarField::fixedUnburntEnthalpyFvPatchScalarField
00062 (
00063 const fvPatch& p,
00064 const DimensionedField<scalar, volMesh>& iF,
00065 const dictionary& dict
00066 )
00067 :
00068 fixedValueFvPatchScalarField(p, iF, dict)
00069 {}
00070
00071
00072 fixedUnburntEnthalpyFvPatchScalarField::fixedUnburntEnthalpyFvPatchScalarField
00073 (
00074 const fixedUnburntEnthalpyFvPatchScalarField& tppsf
00075 )
00076 :
00077 fixedValueFvPatchScalarField(tppsf)
00078 {}
00079
00080
00081 fixedUnburntEnthalpyFvPatchScalarField::fixedUnburntEnthalpyFvPatchScalarField
00082 (
00083 const fixedUnburntEnthalpyFvPatchScalarField& tppsf,
00084 const DimensionedField<scalar, volMesh>& iF
00085 )
00086 :
00087 fixedValueFvPatchScalarField(tppsf, iF)
00088 {}
00089
00090
00091
00092
00093 void fixedUnburntEnthalpyFvPatchScalarField::updateCoeffs()
00094 {
00095 if (updated())
00096 {
00097 return;
00098 }
00099
00100 const hhuCombustionThermo& thermo = db().lookupObject<hhuCombustionThermo>
00101 (
00102 "thermophysicalProperties"
00103 );
00104
00105 const label patchi = patch().index();
00106
00107 fvPatchScalarField& Tw =
00108 const_cast<fvPatchScalarField&>(thermo.Tu().boundaryField()[patchi]);
00109 Tw.evaluate();
00110 operator==(thermo.hu(Tw, patchi));
00111
00112 fixedValueFvPatchScalarField::updateCoeffs();
00113 }
00114
00115
00116
00117
00118 makePatchTypeField(fvPatchScalarField, fixedUnburntEnthalpyFvPatchScalarField);
00119
00120
00121
00122 }
00123
00124