This repository has been archived on 2023-11-20. You can view files and clone it, but cannot push or open issues or pull requests.
foam-extend4.1-coherent-io/applications/solvers/multiphase/interMixingFoam/threePhaseInterfaceProperties/threePhaseInterfaceProperties.H

157 lines
4.4 KiB
C++
Raw Normal View History

2010-09-22 18:13:13 +00:00
/*---------------------------------------------------------------------------*\
========= |
2013-12-11 16:09:41 +00:00
\\ / F ield | foam-extend: Open Source CFD
2016-06-20 15:00:40 +00:00
\\ / O peration | Version: 4.0
\\ / A nd | Web: http://www.foam-extend.org
\\/ M anipulation | For copyright notice see file Copyright
2010-09-22 18:13:13 +00:00
-------------------------------------------------------------------------------
License
2013-12-11 16:09:41 +00:00
This file is part of foam-extend.
2010-09-22 18:13:13 +00:00
2013-12-11 16:09:41 +00:00
foam-extend is free software: you can redistribute it and/or modify it
2010-09-22 18:13:13 +00:00
under the terms of the GNU General Public License as published by the
2013-12-11 16:09:41 +00:00
Free Software Foundation, either version 3 of the License, or (at your
2010-09-22 18:13:13 +00:00
option) any later version.
2013-12-11 16:09:41 +00:00
foam-extend is distributed in the hope that it will be useful, but
WITHOUT ANY WARRANTY; without even the implied warranty of
MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
General Public License for more details.
2010-09-22 18:13:13 +00:00
You should have received a copy of the GNU General Public License
2013-12-11 16:09:41 +00:00
along with foam-extend. If not, see <http://www.gnu.org/licenses/>.
2010-09-22 18:13:13 +00:00
Class
threePhaseInterfaceProperties
Description
Properties to aid interFoam :
1. Correct the alpha boundary condition for dynamic contact angle.
2. Calculate interface curvature.
SourceFiles
threePhaseInterfaceProperties.C
\*---------------------------------------------------------------------------*/
#ifndef threePhaseInterfaceProperties_H
#define threePhaseInterfaceProperties_H
#include "threePhaseMixture.H"
#include "surfaceFields.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
/*---------------------------------------------------------------------------*\
Class threePhaseInterfaceProperties Declaration
\*---------------------------------------------------------------------------*/
class threePhaseInterfaceProperties
{
// Private data
const threePhaseMixture& mixture_;
//- Compression coefficient
scalar cAlpha_;
//- Surface tension 1-2
dimensionedScalar sigma12_;
//- Surface tension 1-3
dimensionedScalar sigma13_;
//- Stabilisation for normalisation of the interface normal
const dimensionedScalar deltaN_;
surfaceScalarField nHatf_;
volScalarField K_;
// Private Member Functions
//- Disallow default bitwise copy construct and assignment
threePhaseInterfaceProperties(const threePhaseInterfaceProperties&);
void operator=(const threePhaseInterfaceProperties&);
//- Correction for the boundary condition on the unit normal nHat on
// walls to produce the correct contact dynamic angle
// calculated from the component of U parallel to the wall
void correctContactAngle
(
surfaceVectorField::GeometricBoundaryField& nHat
) const;
//- Re-calculate the interface curvature
void calculateK();
public:
//- Conversion factor for degrees into radians
static const scalar convertToRad;
// Constructors
//- Construct from volume fraction field alpha and IOdictionary
threePhaseInterfaceProperties(const threePhaseMixture& mixture);
// Member Functions
scalar cAlpha() const
{
return cAlpha_;
}
const dimensionedScalar& deltaN() const
{
return deltaN_;
}
const surfaceScalarField& nHatf() const
{
return nHatf_;
}
const volScalarField& K() const
{
return K_;
}
tmp<volScalarField> sigma() const
{
volScalarField limitedAlpha2 = max(mixture_.alpha2(), scalar(0));
volScalarField limitedAlpha3 = max(mixture_.alpha3(), scalar(0));
return
(limitedAlpha2*sigma12_ + limitedAlpha3*sigma13_)
/(limitedAlpha2 + limitedAlpha3 + SMALL);
}
tmp<volScalarField> sigmaK() const
{
return sigma()*K_;
}
void correct()
{
calculateK();
}
};
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace Foam
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#endif
// ************************************************************************* //