surfaceFilmModel.H
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2017 OpenFOAM Foundation
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26Class
27 Foam::regionModels::surfaceFilmModel
28
29Description
30 Base class for surface film models
31
32SourceFiles
33 surfaceFilmModelI.H
34 surfaceFilmModel.C
35 surfaceFilmModelNew.C
36
37\*---------------------------------------------------------------------------*/
38
39#ifndef surfaceFilmModel_H
40#define surfaceFilmModel_H
41
43#include "volFields.H"
44
45// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
46
47namespace Foam
48{
49namespace regionModels
50{
51
52/*---------------------------------------------------------------------------*\
53 Class surfaceFilmModel Declaration
54\*---------------------------------------------------------------------------*/
57{
58 // Private Member Functions
59
60 //- No copy construct
61 surfaceFilmModel(const surfaceFilmModel&) = delete;
62
63 //- No copy assignment
64 void operator=(const surfaceFilmModel&) = delete;
65
66
67public:
68
69 //- Runtime type information
70 TypeName("surfaceFilmModel");
71
72 //- Reference temperature for enthalpy
73 static const dimensionedScalar Tref;
74
75
76 // Declare runtime constructor selection table
79 (
80 autoPtr,
82 mesh,
83 (
84 const word& modelType,
85 const fvMesh& mesh,
86 const dimensionedVector& g,
87 const word& regionType
88 ),
89 (modelType, mesh, g, regionType)
90 );
91
92
93 //- Constructor
95
96
97 // Selectors
98
99 //- Return a reference to the selected surface film model
101 (
102 const fvMesh& mesh,
103 const dimensionedVector& g,
104 const word& regionType = "surfaceFilm"
105 );
106
107
108 //- Destructor
109 virtual ~surfaceFilmModel();
110
111
112 // Member Functions
113
114 // Solution parameters
115
116 //- Courant number evaluation
117 virtual scalar CourantNumber() const = 0;
118
119
120 // Primary region source fields
121
122 //- Return total mass source - Eulerian phase only
123 virtual tmp<volScalarField::Internal> Srho() const = 0;
124
125 //- Return mass source for specie i - Eulerian phase only
127 (
128 const label i
129 ) const = 0;
130
131 //- Return enthalpy source - Eulerian phase only
132 virtual tmp<volScalarField::Internal> Sh() const = 0;
133
134
135 // Evolution
136
137 //- Main driver routing to evolve the region - calls other evolves
138 virtual void evolve() = 0;
139};
140
141
142// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
143
144} // End namespace regionModels
145} // End namespace Foam
146
147// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
148
149#endif
150
151// ************************************************************************* //
const uniformDimensionedVectorField & g
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
Base class for surface film models.
static autoPtr< surfaceFilmModel > New(const fvMesh &mesh, const dimensionedVector &g, const word &regionType="surfaceFilm")
Return a reference to the selected surface film model.
static const dimensionedScalar Tref
Reference temperature for enthalpy.
TypeName("surfaceFilmModel")
Runtime type information.
virtual scalar CourantNumber() const =0
Courant number evaluation.
declareRunTimeSelectionTable(autoPtr, surfaceFilmModel, mesh,(const word &modelType, const fvMesh &mesh, const dimensionedVector &g, const word &regionType),(modelType, mesh, g, regionType))
virtual void evolve()=0
Main driver routing to evolve the region - calls other evolves.
virtual tmp< volScalarField::Internal > Srho() const =0
Return total mass source - Eulerian phase only.
virtual tmp< volScalarField::Internal > Srho(const label i) const =0
Return mass source for specie i - Eulerian phase only.
virtual tmp< volScalarField::Internal > Sh() const =0
Return enthalpy source - Eulerian phase only.
A class for managing temporary objects.
Definition: tmp.H:65
A class for handling words, derived from Foam::string.
Definition: word.H:68
dynamicFvMesh & mesh
Namespace for OpenFOAM.
Macros to ease declaration of run-time selection tables.
#define declareRunTimeSelectionTable(ptrWrapper, baseType, argNames, argList, parList)
Declare a run-time selection (variables and adder classes)
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition: typeInfo.H:73