incompressibleTwoPhaseVoFMixture.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-2025 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 
27 #include "surfaceInterpolate.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
35 }
36 
37 
38 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
39 
41 (
42  const fvMesh& mesh
43 )
44 :
46 
47  nuModel1_(viscosityModel::New(mesh, phase1Name())),
48  nuModel2_(viscosityModel::New(mesh, phase2Name())),
49 
50  rho1_("rho", dimDensity, nuModel1_()),
51  rho2_("rho", dimDensity, nuModel2_()),
52 
53  rho_
54  (
55  IOobject
56  (
57  "rho",
58  mesh.time().name(),
59  mesh,
60  IOobject::NO_READ,
61  IOobject::AUTO_WRITE
62  ),
63  mesh,
64  dimensionedScalar("rho", dimDensity, 0)
65  ),
66 
67  nu_
68  (
69  IOobject
70  (
71  "nu",
72  mesh.time().name(),
73  mesh
74  ),
75  mesh,
77  calculatedFvPatchScalarField::typeName
78  )
79 {
80  correct();
81 }
82 
83 
84 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
85 
87 {
89  {
90  nuModel1_->lookup("rho") >> rho1_;
91  nuModel2_->lookup("rho") >> rho2_;
92 
93  return true;
94  }
95  else
96  {
97  return false;
98  }
99 }
100 
101 
103 {
104  rho_ = alpha1()*rho1_ + alpha2()*rho2_;
105 
106  nuModel1_->correct();
107  nuModel2_->correct();
108 
109  const volScalarField limitedAlpha1
110  (
111  "limitedAlpha1",
112  min(max(alpha1(), scalar(0)), scalar(1))
113  );
114 
115  // Average kinematic viscosity calculated from dynamic viscosity
116  nu_ =
117  (
118  limitedAlpha1*rho1_*nuModel1_->nu()
119  + (scalar(1) - limitedAlpha1)*rho2_*nuModel2_->nu()
120  )/(limitedAlpha1*rho1_ + (scalar(1) - limitedAlpha1)*rho2_);
121 }
122 
123 
124 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Generic GeometricField class.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:96
Class to represent a mixture of two constant density phases.
virtual void correct()
Correct the mixture density and laminar viscosity.
incompressibleTwoPhaseVoFMixture(const fvMesh &mesh)
Construct from a mesh.
virtual bool read()
Read base phaseProperties dictionary.
Class to represent a VoF mixture.
virtual bool read()=0
Read base phaseProperties dictionary.
An abstract base class for Newtonian viscosity models.
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
Namespace for OpenFOAM.
const dimensionSet dimKinematicViscosity
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
const dimensionSet dimDensity
defineTypeNameAndDebug(combustionModel, 0)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.