fixedUnburntEnthalpyFvPatchScalarField.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-2019 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 
33 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
34 
37 (
38  const fvPatch& p,
40 )
41 :
42  fixedValueFvPatchScalarField(p, iF)
43 {}
44 
45 
48 (
50  const fvPatch& p,
52  const fvPatchFieldMapper& mapper
53 )
54 :
55  fixedValueFvPatchScalarField(ptf, p, iF, mapper)
56 {}
57 
58 
61 (
62  const fvPatch& p,
64  const dictionary& dict
65 )
66 :
67  fixedValueFvPatchScalarField(p, iF, dict)
68 {}
69 
70 
73 (
75 )
76 :
77  fixedValueFvPatchScalarField(tppsf)
78 {}
79 
80 
83 (
86 )
87 :
88  fixedValueFvPatchScalarField(tppsf, iF)
89 {}
90 
91 
92 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
93 
95 {
96  if (updated())
97  {
98  return;
99  }
100 
101  const psiuReactionThermo& thermo = db().lookupObject<psiuReactionThermo>
102  (
104  );
105 
106  const label patchi = patch().index();
107 
108  fvPatchScalarField& Tw =
109  const_cast<fvPatchScalarField&>(thermo.Tu().boundaryField()[patchi]);
110  Tw.evaluate();
111  operator==(thermo.heu(Tw, patchi));
112 
113  fixedValueFvPatchScalarField::updateCoeffs();
114 }
115 
116 
117 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
118 
119 namespace Foam
120 {
122  (
125  );
126 }
127 
128 // ************************************************************************* //
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:158
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:123
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
Foam::fvPatchFieldMapper.
Foam::psiuReactionThermo.
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:249
label patchi
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
fixedUnburntEnthalpyFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Namespace for OpenFOAM.