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-2020 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::alphatPhaseJayatillekeWallFunctionFvPatchScalarField
32 
33 SourceFiles
34  alphatPhaseChangeWallFunctionFvPatchScalarField.C
35 
36 \*---------------------------------------------------------------------------*/
37 
38 #ifndef alphatPhaseChangeWallFunctionFvPatchScalarField_H
39 #define alphatPhaseChangeWallFunctionFvPatchScalarField_H
40 
42 #include "phasePairKey.H"
43 
44 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
45 
46 namespace Foam
47 {
48 namespace compressible
49 {
50 
51 /*---------------------------------------------------------------------------*\
52  Class alphatPhaseChangeWallFunctionFvPatchScalarField Declaration
53 \*---------------------------------------------------------------------------*/
54 
56 :
58 {
59 protected:
60 
61  // Protected data
62 
63  //- Name of the other phase
64  const word otherPhaseName_;
65 
66  //- Relaxation factor
67  const scalar relax_;
68 
69  //- Rate of phase-change
71 
72  //- Latent heat of the phase-change
74 
75 
76 public:
77 
78  //- Runtime type information
79  TypeName("compressible::alphatPhaseChangeWallFunction");
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  // alphatPhaseChangeWallFunctionFvPatchScalarField
101  // onto a new patch
103  (
105  const fvPatch&,
107  const fvPatchFieldMapper&
108  );
109 
110  //- Copy constructor
112  (
114  );
115 
116  //- Copy constructor setting internal field reference
118  (
121  );
122 
123 
124  // Member Functions
125 
126  //- Is there phase change mass transfer for this phasePair
127  bool activePhasePair(const phasePairKey&) const;
128 
129  //- Return the rate of phase-change
130  const scalarField& dmdtf() const;
131 
132  //- Return the rate of phase-change for specific phase pair
133  const scalarField& dmdtf(const phasePairKey&) const;
134 
135  //- Return the enthalpy source due to phase-change
136  const scalarField& dmdtLf() const;
137 
138  //- Return the rate of phase-change for specific phase pair
139  const scalarField& dmdtLf(const phasePairKey&) const;
140 
141 
142  // Mapping functions
143 
144  //- Map (and resize as needed) from self given a mapping object
145  // Used to update fields following mesh topology change
146  virtual void autoMap(const fvPatchFieldMapper&);
147 
148  //- Reverse map the given fvPatchField onto this fvPatchField
149  // Used to reconstruct fields
150  virtual void rmap(const fvPatchScalarField&, const labelList&);
151 
152 
153  // Evaluation functions
154 
155  //- Update the coefficients associated with the patch field
156  virtual void updateCoeffs() = 0;
157 
158 
159  // I-O
160 
161  //- Write
162  virtual void write(Ostream&) const;
163 };
164 
165 
166 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
167 
168 } // End namespace compressible
169 } // End namespace Foam
170 
171 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
172 
173 #endif
174 
175 // ************************************************************************* //
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
virtual void rmap(const fvPatchScalarField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
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:54
alphatPhaseChangeWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
const scalarField & dmdtLf() const
Return the enthalpy source due to phase-change.
TypeName("compressible::alphatPhaseChangeWallFunction")
Runtime type information.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
This boundary condition provides a thermal wall function for turbulent thermal diffusivity (usuallyal...
bool activePhasePair(const phasePairKey &) const
Is there phase change mass transfer for this phasePair.
const scalarField & dmdtf() const
Return the rate of phase-change.
Namespace for OpenFOAM.