mixedUnburntEnthalpyFvPatchScalarField.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  mixedFvPatchScalarField(p, iF)
42 {
43  valueFraction() = 0.0;
44  refValue() = 0.0;
45  refGrad() = 0.0;
46 }
47 
48 
51 (
53  const fvPatch& p,
55  const fvPatchFieldMapper& mapper
56 )
57 :
58  mixedFvPatchScalarField(ptf, p, iF, mapper)
59 {}
60 
61 
64 (
65  const fvPatch& p,
67  const dictionary& dict
68 )
69 :
70  mixedFvPatchScalarField(p, iF, dict)
71 {}
72 
73 
76 (
78 )
79 :
80  mixedFvPatchScalarField(tppsf)
81 {}
82 
83 
86 (
89 )
90 :
91  mixedFvPatchScalarField(tppsf, iF)
92 {}
93 
94 
95 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
96 
98 {
99  if (updated())
100  {
101  return;
102  }
103 
104  const psiuReactionThermo& thermo = db().lookupObject<psiuReactionThermo>
105  (
107  );
108 
109  const label patchi = patch().index();
110 
111  const scalarField& pw = thermo.p().boundaryField()[patchi];
112  mixedFvPatchScalarField& Tw = refCast<mixedFvPatchScalarField>
113  (
114  const_cast<fvPatchScalarField&>(thermo.Tu().boundaryField()[patchi])
115  );
116 
117  Tw.evaluate();
118 
119  valueFraction() = Tw.valueFraction();
120  refValue() = thermo.heu(pw, Tw.refValue(), patchi);
121  refGrad() = thermo.Cp(pw, Tw, patchi)*Tw.refGrad()
122  + patch().deltaCoeffs()*
123  (
124  thermo.heu(pw, Tw, patchi)
125  - thermo.heu(pw, Tw, patch().faceCells())
126  );
127 
128  mixedFvPatchScalarField::updateCoeffs();
129 }
130 
131 
132 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
133 
134 namespace Foam
135 {
137  (
140  );
141 }
142 
143 // ************************************************************************* //
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
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
Foam::fvPatchFieldMapper.
Foam::psiuReactionThermo.
virtual tmp< volScalarField > Cp() const =0
Heat capacity at constant pressure [J/kg/K].
virtual volScalarField & heu()=0
Unburnt gas enthalpy [J/kg].
mixedUnburntEnthalpyFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
label patchi
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Namespace for OpenFOAM.