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-2024 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"
34 #include "fvcGrad.H"
35 
36 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
37 
40 (
41  const fvPatch& p,
43  const dictionary& dict,
44  const bool readValue
45 )
46 :
48  phiName_(dict.lookupOrDefault<word>("phi", "phi")),
49  UName_(dict.lookupOrDefault<word>("U", "U"))
50 {
51  if (readValue)
52  {
54  (
55  scalarField("value", iF.dimensions(), dict, p.size())
56  );
57  }
58 }
59 
60 
63 (
65  const fvPatch& p,
67  const fieldMapper& mapper
68 )
69 :
71  phiName_(ptf.phiName_),
72  UName_(ptf.UName_)
73 {}
74 
75 
78 (
81 )
82 :
84  phiName_(ptf.phiName_),
85  UName_(ptf.UName_)
86 {}
87 
88 
89 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
90 
93 {
96  const PtrList<volScalarField>& Y = thermo.Y();
97 
98  const fvPatchScalarField& Tp = *this;
99  const fvPatchScalarField& pp = thermo.p().boundaryField()[patch().index()];
100 
101  // Sum up the phiHep from all the species
102  tmp<scalarField> tPhiHep(new scalarField(this->size(), 0));
103  scalarField& phiHep = tPhiHep.ref();
104  forAll(Y, i)
105  {
106  const fvPatchScalarField& Yp = Y[i].boundaryField()[patch().index()];
107 
108  if (!isA<YBCType>(Yp))
109  {
111  << "The mass-fraction condition on patch " << patch().name()
112  << " is not of type " << YBCType::typeName << "."
113  << exit(FatalError);
114  }
115 
116  phiHep += refCast<const YBCType>(Yp).phiYp()*thermo.hei(i, pp, Tp);
117  }
118 
119  return tPhiHep;
120 }
121 
122 
124 {
125  if (updated())
126  {
127  return;
128  }
129 
130  // Get the fluxes
131  const scalarField& phip =
132  patch().lookupPatchField<surfaceScalarField, scalar>(phiName_);
133  const fvPatchVectorField& Up =
134  patch().lookupPatchField<volVectorField, vector>(UName_);
135  tmp<scalarField> uPhip =
136  refCast<const specieTransferVelocityFvPatchVectorField>(Up).phip();
137 
139  db().lookupType<fluidThermophysicalTransportModel>();
140 
141  // Get the current energy to linearise around
142  const fvPatchScalarField& hep =
143  ttm.thermo().he().boundaryField()[patch().index()];
144 
145  // Get the energy diffusivity
146  const scalarField AAlphaEffp
147  (
148  patch().magSf()*ttm.alphaEff(patch().index())
149  );
150 
151  // Compute the flux that we need to recover
152  const scalarField phiHep
153  (
154  this->phiHep()
155  - AAlphaEffp*hep.snGrad()
156  - patch().magSf()*ttm.q(patch().index())
157  );
158 
159  // Set the gradient and value so that the transport and diffusion combined
160  // result in the desired energy flux
161  heValueFraction() = phip/(phip - patch().deltaCoeffs()*AAlphaEffp);
162  heRefValue() = hep;
163  heRefGrad() = phip*(hep - phiHep/uPhip)/AAlphaEffp;
164 
165  mixedEnergyCalculatedTemperatureFvPatchScalarField::updateCoeffs();
166 }
167 
168 
170 (
171  Ostream& os
172 ) const
173 {
175  writeEntryIfDifferent<word>(os, "phi", "phi", phiName_);
176  writeEntryIfDifferent<word>(os, "U", "U", UName_);
177  writeEntry(os, "value", *this);
178 }
179 
180 
181 // * * * * * * * * * * * * * * Build Macro Function * * * * * * * * * * * * //
182 
183 namespace Foam
184 {
186  (
189  );
190 }
191 
192 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
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 dimensionSet & dimensions() const
Return dimensions.
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
virtual const volScalarField & he() const =0
Enthalpy/Internal energy [J/kg].
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
Abstract base class for field mapping.
Definition: fieldMapper.H:48
Base-class for multi-component fluid thermodynamic properties.
virtual const volScalarField & p() const =0
Pressure [Pa].
Abstract base class for fluid thermophysical transport models RAS, LES and laminar.
virtual tmp< scalarField > alphaEff(const label patchi) const =0
Effective thermal turbulent diffusivity.
virtual const fluidThermo & thermo() const =0
Access function to fluid thermophysical properties.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:91
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:237
friend Ostream & operator(Ostream &, const fvPatchField< Type > &)
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 ...
virtual scalar hei(const label speciei, const scalar p, const scalar T) const =0
Enthalpy/Internal energy [J/kg].
virtual PtrList< volScalarField > & Y()=0
Access the mass-fraction fields.
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< surfaceScalarField > q() const =0
Return the heat flux [W/m^2].
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:197
A class for handling words, derived from string.
Definition: word.H:62
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
Calculate the gradient of the given field.
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
VolField< vector > volVectorField
Definition: volFieldsFwd.H:63
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
volScalarField & p
PtrList< volScalarField > & Y
fluidMulticomponentThermo & thermo
Definition: createFields.H:31
Foam::surfaceFields.