multiphaseExternalTemperatureFvPatchScalarField.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) 2023-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 "phaseSystem.H"
29 
30 // * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
31 
33 (
35  tmp<scalarField>& sumKappaTcByDelta,
36  tmp<scalarField>& sumKappaByDelta,
38  tmp<scalarField>& sumq
39 ) const
40 {
41  const phaseSystem& fluid =
42  db().lookupObject<phaseSystem>(phaseSystem::propertiesName);
43 
44  const label patchi = patch().index();
45 
46  const phaseModel& thisPhase = fluid.phases()[internalField().group()];
47 
48  // Get alpha and kappa for this phase
49  const scalarField& alpha = thisPhase.boundaryField()[patchi];
50  const scalarField alphaKappaEff(alpha*thisPhase.kappaEff(patchi));
51 
52  // Get alpha and kappa sums for all other phases
53  scalarField sumKappa(size(), scalar(0));
54  scalarField sumKappaTc(size(), scalar(0));
55  scalarField sumKappaT(size(), scalar(0));
56  forAll(fluid.phases(), phasei)
57  {
58  const phaseModel& phase = fluid.phases()[phasei];
59 
60  if (&phase != &thisPhase)
61  {
62  const scalarField& alpha = phase.boundaryField()[patchi];
63  const scalarField alphaKappaEff(alpha*phase.kappaEff(patchi));
64 
65  const fvPatchScalarField& T =
66  phase.thermo().T().boundaryField()[patchi];
67 
68  sumKappa += alphaKappaEff;
69  sumKappaTc += alphaKappaEff*T.patchInternalField();
70  sumKappaT += alphaKappaEff*T;
71  }
72  }
73 
74  kappa = alphaKappaEff;
75 
76  if (fluid.phases().size() > 1)
77  {
78  sumKappaByDelta = sumKappa*patch().deltaCoeffs();
79  sumKappaTcByDelta = sumKappaTc*patch().deltaCoeffs();
80  }
81 
82  T = (sumKappaT + alphaKappaEff*(*this))/(sumKappa + alphaKappaEff);
83 }
84 
85 
86 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
87 
90 (
91  const fvPatch& p,
93  const dictionary& dict
94 )
95 :
97 {}
98 
99 
102 (
104  const fvPatch& p,
106  const fieldMapper& mapper
107 )
108 :
109  externalTemperatureFvPatchScalarField(psf, p, iF, mapper)
110 {}
111 
112 
115 (
118 )
119 :
121 {}
122 
123 
124 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
125 
126 namespace Foam
127 {
129  (
132  );
133 }
134 
135 
136 // ************************************************************************* //
#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].
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
This boundary condition applies a heat flux condition to temperature on an external wall....
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 getKappa(scalarField &kappa, tmp< scalarField > &sumKappaTcByDelta, tmp< scalarField > &sumKappaByDelta, tmp< scalarField > &T, tmp< scalarField > &sumq) const
Get the patch kappa, sum kappa*Tc/delta, kappa/delta and.
multiphaseExternalTemperatureFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &)
Construct from patch, internal field and dictionary.
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
A class for managing temporary objects.
Definition: tmp.H:55
label patchi
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.
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
void T(LagrangianPatchField< Type > &f, const LagrangianPatchField< Type > &f1)
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
dictionary dict
volScalarField & p