alphatPhaseChangeWallFunctionFvPatchScalarField.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::alphatPhaseChangeWallFunctionFvPatchScalarField
26 
27 Description
28  Abstract base-class for all alphatWallFunctions supporting phase-change.
29 
30 See also
31  Foam::fixedValueFvPatchScalarField
32  Foam::alphatWallFunctionFvPatchScalarField
33 
34 SourceFiles
35  alphatPhaseChangeWallFunctionFvPatchScalarField.C
36 
37 \*---------------------------------------------------------------------------*/
38 
39 #ifndef alphatPhaseChangeWallFunctionFvPatchScalarField_H
40 #define alphatPhaseChangeWallFunctionFvPatchScalarField_H
41 
43 #include "phasePairKey.H"
44 
45 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
46 
47 namespace Foam
48 {
49 namespace compressible
50 {
51 
52 /*---------------------------------------------------------------------------*\
53  Class alphatPhaseChangeWallFunctionFvPatchScalarField Declaration
54 \*---------------------------------------------------------------------------*/
55 
57 :
58  public fixedValueFvPatchScalarField
59 {
60 protected:
61 
62  // Protected data
63 
64  //- Rate of phase-change
66 
67  //- Latent heat of the phase-change
69 
70 
71 public:
72 
73  //- Runtime type information
74  TypeName("compressible::alphatPhaseChangeWallFunction");
75 
76 
77  // Constructors
78 
79  //- Construct from patch and internal field
81  (
82  const fvPatch&,
84  );
85 
86  //- Construct from patch, internal field and dictionary
88  (
89  const fvPatch&,
91  const dictionary&
92  );
93 
94  //- Construct by mapping given
95  // alphatPhaseChangeWallFunctionFvPatchScalarField
96  // onto a new patch
98  (
100  const fvPatch&,
102  const fvPatchFieldMapper&
103  );
104 
105  //- Construct as copy
107  (
109  );
110 
111  //- Construct as copy setting internal field reference
113  (
116  );
117 
118 
119  // Member functions
120 
121  //- Return the rate of phase-change
122  virtual const scalarField& dmdt() const
123  {
124  return dmdt_;
125  }
126 
127  //- Return the enthalpy source due to phase-change
128  virtual const scalarField& mDotL() const
129  {
130  return mDotL_;
131  }
132 
133  //- Is there phase change mass transfer for this phasePair
134  virtual bool activePhasePair(const phasePairKey&) const
135  {
136  return false;
137  }
138 
139  //- Return the rate of phase-change for specific phase pair
140  virtual const scalarField& dmdt(const phasePairKey&) const
141  {
142  return dmdt_;
143  }
144 
145  //- Return the rate of phase-change for specific phase pair
146  virtual const scalarField& mDotL(const phasePairKey&) const
147  {
148  return mDotL_;
149  }
150 
151  //- Return the rate of phase-change for specific phase
152  virtual scalarField dmdt(const word&) const
153  {
154  return dmdt_;
155  }
156 
157  //- Return the enthalpy source due to phase-change for specific phase
158  virtual scalarField mDotL(const word&) const
159  {
160  return mDotL_;
161  }
162 
163  // Evaluation functions
164 
165  //- Update the coefficients associated with the patch field
166  virtual void updateCoeffs() = 0;
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.
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.
virtual bool activePhasePair(const phasePairKey &) const
Is there phase change mass transfer for this phasePair.
virtual void updateCoeffs()=0
Update the coefficients associated with the patch field.
A class for handling words, derived from string.
Definition: word.H:59
Foam::fvPatchFieldMapper.
Abstract base-class for all alphatWallFunctions supporting phase-change.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
alphatPhaseChangeWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
TypeName("compressible::alphatPhaseChangeWallFunction")
Runtime type information.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Namespace for OpenFOAM.