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-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 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 
43 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
44 
45 namespace Foam
46 {
47 
48 class phaseInterface;
49 
50 namespace compressible
51 {
52 
53 /*---------------------------------------------------------------------------*\
54  Class alphatPhaseChangeWallFunctionFvPatchScalarField Declaration
55 \*---------------------------------------------------------------------------*/
56 
58 :
60 {
61 protected:
62 
63  // Protected data
64 
65  //- Name of the other phase
66  const word otherPhaseName_;
67 
68  //- Relaxation factor
69  const scalar relax_;
70 
71  //- Rate of phase-change
73 
74 
75 public:
76 
77  //- Runtime type information
78  TypeName("compressible::alphatPhaseChangeWallFunction");
79 
80 
81  // Constructors
82 
83  //- Construct from patch and internal field
85  (
86  const fvPatch&,
88  );
89 
90  //- Construct from patch, internal field and dictionary
92  (
93  const fvPatch&,
95  const dictionary&
96  );
97 
98  //- Construct by mapping given
99  // alphatPhaseChangeWallFunctionFvPatchScalarField
100  // onto a new patch
102  (
104  const fvPatch&,
106  const fvPatchFieldMapper&
107  );
108 
109  //- Disallow copy without setting internal field reference
111  (
113  ) = delete;
114 
115  //- Copy constructor setting internal field reference
117  (
120  );
121 
122 
123  // Member Functions
124 
125  //- Is there phase change mass transfer for this interface?
126  bool activeInterface(const phaseInterface&) const;
127 
128  //- Return the rate of phase-change
129  const scalarField& dmdtf() const;
130 
131  //- Return the rate of phase-change for an interface
132  const scalarField& dmdtf(const phaseInterface&) const;
133 
134 
135  // Mapping functions
136 
137  //- Map (and resize as needed) from self given a mapping object
138  // Used to update fields following mesh topology change
139  virtual void autoMap(const fvPatchFieldMapper&);
140 
141  //- Reverse map the given fvPatchField onto this fvPatchField
142  // Used to reconstruct fields
143  virtual void rmap(const fvPatchScalarField&, const labelList&);
144 
145  //- Reset the fvPatchField to the given fvPatchField
146  // Used for mesh to mesh mapping
147  virtual void reset(const fvPatchScalarField&);
148 
149 
150  // Evaluation functions
151 
152  //- Update the coefficients associated with the patch field
153  virtual void updateCoeffs() = 0;
154 
155 
156  // I-O
157 
158  //- Write
159  virtual void write(Ostream&) const;
160 };
161 
162 
163 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
164 
165 } // End namespace compressible
166 } // End namespace Foam
167 
168 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
169 
170 #endif
171 
172 // ************************************************************************* //
virtual void reset(const fvPatchScalarField &)
Reset the fvPatchField to the given fvPatchField.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
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:63
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.
Class to represent an interface between phases. Derivations can further specify the configuration of ...
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.
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 activeInterface(const phaseInterface &) const
Is there phase change mass transfer for this interface?
const scalarField & dmdtf() const
Return the rate of phase-change.
Namespace for OpenFOAM.