gradientEnergyFvPatchScalarField.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2012 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
30 #include "fvPatchFieldMapper.H"
31 #include "volFields.H"
32 #include "basicThermo.H"
33 
34 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
35 
38 (
39  const fvPatch& p,
41 )
42 :
43  fixedGradientFvPatchScalarField(p, iF)
44 {}
45 
46 
49 (
51  const fvPatch& p,
53  const fvPatchFieldMapper& mapper
54 )
55 :
56  fixedGradientFvPatchScalarField(ptf, p, iF, mapper)
57 {}
58 
59 
62 (
63  const fvPatch& p,
65  const dictionary& dict
66 )
67 :
68  fixedGradientFvPatchScalarField(p, iF, dict)
69 {}
70 
71 
74 (
76 )
77 :
78  fixedGradientFvPatchScalarField(tppsf)
79 {}
80 
81 
84 (
87 )
88 :
89  fixedGradientFvPatchScalarField(tppsf, iF)
90 {}
91 
92 
93 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
94 
96 {
97  if (updated())
98  {
99  return;
100  }
101 
102  const basicThermo& thermo = basicThermo::lookupThermo(*this);
103  const label patchi = patch().index();
104 
105  const scalarField& pw = thermo.p().boundaryField()[patchi];
106  fvPatchScalarField& Tw =
107  const_cast<fvPatchScalarField&>(thermo.T().boundaryField()[patchi]);
108 
109  Tw.evaluate();
110 
111  gradient() = thermo.Cpv(pw, Tw, patchi)*Tw.snGrad()
112  + patch().deltaCoeffs()*
113  (
114  thermo.he(pw, Tw, patchi)
115  - thermo.he(pw, Tw, patch().faceCells())
116  );
117 
118  fixedGradientFvPatchScalarField::updateCoeffs();
119 }
120 
121 
123 {
126 }
127 
128 
129 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
130 
131 namespace Foam
132 {
134  (
136  gradientEnergyFvPatchScalarField
137  );
138 }
139 
140 // ************************************************************************* //
dictionary dict
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:70
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
void writeValueEntry(Ostream &os) const
Write *this field as a "value" entry.
Definition: fvPatchField.H:375
Smooth ATC in cells next to a set of patches supplied by type.
Definition: faceCells.H:52
Macros for easy insertion into run-time selection tables.
This boundary condition provides a gradient condition for internal energy, where the gradient is calc...
psiReactionThermo & thermo
Definition: createFields.H:28
virtual void write(Ostream &) const
Write.
fvPatchField< scalar > fvPatchScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
A FieldMapper for finite-volume patch fields.
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
OBJstream os(runTime.globalPath()/outputName)
virtual void evaluate(const Pstream::commsTypes commsType=Pstream::commsTypes::blocking)
Evaluate the patch field, sets updated() to false.
Definition: fvPatchField.C:329
gradientEnergyFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: areaFieldsFwd.H:42
const std::string patch
OpenFOAM patch number as a std::string.
static const basicThermo & lookupThermo(const fvPatchScalarField &pf)
Definition: basicThermo.C:448
volScalarField & p
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Namespace for OpenFOAM.