specieTransferTemperatureFvPatchScalarField.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) 2019-2023 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 
29 #include "volFields.H"
30 #include "surfaceFields.H"
32 #include "basicSpecieMixture.H"
34 
35 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
36 
39 (
40  const fvPatch& p,
42  const dictionary& dict,
43  const bool readValue
44 )
45 :
47  phiName_(dict.lookupOrDefault<word>("phi", "phi")),
48  UName_(dict.lookupOrDefault<word>("U", "U"))
49 {
50  if (readValue)
51  {
53  }
54 }
55 
56 
59 (
61  const fvPatch& p,
63  const fvPatchFieldMapper& mapper
64 )
65 :
67  phiName_(ptf.phiName_),
68  UName_(ptf.UName_)
69 {}
70 
71 
74 (
77 )
78 :
80  phiName_(ptf.phiName_),
81  UName_(ptf.UName_)
82 {}
83 
84 
85 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
86 
89 {
92  const PtrList<volScalarField>& Y = mixture.Y();
93 
94  // Get thermodynamic properties
95  const fluidThermo& thermo =
96  db().lookupObject<fluidThermo>(physicalProperties::typeName);
97  const fvPatchScalarField& Tp = *this;
98  const fvPatchScalarField& pp = thermo.p().boundaryField()[patch().index()];
99 
100  // Sum up the phiHep from all the species
101  tmp<scalarField> tPhiHep(new scalarField(this->size(), 0));
102  scalarField& phiHep = tPhiHep.ref();
103  forAll(Y, i)
104  {
105  const fvPatchScalarField& Yp = Y[i].boundaryField()[patch().index()];
106 
107  if (!isA<YBCType>(Yp))
108  {
110  << "The mass-fraction condition on patch " << patch().name()
111  << " is not of type " << YBCType::typeName << "."
112  << exit(FatalError);
113  }
114 
115  phiHep += refCast<const YBCType>(Yp).phiYp()*mixture.HE(i, pp, Tp);
116  }
117 
118  return tPhiHep;
119 }
120 
121 
123 {
124  if (updated())
125  {
126  return;
127  }
128 
129  // Get the fluxes
130  const scalarField& phip =
131  patch().lookupPatchField<surfaceScalarField, scalar>(phiName_);
132  const fvPatchVectorField& Up =
133  patch().lookupPatchField<volVectorField, vector>(UName_);
134  tmp<scalarField> uPhip =
135  refCast<const specieTransferVelocityFvPatchVectorField>(Up).phip();
136 
138  db().lookupType<fluidThermophysicalTransportModel>();
139 
140  // Get the diffusivity
141  const scalarField AAlphaEffp
142  (
143  patch().magSf()
144  *ttm.kappaEff(patch().index())
145  /ttm.thermo().Cp().boundaryField()[patch().index()]
146  );
147 
148  // Get the current energy to linearise around
149  const fluidThermo& thermo =
150  db().lookupObject<fluidThermo>(physicalProperties::typeName);
151  const scalarField& hep = thermo.he().boundaryField()[patch().index()];
152 
153  heValueFraction() = phip/(phip - patch().deltaCoeffs()*AAlphaEffp);
154  heRefValue() = hep;
155  heRefGrad() = phip*(hep - phiHep()/uPhip)/AAlphaEffp;
156 
157  mixedEnergyCalculatedTemperatureFvPatchScalarField::updateCoeffs();
158 }
159 
160 
162 (
163  Ostream& os
164 ) const
165 {
167  writeEntryIfDifferent<word>(os, "phi", "phi", phiName_);
168  writeEntryIfDifferent<word>(os, "U", "U", UName_);
169  writeEntry(os, "value", *this);
170 }
171 
172 
173 // * * * * * * * * * * * * * * Build Macro Function * * * * * * * * * * * * //
174 
175 namespace Foam
176 {
178  (
181  );
182 }
183 
184 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const Boundary & boundaryField() const
Return const-reference to the boundary field.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: PtrList.H:75
Specialisation of basicMixture for a mixture consisting of a number for molecular species.
virtual volScalarField & he()=0
Enthalpy/Internal energy [J/kg].
virtual const volScalarField & Cp() const =0
Heat capacity at constant pressure [J/kg/K].
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
Base-class for fluid thermodynamic properties.
Definition: fluidThermo.H:57
virtual volScalarField & p()=0
Pressure [Pa].
Abstract base class for fluid thermophysical transport models RAS, LES and laminar.
virtual const fluidThermo & thermo() const =0
Access function to fluid thermophysical properties.
Foam::fvPatchFieldMapper.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:87
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:231
virtual void operator==(const fvPatchField< Type > &)
Definition: fvPatchField.C:417
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
Base class for temperature boundary conditions in which the parameters of the mixed energy condition ...
Abstract base class for specie-transferring mass fraction boundary conditions.
This is a temperature boundary condition for a specie-transferring wall.
const tmp< scalarField > phiHep() const
Return the flux of energy.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
specieTransferTemperatureFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &, const bool readValue=true)
Construct from patch, internal field and dictionary.
virtual tmp< volScalarField > kappaEff() const =0
Effective thermal turbulent conductivity.
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
A class for handling words, derived from string.
Definition: word.H:62
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
VolField< vector > volVectorField
Definition: volFieldsFwd.H:62
SurfaceField< scalar > surfaceScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
error FatalError
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
dictionary dict
basicSpecieMixture & composition
volScalarField & p
PtrList< volScalarField > & Y
fluidMulticomponentThermo & thermo
Definition: createFields.H:31
Foam::surfaceFields.