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 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  scalarField& sumKappaTByDelta,
36  scalarField& sumKappaByDelta,
37  scalarField& Tref,
38  scalarField& Tw,
39  scalarField& sumq,
40  scalarField& qByKappa
41 ) const
42 {
43  // Lookup the fluid model
44  const phaseSystem& fluid =
45  db().lookupObject<phaseSystem>(phaseSystem::propertiesName);
46 
47  const phaseModel& thisPhase = fluid.phases()[internalField().group()];
48 
49  scalarField sumKappa(size(), scalar(0));
50  scalarField sumKappaT(size(), scalar(0));
51  scalarField sumKappaTw(size(), scalar(0));
52 
53  const label patchi = patch().index();
54 
55  forAll(fluid.phases(), phasei)
56  {
57  const phaseModel& phase = fluid.phases()[phasei];
58 
59  if (&phase != &thisPhase)
60  {
61  const scalarField& alpha = phase.boundaryField()[patchi];
62  const scalarField kappaEff(phase.kappaEff(patchi));
63  const scalarField alphaKappaEff(alpha*kappaEff);
64 
65  const fvPatchScalarField& T =
66  phase.thermo().T().boundaryField()[patchi];
67 
68  sumKappa += alphaKappaEff;
69  sumKappaT += alphaKappaEff*T.patchInternalField();
70  sumKappaTw += alphaKappaEff*T;
71  }
72  }
73 
74  const scalarField& alpha = thisPhase.boundaryField()[patchi];
75  const scalarField kappaEff(thisPhase.kappaEff(patchi));
76  const scalarField alphaKappaEff(alpha*kappaEff);
77 
78  kappa = alphaKappaEff;
79  qByKappa = sumq/(max(alpha, rootSmall)*kappaEff);
80  // sumq -= alpha*sumq;
81 
82  const scalarField& T = *this;
83  Tw = (sumKappaTw + alphaKappaEff*T)/(sumKappa + alphaKappaEff);
84 
85  Tref =
86  (sumKappaT + rootSmall*kappaEff*patchInternalField())
87  /(sumKappa + rootSmall*kappaEff);
88 
89  sumKappaByDelta = sumKappa*patch().deltaCoeffs();
90  sumKappaTByDelta = sumKappaT*patch().deltaCoeffs();
91 }
92 
93 
94 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
95 
98 (
99  const fvPatch& p,
101  const dictionary& dict
102 )
103 :
105 {}
106 
107 
110 (
112  const fvPatch& p,
114  const fieldMapper& mapper
115 )
116 :
117  externalTemperatureFvPatchScalarField(psf, p, iF, mapper)
118 {}
119 
120 
123 (
126 )
127 :
129 {}
130 
131 
132 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
133 
134 namespace Foam
135 {
137  (
140  );
141 }
142 
143 
144 // ************************************************************************* //
#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.
virtual const volScalarField & T() const =0
Temperature [K].
A list of keyword definitions, which are a keyword followed by any number of values (e....
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:88
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, scalarField &sumKappaTByDelta, scalarField &sumKappaByDelta, scalarField &Tref, scalarField &Tw, scalarField &sumq, scalarField &qByKappa) 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 and model interfacial transfers between them.
Definition: phaseSystem.H:73
static const word propertiesName
Default name of the phase properties dictionary.
Definition: phaseSystem.H:226
const phaseModelList & phases() const
Return the phase models.
Definition: phaseSystemI.H:102
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
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
dictionary dict
volScalarField & p