incompressibleThreePhaseMixture.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 | Copyright (C) 2011-2015 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::incompressibleThreePhaseMixture
26 
27 Description
28 
29 SourceFiles
30  incompressibleThreePhaseMixture.C
31 
32 \*---------------------------------------------------------------------------*/
33 
34 #ifndef incompressibleThreePhaseMixture_H
35 #define incompressibleThreePhaseMixture_H
36 
38 #include "IOdictionary.H"
40 #include "dimensionedScalar.H"
41 #include "volFields.H"
42 
43 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
44 
45 namespace Foam
46 {
47 
48 /*---------------------------------------------------------------------------*\
49  Class incompressibleThreePhaseMixture Declaration
50 \*---------------------------------------------------------------------------*/
51 
53 :
54  public IOdictionary,
55  public transportModel
56 {
57  // Private data
58 
59  word phase1Name_;
60  word phase2Name_;
61  word phase3Name_;
62 
63  volScalarField alpha1_;
64  volScalarField alpha2_;
65  volScalarField alpha3_;
66 
67  const volVectorField& U_;
68  const surfaceScalarField& phi_;
69 
70  volScalarField nu_;
71 
72  autoPtr<viscosityModel> nuModel1_;
73  autoPtr<viscosityModel> nuModel2_;
74  autoPtr<viscosityModel> nuModel3_;
75 
76  dimensionedScalar rho1_;
77  dimensionedScalar rho2_;
78  dimensionedScalar rho3_;
79 
80 
81  // Private Member Functions
82 
83  //- Calculate and return the laminar viscosity
84  void calcNu();
85 
86 
87 public:
88 
89  // Constructors
90 
91  //- Construct from components
93  (
94  const volVectorField& U,
95  const surfaceScalarField& phi
96  );
97 
98 
99  //- Destructor
101  {}
102 
103 
104  // Member Functions
106  const word phase1Name() const
107  {
108  return phase1Name_;
109  }
111  const word phase2Name() const
112  {
113  return phase2Name_;
114  }
116  const word phase3Name() const
117  {
118  return phase3Name_;
119  }
121  const volScalarField& alpha1() const
122  {
123  return alpha1_;
124  }
127  {
128  return alpha1_;
129  }
131  const volScalarField& alpha2() const
132  {
133  return alpha2_;
134  }
137  {
138  return alpha2_;
139  }
141  const volScalarField& alpha3() const
142  {
143  return alpha3_;
144  }
147  {
148  return alpha3_;
149  }
150 
151  //- Return const-access to phase1 density
152  const dimensionedScalar& rho1() const
153  {
154  return rho1_;
155  }
156 
157  //- Return const-access to phase2 density
158  const dimensionedScalar& rho2() const
159  {
160  return rho2_;
161  };
162 
163  //- Return const-access to phase3 density
164  const dimensionedScalar& rho3() const
165  {
166  return rho3_;
167  };
168 
169  //- Return the velocity
170  const volVectorField& U() const
171  {
172  return U_;
173  }
174 
175  //- Return the flux
176  const surfaceScalarField& phi() const
177  {
178  return phi_;
179  }
180 
181  //- Return const-access to phase1 viscosityModel
182  const viscosityModel& nuModel1() const
183  {
184  return nuModel1_();
185  }
186 
187  //- Return const-access to phase2 viscosityModel
188  const viscosityModel& nuModel2() const
189  {
190  return nuModel2_();
191  }
192 
193  //- Return const-access to phase3 viscosityModel
194  const viscosityModel& nuModel3() const
195  {
196  return nuModel3_();
197  }
198 
199  //- Return the dynamic laminar viscosity
200  tmp<volScalarField> mu() const;
201 
202  //- Return the face-interpolated dynamic laminar viscosity
204 
205  //- Return the kinematic laminar viscosity
206  tmp<volScalarField> nu() const
207  {
208  return nu_;
209  }
210 
211  //- Return the laminar viscosity for patch
212  tmp<scalarField> nu(const label patchi) const
213  {
214  return nu_.boundaryField()[patchi];
215  }
216 
217  //- Return the face-interpolated dynamic laminar viscosity
219 
220  //- Correct the laminar viscosity
221  void correct()
222  {
223  calcNu();
224  }
225 
226  //- Read base transportProperties dictionary
227  bool read();
228 };
229 
230 
231 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
232 
233 } // End namespace Foam
234 
235 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
236 
237 #endif
238 
239 // ************************************************************************* //
const surfaceScalarField & phi() const
Return the flux.
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
incompressibleThreePhaseMixture(const volVectorField &U, const surfaceScalarField &phi)
Construct from components.
tmp< volScalarField > nu() const
Return the kinematic laminar viscosity.
tmp< volScalarField > mu() const
Return the dynamic laminar viscosity.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
const viscosityModel & nuModel2() const
Return const-access to phase2 viscosityModel.
const dimensionedScalar & rho2() const
Return const-access to phase2 density.
const viscosityModel & nuModel1() const
Return const-access to phase1 viscosityModel.
An abstract base class for incompressible viscosityModels.
const dimensionedScalar & rho3() const
Return const-access to phase3 density.
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:52
const dimensionedScalar & rho1() const
Return const-access to phase1 density.
const viscosityModel & nuModel3() const
Return const-access to phase3 viscosityModel.
A class for handling words, derived from string.
Definition: word.H:59
tmp< surfaceScalarField > muf() const
Return the face-interpolated dynamic laminar viscosity.
const volVectorField & U() const
Return the velocity.
label patchi
void correct()
Correct the laminar viscosity.
Base-class for all transport models used by the incompressible turbulence models. ...
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
A class for managing temporary objects.
Definition: PtrList.H:53
tmp< surfaceScalarField > nuf() const
Return the face-interpolated dynamic laminar viscosity.
bool read()
Read base transportProperties dictionary.
Namespace for OpenFOAM.