alphatFixedDmdtWallBoilingWallFunctionFvPatchScalarField.H
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) 2015-2018 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 Class
25  Foam::compressible::
26  alphatFixedDmdtWallBoilingWallFunctionFvPatchScalarField
27 
28 Description
29  A simple alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField with
30  a fixed volumetric phase-change mass flux.
31 
32 See also
33  Foam::compressible::
34  alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField
35 
36 SourceFiles
37  alphatFixedDmdtWallBoilingWallFunctionFvPatchScalarField.C
38 
39 \*---------------------------------------------------------------------------*/
40 
41 #ifndef alphatFixedDmdtWallBoilingWallFunctionFvPatchScalarField_H
42 #define alphatFixedDmdtWallBoilingWallFunctionFvPatchScalarField_H
43 
45 
46 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
47 
48 namespace Foam
49 {
50 namespace compressible
51 {
52 
53 /*---------------------------------------------------------------------------*\
54  Class alphatFixedDmdtWallBoilingWallFunctionFvPatchScalarField Declaration
55 \*---------------------------------------------------------------------------*/
56 
58 :
60 {
61  // Private data
62 
63  //- name on the phase
64  word vaporPhaseName_;
65 
66  //- dmdt relaxationFactor
67  scalar relax_;
68 
69  //- Volumetric phase-change mass flux in near wall cells
70  scalar fixedDmdt_;
71 
72  //- Latent heat
73  scalar L_;
74 
75 
76 public:
77 
78  //- Runtime type information
79  TypeName("compressible::alphatFixedDmdtWallBoilingWallFunction");
80 
81 
82  // Constructors
83 
84  //- Construct from patch and internal field
86  (
87  const fvPatch&,
89  );
90 
91  //- Construct from patch, internal field and dictionary
93  (
94  const fvPatch&,
96  const dictionary&
97  );
98 
99  //- Construct by mapping given
100  // alphatFixedDmdtWallBoilingWallFunctionFvPatchScalarField
101  // onto a new patch
103  (
105  const fvPatch&,
107  const fvPatchFieldMapper&
108  );
109 
110  //- Construct as copy
112  (
114  );
116  //- Construct and return a clone
117  virtual tmp<fvPatchScalarField> clone() const
118  {
120  (
122  (
123  *this
124  )
125  );
126  }
127 
128  //- Construct as copy setting internal field reference
130  (
133  );
134 
135  //- Construct and return a clone setting internal field reference
137  (
139  ) const
140  {
142  (
144  (
145  *this,
146  iF
147  )
148  );
149  }
150 
151 
152  // Member functions
153 
154  //- Is there phase change mass transfer for this phasePair
155  virtual bool activePhasePair(const phasePairKey&) const;
156 
157  //- Return the rate of phase-change for specific phase pair
158  virtual const scalarField& dmdt(const phasePairKey&) const;
159 
160  //- Return the rate of phase-change for specific phase pair
161  virtual const scalarField& mDotL(const phasePairKey&) const;
162 
163  // Evaluation functions
164 
165  //- Update the coefficients associated with the patch field
166  virtual void updateCoeffs();
167 
168 
169  // I-O
170 
171  //- Write
172  virtual void write(Ostream&) const;
173 };
174 
175 
176 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
177 
178 } // End namespace compressible
179 } // End namespace Foam
180 
181 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
182 
183 #endif
184 
185 // ************************************************************************* //
virtual const scalarField & dmdt() const
Return the rate of phase-change.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
virtual const scalarField & mDotL() const
Return the enthalpy source due to phase-change.
TypeName("compressible::alphatFixedDmdtWallBoilingWallFunction")
Runtime type information.
virtual bool activePhasePair(const phasePairKey &) const
Is there phase change mass transfer for this phasePair.
A class for handling words, derived from string.
Definition: word.H:59
Foam::fvPatchFieldMapper.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
alphatFixedDmdtWallBoilingWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
A simple alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField with a fixed volumetric phase-cha...
This boundary condition provides a thermal wall function for turbulent thermal diffusivity (usuallyal...
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
A class for managing temporary objects.
Definition: PtrList.H:53
Namespace for OpenFOAM.