noFilm.C
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) 2011-2021 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 \*---------------------------------------------------------------------------*/
25 
26 #include "noFilm.H"
28 #include "volFields.H"
29 
30 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace regionModels
35 {
36 namespace surfaceFilmModels
37 {
38 
39 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
40 
41 defineTypeNameAndDebug(noFilm, 0);
42 addToRunTimeSelectionTable(surfaceFilmModel, noFilm, mesh);
43 
44 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
45 
47 (
48  const word& modelType,
49  const fvMesh& mesh,
50  const dimensionedVector& g,
51  const word& regionType
52 )
53 :
55  mesh_(mesh)
56 {}
57 
58 
59 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
60 
62 {}
63 
64 
65 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
66 
67 Foam::scalar noFilm::CourantNumber() const
68 {
69  return 0;
70 }
71 
72 
74 {
76  (
77  "noFilm::Srho",
78  mesh_,
80  );
81 }
82 
83 
85 {
87  (
88  "noFilm::SY(" + Foam::name(i) + ")",
89  mesh_,
91  );
92 }
93 
94 
96 {
98  (
99  IOobject::modelName("SU", typeName),
100  mesh_,
102  );
103 }
104 
105 
107 {
109  (
110  "noFilm::Sh",
111  mesh_,
113  );
114 }
115 
116 
118 {}
119 
120 
121 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
122 
123 } // End namespace surfaceFilmModels
124 } // End namespace regionModels
125 } // End namespace Foam
126 
127 // ************************************************************************* //
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
static tmp< DimensionedField< Type, GeoMesh > > New(const word &name, const Mesh &mesh, const dimensionSet &)
Return a temporary field constructed from name, mesh.
Base class for surface film models.
virtual tmp< volScalarField::Internal > SYi(const label i) const
Return mass source for specie i - Eulerian phase only.
Definition: noFilm.C:84
virtual tmp< volScalarField::Internal > Srho() const
Return total mass source - Eulerian phase only.
Definition: noFilm.C:73
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
noFilm(const word &modelType, const fvMesh &mesh, const dimensionedVector &g, const word &regionType)
Construct from components.
Definition: noFilm.C:47
virtual void evolve()
Main driver routing to evolve the region - calls other evolves.
Definition: noFilm.C:117
Macros for easy insertion into run-time selection tables.
const dimensionSet dimTime
dynamicFvMesh & mesh
A class for handling words, derived from string.
Definition: word.H:59
static const zero Zero
Definition: zero.H:97
const dimensionSet dimEnergy
const dimensionSet dimMass
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
const dimensionSet dimVolume
virtual tmp< volVectorField::Internal > SU() const
Return momentum source - Eulerian phase only.
Definition: noFilm.C:95
static word modelName(Name name, const word &model)
Return the name of the object within the given model.
A class for managing temporary objects.
Definition: PtrList.H:53
const dimensionedVector & g
virtual scalar CourantNumber() const
Courant number evaluation.
Definition: noFilm.C:67
defineTypeNameAndDebug(kinematicSingleLayer, 0)
virtual tmp< volScalarField::Internal > Sh() const
Return enthalpy source - Eulerian phase only.
Definition: noFilm.C:106
addToRunTimeSelectionTable(surfaceFilmModel, noFilm, mesh)
Namespace for OpenFOAM.