gradientUnburntEnthalpyFvPatchScalarField.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2018 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
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.
15 
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
19  for more details.
20 
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/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
28 #include "fvPatchFieldMapper.H"
29 #include "volFields.H"
30 #include "psiuReactionThermo.H"
31 
32 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
33 
36 (
37  const fvPatch& p,
39 )
40 :
41  fixedGradientFvPatchScalarField(p, iF)
42 {}
43 
44 
47 (
49  const fvPatch& p,
51  const fvPatchFieldMapper& mapper
52 )
53 :
54  fixedGradientFvPatchScalarField(ptf, p, iF, mapper)
55 {}
56 
57 
60 (
61  const fvPatch& p,
63  const dictionary& dict
64 )
65 :
66  fixedGradientFvPatchScalarField(p, iF, dict)
67 {}
68 
69 
72 (
74 )
75 :
76  fixedGradientFvPatchScalarField(tppsf)
77 {}
78 
79 
82 (
85 )
86 :
87  fixedGradientFvPatchScalarField(tppsf, iF)
88 {}
89 
90 
91 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
92 
94 {
95  if (updated())
96  {
97  return;
98  }
99 
100  const psiuReactionThermo& thermo = db().lookupObject<psiuReactionThermo>
101  (
103  );
104 
105  const label patchi = patch().index();
106 
107  const scalarField& pw = thermo.p().boundaryField()[patchi];
108  fvPatchScalarField& Tw =
109  const_cast<fvPatchScalarField&>(thermo.Tu().boundaryField()[patchi]);
110 
111  Tw.evaluate();
112 
113  gradient() = thermo.Cp(pw, Tw, patchi)*Tw.snGrad()
114  + patch().deltaCoeffs()*
115  (
116  thermo.heu(pw, Tw, patchi)
117  - thermo.heu(pw, Tw, patch().faceCells())
118  );
119 
120  fixedGradientFvPatchScalarField::updateCoeffs();
121 }
122 
123 
124 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
125 
126 namespace Foam
127 {
129  (
132  );
133 }
134 
135 // ************************************************************************* //
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
virtual tmp< Field< Type > > snGrad() const
Return patch-normal gradient.
Definition: fvPatchField.C:216
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
const Boundary & boundaryField() const
Return const-reference to the boundary field.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
virtual const volScalarField & Tu() const =0
Unburnt gas temperature [K].
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
rhoReactionThermo & thermo
Definition: createFields.H:28
Macros for easy insertion into run-time selection tables.
const word dictName() const
Return the local dictionary name (final part of scoped name)
Definition: dictionary.H:115
virtual volScalarField & p()
Pressure [Pa].
Definition: basicThermo.C:487
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Foam::fvPatchFieldMapper.
Foam::psiuReactionThermo.
gradientUnburntEnthalpyFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
virtual tmp< volScalarField > Cp() const =0
Heat capacity at constant pressure [J/kg/K].
virtual volScalarField & heu()=0
Unburnt gas enthalpy [J/kg].
virtual void evaluate(const Pstream::commsTypes commsType=Pstream::commsTypes::blocking)
Evaluate the patch field, sets Updated to false.
Definition: fvPatchField.C:331
label patchi
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Namespace for OpenFOAM.