multiphaseCoupledTemperatureFvPatchScalarField.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) 2022-2025 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 
27 #include "fieldMapper.H"
28 #include "phaseSystem.H"
32 
33 // * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
34 
36 (
38  tmp<scalarField>& sumKappaTByDelta,
39  tmp<scalarField>& sumKappaByDelta,
40  scalarField& sumq,
41  tmp<scalarField>& qByKappa
42 ) const
43 {
44  // Lookup the fluid model
45  const phaseSystem& fluid =
46  db().lookupObject<phaseSystem>(phaseSystem::propertiesName);
47 
48  scalarField sumKappa(size(), scalar(0));
49  scalarField sumKappaT(size(), scalar(0));
50 
51  forAll(fluid.phases(), phasei)
52  {
53  const phaseModel& phase = fluid.phases()[phasei];
54  const rhoThermo& thermo = phase.thermo();
55 
56  const fvPatchScalarField& Tw =
57  thermo.T().boundaryField()[patch().index()];
58 
59  const fvPatchScalarField& alpha =
60  phase.boundaryField()[patch().index()];
61 
62  tmp<scalarField> kappaEff(phase.kappaEff(patch().index()));
63  tmp<scalarField> alphaKappaEff(alpha*kappaEff());
64 
65  if (&Tw == this)
66  {
67  kappa = alphaKappaEff;
68  qByKappa = sumq/kappaEff;
69  sumq -= alpha*sumq;
70  }
71  else
72  {
73  const scalarField T
74  (
75  thermo.T().boundaryField()[patch().index()]
76  .patchInternalField()
77  );
78 
79  sumKappa += alphaKappaEff();
80  sumKappaT += alphaKappaEff*T;
81  }
82  }
83 
84  sumKappaByDelta = sumKappa*patch().deltaCoeffs();
85  sumKappaTByDelta = sumKappaT*patch().deltaCoeffs();
86 }
87 
88 
90 (
91  tmp<scalarField>& sumKappaTByDeltaNbr,
92  tmp<scalarField>& sumKappaByDeltaNbr,
93  tmp<scalarField>& qNbr
94 ) const
95 {
96  // Lookup the fluid model
97  const phaseSystem& fluid =
98  patch().boundaryMesh().mesh()
100 
101  scalarField sumKappa(size(), scalar(0));
102  scalarField sumKappaT(size(), scalar(0));
103 
104  forAll(fluid.phases(), phasei)
105  {
106  const phaseModel& phase = fluid.phases()[phasei];
107  const rhoThermo& thermo = phase.thermo();
108 
109  const fvPatchScalarField& alpha =
110  phase.boundaryField()[patch().index()];
111 
112  const scalarField T
113  (
114  thermo.T().boundaryField()[patch().index()].patchInternalField()
115  );
116 
117  const scalarField alphaKappaEff(alpha*phase.kappaEff(patch().index()));
118 
119  sumKappa += alphaKappaEff;
120  sumKappaT += alphaKappaEff*T;
121  }
122 
123  sumKappaByDeltaNbr = sumKappa*patch().deltaCoeffs();
124  sumKappaTByDeltaNbr = sumKappaT*patch().deltaCoeffs();
125 }
126 
127 
129 (
130  tmp<scalarField>& TrefNbr,
131  tmp<scalarField>& qNbr
132 ) const
133 {
134  // Lookup the fluid model
135  const phaseSystem& fluid =
136  patch().boundaryMesh().mesh()
138 
139  TrefNbr = new scalarField(size(), scalar(0));
140  scalarField& Tw = TrefNbr.ref();
141 
142  forAll(fluid.phases(), phasei)
143  {
144  const phaseModel& phase = fluid.phases()[phasei];
145  const rhoThermo& thermo = phase.thermo();
146 
147  const fvPatchScalarField& alpha =
148  phase.boundaryField()[patch().index()];
149 
150  Tw += alpha*thermo.T().boundaryField()[patch().index()];
151 
152  /*
153  // Pending the addition of phase anisotropic thermal transport
154  tmp<scalarField> qCorr = ttm.qCorr(patch().index());
155 
156  if (qCorr.valid())
157  {
158  if (qCorr.valid())
159  {
160  qNbr.ref() += alpha*qCorr;
161  }
162  else
163  {
164  qNbr = alpha*qCorr;
165  }
166  }
167  */
168  }
169 }
170 
171 
172 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
173 
176 (
177  const fvPatch& p,
179  const dictionary& dict
180 )
181 :
183 {}
184 
185 
188 (
190  const fvPatch& p,
192  const fieldMapper& mapper
193 )
194 :
195  coupledTemperatureFvPatchScalarField(psf, p, iF, mapper)
196 {}
197 
198 
201 (
204 )
205 :
207 {}
208 
209 
210 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
211 
212 namespace Foam
213 {
215  (
218  );
219 }
220 
221 
222 // ************************************************************************* //
#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 Boundary & boundaryField() const
Return const-reference to the boundary field.
virtual const volScalarField & T() const =0
Temperature [K].
Mixed boundary condition for temperature, to be used for heat-transfer with another region in a CHT c...
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
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:91
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
Mixed boundary condition for the phase temperature of a phase in an Euler-Euler multiphase simulation...
virtual void getNbr(tmp< scalarField > &sumKappaTByDeltaNbr, tmp< scalarField > &sumKappaByDeltaNbr, tmp< scalarField > &qNbr) const
Get the neighbour patch sum kappa*Tc/delta and kappa/delta.
virtual void getThis(tmp< scalarField > &kappa, tmp< scalarField > &sumKappaTByDelta, tmp< scalarField > &sumKappaByDelta, scalarField &sumq, tmp< scalarField > &qByKappa) const
Get the patch kappa, sum kappa*Tc/delta and kappa/delta for all.
multiphaseCoupledTemperatureFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &)
Construct from patch, internal field and dictionary.
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type and name.
virtual tmp< scalarField > kappaEff(const label patchi) const =0
Effective thermal turbulent conductivity.
virtual const rhoThermo & thermo() const =0
Return the thermophysical model.
Class to represent a system of phases.
Definition: phaseSystem.H:74
static const word propertiesName
Default name of the phase properties dictionary.
Definition: phaseSystem.H:252
const phaseModelList & phases() const
Return the phase models.
Definition: phaseSystemI.H:104
const fvMesh & mesh() const
Return the mesh.
Definition: phaseSystemI.H:91
Base-class for thermodynamic properties based on density.
Definition: rhoThermo.H:54
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
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
Namespace for OpenFOAM.
void T(LagrangianPatchField< Type > &f, const LagrangianPatchField< Type > &f1)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
dictionary dict
volScalarField & p
fluidMulticomponentThermo & thermo
Definition: createFields.H:31