IATEwallBoiling.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) 2016-2022 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 
26 #include "IATEwallBoiling.H"
28 #include "fvmSup.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace diameterModels
36 {
37 namespace IATEsources
38 {
39  defineTypeNameAndDebug(wallBoiling, 0);
40  addToRunTimeSelectionTable(IATEsource, wallBoiling, dictionary);
41 }
42 }
43 }
44 
45 
46 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
47 
49 (
50  const IATE& iate,
51  const dictionary& dict
52 )
53 :
54  IATEsource(iate)
55 {}
56 
57 
58 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
59 
62 (
63  const volScalarField& alphai,
64  volScalarField& kappai
65 ) const
66 {
68  (
69  IOobject
70  (
71  "wallBoiling:R",
72  phase().time().timeName(),
73  phase().mesh()
74  ),
75  phase().mesh(),
77  );
78 
80  (
81  IOobject
82  (
83  "wallBoiling:Rdk",
84  phase().time().timeName(),
85  phase().mesh()
86  ),
87  phase().mesh(),
88  dimensionedScalar(kappai.dimensions()/dimTime, 0)
89  );
90 
91  const volScalarField& alphat =
92  phase().mesh().lookupObject<volScalarField>
93  (
94  IOobject::groupName("alphat", otherPhase().name())
95  );
96 
97  const volScalarField::Boundary& alphatBf = alphat.boundaryField();
98 
99  const scalarField& rho = phase().rho();
100 
101  typedef compressible::alphatWallBoilingWallFunctionFvPatchScalarField
102  alphatWallBoilingWallFunction;
103 
104  forAll(alphatBf, patchi)
105  {
106  if (isA<alphatWallBoilingWallFunction>(alphatBf[patchi]))
107  {
108  const alphatWallBoilingWallFunction& alphatw =
109  refCast<const alphatWallBoilingWallFunction>(alphatBf[patchi]);
110 
111  const scalarField& dmdt = alphatw.dmdtf();
112  const scalarField& dDep = alphatw.dDeparture();
113 
114  const labelList& faceCells = alphatw.patch().faceCells();
115 
116  forAll(alphatw, facei)
117  {
118  if (dmdt[facei] > small)
119  {
120  const label faceCelli = faceCells[facei];
121  R[faceCelli] =
122  dmdt[facei]/(alphai[faceCelli]*rho[faceCelli]);
123  Rdk[faceCelli] = R[faceCelli]*(6.0/dDep[facei]);
124  }
125  }
126  }
127  }
128 
129  return Rdk - fvm::Sp(R, kappai);
130 }
131 
132 
133 // ************************************************************************* //
virtual tmp< fvScalarMatrix > R(const volScalarField &alphai, volScalarField &kappai) const
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
FvWallInfoData< WallInfo, label > label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
const phaseModel & otherPhase() const
Definition: IATEsource.H:148
virtual tmp< volScalarField > rho() const =0
Return the density field.
GeometricBoundaryField< scalar, fvPatchField, volMesh > Boundary
Type of the boundary field.
const dimensionSet dimless
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
wallBoiling(const IATE &iate, const dictionary &dict)
fvMesh & mesh
DimensionedField< scalar, volMesh > Internal
Type of the internal field from which this GeometricField is derived.
Macros for easy insertion into run-time selection tables.
virtual tmp< fvScalarMatrix > R(const volScalarField &alphai, volScalarField &kappai) const
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:58
const dimensionSet dimTime
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
static word groupName(Name name, const word &group)
const phaseModel & phase() const
Definition: IATEsource.H:138
word timeName
Definition: getTimeIndex.H:3
List< label > labelList
A List of labels.
Definition: labelList.H:56
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
const Mesh & mesh() const
Return mesh.
defineTypeNameAndDebug(combustionModel, 0)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
label patchi
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
A class for managing temporary objects.
Definition: PtrList.H:53
Calculate the matrix for implicit and explicit sources.
Namespace for OpenFOAM.