noFilm.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) 2011-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::surfaceFilmModels::noFilm
28
29Description
30 Dummy surfaceFilmModel to allow solvers supporting film simulations to be
31 run without a film region.
32
33SourceFiles
34 noFilm.C
35
36\*---------------------------------------------------------------------------*/
37
38#ifndef noFilm_H
39#define noFilm_H
40
41#include "surfaceFilmModel.H"
42
43// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
44
45namespace Foam
46{
47namespace regionModels
48{
49namespace surfaceFilmModels
50{
51
52/*---------------------------------------------------------------------------*\
53 Class noFilm Declaration
54\*---------------------------------------------------------------------------*/
56class noFilm
57:
58 public surfaceFilmModel
59{
60 // Private member data
61
62 //- Reference to the mesh
63 const fvMesh& mesh_;
64
65
66 // Private member functions
67
68 //- No copy construct
69 noFilm(const noFilm&) = delete;
70
71 //- No copy assignment
72 void operator=(const noFilm&) = delete;
73
74
75public:
76
77 //- Runtime type information
78 TypeName("none");
79
80
81 // Constructors
82
83 //- Construct from components
84 noFilm
85 (
86 const word& modelType,
87 const fvMesh& mesh,
88 const dimensionedVector& g,
89 const word& regionType
90 );
91
92
93 //- Destructor
94 virtual ~noFilm();
95
96
97 // Member Functions
98
99 // Solution parameters
100
101 //- Courant number evaluation
102 virtual scalar CourantNumber() const;
103
104
105 // Primary region source fields
106
107 //- Return total mass source - Eulerian phase only
108 virtual tmp<volScalarField::Internal> Srho() const;
109
110 //- Return mass source for specie i - Eulerian phase only
112 (
113 const label i
114 ) const;
115
116 //- Return enthalpy source - Eulerian phase only
117 virtual tmp<volScalarField::Internal> Sh() const;
118
119
120 // Evolution
121
122 //- Main driver routing to evolve the region - calls other evolves
123 virtual void evolve();
124};
125
126
127// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
128
129} // End namespace surfaceFilmModels
130} // regionModels
131} // End namespace Foam
132
133// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
134
135#endif
136
137// ************************************************************************* //
const uniformDimensionedVectorField & g
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
Base class for surface film models.
Dummy surfaceFilmModel to allow solvers supporting film simulations to be run without a film region.
Definition: noFilm.H:58
TypeName("none")
Runtime type information.
virtual scalar CourantNumber() const
Courant number evaluation.
Definition: noFilm.C:69
virtual void evolve()
Main driver routing to evolve the region - calls other evolves.
Definition: noFilm.C:132
virtual tmp< volScalarField::Internal > Srho() const
Return total mass source - Eulerian phase only.
Definition: noFilm.C:75
virtual tmp< volScalarField::Internal > Sh() const
Return enthalpy source - Eulerian phase only.
Definition: noFilm.C:113
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.
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition: typeInfo.H:73