1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 2011 OpenFOAM Foundation
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
13 the Free Software Foundation, either version 3 of the License, or
14 (at your 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, see <http://www.gnu.org/licenses/>.
24 \*---------------------------------------------------------------------------*/
26 #include "gradientUnburntEnthalpyFvPatchScalarField.H"
27 #include "addToRunTimeSelectionTable.H"
28 #include "fvPatchFieldMapper.H"
29 #include "volFields.H"
30 #include "hhuCombustionThermo.H"
32 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
34 Foam::gradientUnburntEnthalpyFvPatchScalarField::
35 gradientUnburntEnthalpyFvPatchScalarField
38 const DimensionedField<scalar, volMesh>& iF
41 fixedGradientFvPatchScalarField(p, iF)
45 Foam::gradientUnburntEnthalpyFvPatchScalarField::
46 gradientUnburntEnthalpyFvPatchScalarField
48 const gradientUnburntEnthalpyFvPatchScalarField& ptf,
50 const DimensionedField<scalar, volMesh>& iF,
51 const fvPatchFieldMapper& mapper
54 fixedGradientFvPatchScalarField(ptf, p, iF, mapper)
58 Foam::gradientUnburntEnthalpyFvPatchScalarField::
59 gradientUnburntEnthalpyFvPatchScalarField
62 const DimensionedField<scalar, volMesh>& iF,
63 const dictionary& dict
66 fixedGradientFvPatchScalarField(p, iF, dict)
70 Foam::gradientUnburntEnthalpyFvPatchScalarField::
71 gradientUnburntEnthalpyFvPatchScalarField
73 const gradientUnburntEnthalpyFvPatchScalarField& tppsf
76 fixedGradientFvPatchScalarField(tppsf)
80 Foam::gradientUnburntEnthalpyFvPatchScalarField::
81 gradientUnburntEnthalpyFvPatchScalarField
83 const gradientUnburntEnthalpyFvPatchScalarField& tppsf,
84 const DimensionedField<scalar, volMesh>& iF
87 fixedGradientFvPatchScalarField(tppsf, iF)
91 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
93 void Foam::gradientUnburntEnthalpyFvPatchScalarField::updateCoeffs()
100 const hhuCombustionThermo& thermo = db().lookupObject<hhuCombustionThermo>
102 "thermophysicalProperties"
105 const label patchi = patch().index();
107 fvPatchScalarField& Tw =
108 const_cast<fvPatchScalarField&>(thermo.Tu().boundaryField()[patchi]);
112 gradient() = thermo.Cp(Tw, patchi)*Tw.snGrad()
113 + patch().deltaCoeffs()*
115 thermo.hu(Tw, patchi)
116 - thermo.hu(Tw, patch().faceCells())
119 fixedGradientFvPatchScalarField::updateCoeffs();
123 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
130 gradientUnburntEnthalpyFvPatchScalarField
134 // ************************************************************************* //