Skip to content
Snippets Groups Projects
Commit 7c4e3a07 authored by Henry's avatar Henry
Browse files

threePhaseInterfaceProperties: Update constructor to be consistent with the...

threePhaseInterfaceProperties: Update constructor to be consistent with the two-phase interfaceProperties
Avoids problem of duplicate registration of K
parent 4b3c77ca
Branches
Tags
No related merge requests found
......@@ -2,7 +2,7 @@
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2011-2014 OpenFOAM Foundation
\\ / A nd | Copyright (C) 2011-2015 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
......@@ -21,14 +21,6 @@ License
You should have received a copy of the GNU General Public License
along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
Application
threePhaseInterfaceProperties
Description
Properties to aid interFoam :
1. Correct the alpha boundary condition for dynamic contact angle.
2. Calculate interface curvature.
\*---------------------------------------------------------------------------*/
#include "threePhaseInterfaceProperties.H"
......@@ -47,12 +39,6 @@ const Foam::scalar Foam::threePhaseInterfaceProperties::convertToRad =
// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
// Correction for the boundary condition on the unit normal nHat on
// walls to produce the correct contact angle.
// The dynamic contact angle is calculated from the component of the
// velocity on the direction of the interface, parallel to the wall.
void Foam::threePhaseInterfaceProperties::correctContactAngle
(
surfaceVectorField::GeometricBoundaryField& nHatb
......@@ -147,6 +133,7 @@ void Foam::threePhaseInterfaceProperties::calculateK()
// Face unit interface normal
surfaceVectorField nHatfv(gradAlphaf/(mag(gradAlphaf) + deltaN_));
correctContactAngle(nHatfv.boundaryField());
// Face unit interface normal flux
......@@ -157,9 +144,9 @@ void Foam::threePhaseInterfaceProperties::calculateK()
// Complex expression for curvature.
// Correction is formally zero but numerically non-zero.
//volVectorField nHat = gradAlpha/(mag(gradAlpha) + deltaN_);
//nHat.boundaryField() = nHatfv.boundaryField();
//K_ = -fvc::div(nHatf_) + (nHat & fvc::grad(nHatfv) & nHat);
// volVectorField nHat = gradAlpha/(mag(gradAlpha) + deltaN_);
// nHat.boundaryField() = nHatfv.boundaryField();
// K_ = -fvc::div(nHatf_) + (nHat & fvc::grad(nHatfv) & nHat);
}
......@@ -192,21 +179,26 @@ Foam::threePhaseInterfaceProperties::threePhaseInterfaceProperties
nHatf_
(
IOobject
(
fvc::interpolate(fvc::grad(mixture.alpha1()))
/(mag(fvc::interpolate(fvc::grad(mixture.alpha1()))) + deltaN_)
) & mixture.alpha1().mesh().Sf()
"nHatf",
mixture.alpha1().time().timeName(),
mixture.alpha1().mesh()
),
mixture.alpha1().mesh(),
dimensionedScalar("nHatf", dimArea, 0.0)
),
K_
(
IOobject
(
"K",
"interfaceProperties:K",
mixture.alpha1().time().timeName(),
mixture.alpha1().mesh()
),
-fvc::div(nHatf_)
mixture.alpha1().mesh(),
dimensionedScalar("K", dimless/dimLength, 0.0)
)
{
calculateK();
......@@ -233,5 +225,4 @@ Foam::threePhaseInterfaceProperties::nearInterface() const
}
// ************************************************************************* //
......@@ -78,8 +78,8 @@ class 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
// walls to produce the correct contact dynamic angle.
// Calculated from the component of U parallel to the wall
void correctContactAngle
(
surfaceVectorField::GeometricBoundaryField& nHat
......
0% or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment