All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
incompressibleTwoPhaseMixture.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) 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 Class
25  Foam::incompressibleTwoPhaseMixture
26 
27 Description
28  A two-phase incompressible transportModel
29 
30 SourceFiles
31  incompressibleTwoPhaseMixture.C
32 
33 \*---------------------------------------------------------------------------*/
34 
35 #ifndef incompressibleTwoPhaseMixture_H
36 #define incompressibleTwoPhaseMixture_H
37 
39 #include "viscosityModel.H"
40 #include "twoPhaseMixture.H"
41 #include "IOdictionary.H"
42 
43 
44 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
45 
46 namespace Foam
47 {
48 
49 /*---------------------------------------------------------------------------*\
50  Class incompressibleTwoPhaseMixture Declaration
51 \*---------------------------------------------------------------------------*/
52 
54 :
55  public IOdictionary,
57  public twoPhaseMixture
58 {
59 protected:
60 
61  // Protected data
62 
65 
68 
70  const surfaceScalarField& phi_;
71 
73 
74 
75  // Protected Member Functions
76 
77  //- Calculate and return the laminar viscosity
78  void calcNu();
79 
80 
81 public:
82 
83  TypeName("incompressibleTwoPhaseMixture");
84 
85 
86  // Constructors
87 
88  //- Construct from components
90  (
91  const volVectorField& U,
92  const surfaceScalarField& phi
93  );
94 
95 
96  //- Destructor
98  {}
99 
100 
101  // Member Functions
102 
103  //- Return const-access to phase1 viscosityModel
104  const viscosityModel& nuModel1() const
105  {
106  return nuModel1_();
107  }
108 
109  //- Return const-access to phase2 viscosityModel
110  const viscosityModel& nuModel2() const
111  {
112  return nuModel2_();
113  }
114 
115  //- Return const-access to phase1 density
116  const dimensionedScalar& rho1() const
117  {
118  return rho1_;
119  }
120 
121  //- Return const-access to phase2 density
122  const dimensionedScalar& rho2() const
123  {
124  return rho2_;
125  };
126 
127  //- Return const-access to the mixture velocity
128  const volVectorField& U() const
129  {
130  return U_;
131  }
132 
133  //- Return the dynamic laminar viscosity
134  tmp<volScalarField> mu() const;
135 
136  //- Return the face-interpolated dynamic laminar viscosity
138 
139  //- Return the kinematic laminar viscosity
140  virtual tmp<volScalarField> nu() const
141  {
142  return nu_;
143  }
144 
145  //- Return the laminar viscosity for patch
146  virtual tmp<scalarField> nu(const label patchi) const
147  {
148  return nu_.boundaryField()[patchi];
149  }
150 
151  //- Return the face-interpolated kinematic laminar viscosity
153 
154  //- Correct the laminar viscosity
155  virtual void correct()
156  {
157  calcNu();
158  }
159 
160  //- Read base transportProperties dictionary
161  virtual bool read();
162 };
163 
164 
165 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
166 
167 } // End namespace Foam
168 
169 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
170 
171 #endif
172 
173 // ************************************************************************* //
const dimensionedScalar & rho2() const
Return const-access to phase2 density.
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
tmp< surfaceScalarField > nuf() const
Return the face-interpolated kinematic laminar viscosity.
A two-phase incompressible transportModel.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
virtual tmp< volScalarField > nu() const
Return the kinematic laminar viscosity.
incompressibleTwoPhaseMixture(const volVectorField &U, const surfaceScalarField &phi)
Construct from components.
Generic GeometricField class.
An abstract base class for incompressible viscosityModels.
tmp< surfaceScalarField > muf() const
Return the face-interpolated dynamic laminar viscosity.
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:53
const volVectorField & U() const
Return const-access to the mixture velocity.
Base-class for all transport models used by the incompressible turbulence models. ...
phi
Definition: correctPhi.H:3
tmp< volScalarField > mu() const
Return the dynamic laminar viscosity.
virtual void correct()
Correct the laminar viscosity.
const viscosityModel & nuModel1() const
Return const-access to phase1 viscosityModel.
label patchi
const viscosityModel & nuModel2() const
Return const-access to phase2 viscosityModel.
TypeName("incompressibleTwoPhaseMixture")
A two-phase mixture model.
virtual bool read()
Read base transportProperties dictionary.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
void calcNu()
Calculate and return the laminar viscosity.
A class for managing temporary objects.
Definition: PtrList.H:53
Namespace for OpenFOAM.
const dimensionedScalar & rho1() const
Return const-access to phase1 density.