FreeFOAM The Cross-Platform CFD Toolkit
Hosted by SourceForge:
Get FreeFOAM at SourceForge.net.
            Fast, secure and Free Open Source software downloads

gradientUnburntEnthalpyFvPatchScalarField.C

Go to the documentation of this file.
00001 /*---------------------------------------------------------------------------*\
00002   =========                 |
00003   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
00004    \\    /   O peration     |
00005     \\  /    A nd           | Copyright (C) 1991-2010 OpenCFD Ltd.
00006      \\/     M anipulation  |
00007 -------------------------------------------------------------------------------
00008 License
00009     This file is part of OpenFOAM.
00010 
00011     OpenFOAM is free software: you can redistribute it and/or modify it
00012     under the terms of the GNU General Public License as published by
00013     the Free Software Foundation, either version 3 of the License, or
00014     (at your option) any later version.
00015 
00016     OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
00017     ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
00018     FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
00019     for more details.
00020 
00021     You should have received a copy of the GNU General Public License
00022     along with OpenFOAM.  If not, see <http://www.gnu.org/licenses/>.
00023 
00024 \*---------------------------------------------------------------------------*/
00025 
00026 #include "gradientUnburntEnthalpyFvPatchScalarField.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 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
00038 
00039 gradientUnburntEnthalpyFvPatchScalarField::
00040 gradientUnburntEnthalpyFvPatchScalarField
00041 (
00042     const fvPatch& p,
00043     const DimensionedField<scalar, volMesh>& iF
00044 )
00045 :
00046     fixedGradientFvPatchScalarField(p, iF)
00047 {}
00048 
00049 
00050 gradientUnburntEnthalpyFvPatchScalarField::
00051 gradientUnburntEnthalpyFvPatchScalarField
00052 (
00053     const gradientUnburntEnthalpyFvPatchScalarField& ptf,
00054     const fvPatch& p,
00055     const DimensionedField<scalar, volMesh>& iF,
00056     const fvPatchFieldMapper& mapper
00057 )
00058 :
00059     fixedGradientFvPatchScalarField(ptf, p, iF, mapper)
00060 {}
00061 
00062 
00063 gradientUnburntEnthalpyFvPatchScalarField::
00064 gradientUnburntEnthalpyFvPatchScalarField
00065 (
00066     const fvPatch& p,
00067     const DimensionedField<scalar, volMesh>& iF,
00068     const dictionary& dict
00069 )
00070 :
00071     fixedGradientFvPatchScalarField(p, iF, dict)
00072 {}
00073 
00074 
00075 gradientUnburntEnthalpyFvPatchScalarField::
00076 gradientUnburntEnthalpyFvPatchScalarField
00077 (
00078     const gradientUnburntEnthalpyFvPatchScalarField& tppsf
00079 )
00080 :
00081     fixedGradientFvPatchScalarField(tppsf)
00082 {}
00083 
00084 
00085 gradientUnburntEnthalpyFvPatchScalarField::
00086 gradientUnburntEnthalpyFvPatchScalarField
00087 (
00088     const gradientUnburntEnthalpyFvPatchScalarField& tppsf,
00089     const DimensionedField<scalar, volMesh>& iF
00090 )
00091 :
00092     fixedGradientFvPatchScalarField(tppsf, iF)
00093 {}
00094 
00095 
00096 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
00097 
00098 void gradientUnburntEnthalpyFvPatchScalarField::updateCoeffs()
00099 {
00100     if (updated())
00101     {
00102         return;
00103     }
00104 
00105     const hhuCombustionThermo& thermo = db().lookupObject<hhuCombustionThermo>
00106     (
00107         "thermophysicalProperties"
00108     );
00109     
00110     const label patchi = patch().index();
00111 
00112     fvPatchScalarField& Tw = 
00113         const_cast<fvPatchScalarField&>(thermo.Tu().boundaryField()[patchi]);
00114 
00115     Tw.evaluate();
00116 
00117     gradient() = thermo.Cp(Tw, patchi)*Tw.snGrad()
00118       + patch().deltaCoeffs()*
00119         (
00120             thermo.hu(Tw, patchi)
00121           - thermo.hu(Tw, patch().faceCells())
00122         );
00123 
00124     fixedGradientFvPatchScalarField::updateCoeffs();
00125 }
00126 
00127 
00128 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00129 
00130 makePatchTypeField
00131 (
00132     fvPatchScalarField,
00133     gradientUnburntEnthalpyFvPatchScalarField
00134 );
00135 
00136 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00137 
00138 } // End namespace Foam
00139 
00140 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines