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

hhuCombustionThermo.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 "hhuCombustionThermo.H"
00027 #include <finiteVolume/fvMesh.H>
00028 #include <finiteVolume/zeroGradientFvPatchFields.H>
00029 #include <reactionThermophysicalModels/fixedUnburntEnthalpyFvPatchScalarField.H>
00030 #include <reactionThermophysicalModels/gradientUnburntEnthalpyFvPatchScalarField.H>
00031 #include <reactionThermophysicalModels/mixedUnburntEnthalpyFvPatchScalarField.H>
00032 
00033 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00034 
00035 namespace Foam
00036 {
00037 
00038 /* * * * * * * * * * * * * * * private static data * * * * * * * * * * * * * */
00039 
00040 defineTypeNameAndDebug(hhuCombustionThermo, 0);
00041 defineRunTimeSelectionTable(hhuCombustionThermo, fvMesh);
00042 
00043 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
00044 
00045 wordList hhuCombustionThermo::huBoundaryTypes()
00046 {
00047     const volScalarField::GeometricBoundaryField& tbf = Tu_.boundaryField();
00048 
00049     wordList hbt = tbf.types();
00050 
00051     forAll(tbf, patchi)
00052     {
00053         if (isA<fixedValueFvPatchScalarField>(tbf[patchi]))
00054         {
00055             hbt[patchi] = fixedUnburntEnthalpyFvPatchScalarField::typeName;
00056         }
00057         else if
00058         (
00059             isA<zeroGradientFvPatchScalarField>(tbf[patchi])
00060          || isA<fixedGradientFvPatchScalarField>(tbf[patchi])
00061         )
00062         {
00063             hbt[patchi] = gradientUnburntEnthalpyFvPatchScalarField::typeName;
00064         }
00065         else if (isA<mixedFvPatchScalarField>(tbf[patchi]))
00066         {
00067             hbt[patchi] = mixedUnburntEnthalpyFvPatchScalarField::typeName;
00068         }
00069     }
00070 
00071     return hbt;
00072 }
00073 
00074 void hhuCombustionThermo::huBoundaryCorrection(volScalarField& hu)
00075 {
00076     volScalarField::GeometricBoundaryField& hbf = hu.boundaryField();
00077 
00078     forAll(hbf, patchi)
00079     {
00080         if
00081         (
00082             isA<gradientUnburntEnthalpyFvPatchScalarField>(hbf[patchi])
00083         )
00084         {
00085             refCast<gradientUnburntEnthalpyFvPatchScalarField>(hbf[patchi])
00086                 .gradient() = hbf[patchi].fvPatchField::snGrad();
00087         }
00088         else if
00089         (
00090             isA<mixedUnburntEnthalpyFvPatchScalarField>(hbf[patchi])
00091         )
00092         {
00093             refCast<mixedUnburntEnthalpyFvPatchScalarField>(hbf[patchi])
00094                 .refGrad() = hbf[patchi].fvPatchField::snGrad();
00095         }
00096     }
00097 }
00098 
00099 
00100 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
00101 
00102 hhuCombustionThermo::hhuCombustionThermo(const fvMesh& mesh)
00103 :
00104     hCombustionThermo(mesh),
00105 
00106     Tu_
00107     (
00108         IOobject
00109         (
00110             "Tu",
00111             mesh.time().timeName(),
00112             mesh,
00113             IOobject::MUST_READ,
00114             IOobject::AUTO_WRITE
00115         ),
00116         mesh
00117     ),
00118 
00119     hu_
00120     (
00121         IOobject
00122         (
00123             "hu",
00124             mesh.time().timeName(),
00125             mesh,
00126             IOobject::NO_READ,
00127             IOobject::NO_WRITE
00128         ),
00129         mesh,
00130         dimensionSet(0, 2, -2, 0, 0),
00131         huBoundaryTypes()
00132     )
00133 {}
00134 
00135 
00136 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
00137 
00138 hhuCombustionThermo::~hhuCombustionThermo()
00139 {}
00140 
00141 
00142 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00143 
00144 } // End namespace Foam
00145 
00146 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines