FreeFOAM The Cross-Platform CFD Toolkit
Hosted by SourceForge:
Get FreeFOAM at SourceForge.net.
            Fast, secure and Free Open Source software downloads

threePhaseInterfaceProperties.C

Go to the documentation of this file.
00001 /*---------------------------------------------------------------------------*\
00002   =========                 |
00003   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
00004    \\    /   O peration     |
00005     \\  /    A nd           | Copyright (C) 1991-2010 OpenCFD Ltd.
00006      \\/     M anipulation  |
00007 -------------------------------------------------------------------------------
00008 License
00009     This file is part of OpenFOAM.
00010 
00011     OpenFOAM is free software: you can redistribute it and/or modify it
00012     under the terms of the GNU General Public License as published by
00013     the Free Software Foundation, either version 3 of the License, or
00014     (at your option) any later version.
00015 
00016     OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
00017     ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
00018     FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
00019     for more details.
00020 
00021     You should have received a copy of the GNU General Public License
00022     along with OpenFOAM.  If not, see <http://www.gnu.org/licenses/>.
00023     Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
00024 
00025 Application
00026     threePhaseInterfaceProperties
00027 
00028 Description
00029     Properties to aid interFoam :
00030     1. Correct the alpha boundary condition for dynamic contact angle.
00031     2. Calculate interface curvature.
00032 
00033 \*---------------------------------------------------------------------------*/
00034 
00035 #include "threePhaseInterfaceProperties.H"
00036 #include <twoPhaseInterfaceProperties/alphaContactAngleFvPatchScalarField.H>
00037 #include <OpenFOAM/mathematicalConstants.H>
00038 #include <finiteVolume/surfaceInterpolate.H>
00039 #include <finiteVolume/fvcDiv.H>
00040 #include <finiteVolume/fvcGrad.H>
00041 
00042 // * * * * * * * * * * * * * * * Static Member Data  * * * * * * * * * * * * //
00043 
00044 const Foam::scalar Foam::threePhaseInterfaceProperties::convertToRad =
00045     Foam::mathematicalConstant::pi/180.0;
00046 
00047 
00048 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
00049 
00050 // Correction for the boundary condition on the unit normal nHat on
00051 // walls to produce the correct contact angle.
00052 
00053 // The dynamic contact angle is calculated from the component of the
00054 // velocity on the direction of the interface, parallel to the wall.
00055 
00056 void Foam::threePhaseInterfaceProperties::correctContactAngle
00057 (
00058     surfaceVectorField::GeometricBoundaryField& nHatb
00059 ) const
00060 {
00061     const volScalarField::GeometricBoundaryField& alpha1 =
00062         mixture_.alpha1().boundaryField();
00063     const volScalarField::GeometricBoundaryField& alpha2 =
00064         mixture_.alpha2().boundaryField();
00065     const volScalarField::GeometricBoundaryField& alpha3 =
00066         mixture_.alpha3().boundaryField();
00067     const volVectorField::GeometricBoundaryField& U =
00068         mixture_.U().boundaryField();
00069 
00070     const fvMesh& mesh = mixture_.U().mesh();
00071     const fvBoundaryMesh& boundary = mesh.boundary();
00072 
00073     forAll(boundary, patchi)
00074     {
00075         if (isA<alphaContactAngleFvPatchScalarField>(alpha1[patchi]))
00076         {
00077             const alphaContactAngleFvPatchScalarField& a2cap =
00078                 refCast<const alphaContactAngleFvPatchScalarField>
00079                 (alpha2[patchi]);
00080 
00081             const alphaContactAngleFvPatchScalarField& a3cap =
00082                 refCast<const alphaContactAngleFvPatchScalarField>
00083                 (alpha3[patchi]);
00084 
00085             scalarField twoPhaseAlpha2 = max(a2cap, scalar(0));
00086             scalarField twoPhaseAlpha3 = max(a3cap, scalar(0));
00087 
00088             scalarField sumTwoPhaseAlpha =
00089                 twoPhaseAlpha2 + twoPhaseAlpha3 + SMALL;
00090 
00091             twoPhaseAlpha2 /= sumTwoPhaseAlpha;
00092             twoPhaseAlpha3 /= sumTwoPhaseAlpha;
00093 
00094             fvsPatchVectorField& nHatp = nHatb[patchi];
00095 
00096             scalarField theta =
00097                 convertToRad
00098                *(
00099                    twoPhaseAlpha2*(180 - a2cap.theta(U[patchi], nHatp))
00100                  + twoPhaseAlpha3*(180 - a3cap.theta(U[patchi], nHatp))
00101                );
00102 
00103             vectorField nf = boundary[patchi].nf();
00104 
00105             // Reset nHatPatch to correspond to the contact angle
00106 
00107             scalarField a12 = nHatp & nf;
00108 
00109             scalarField b1 = cos(theta);
00110 
00111             scalarField b2(nHatp.size());
00112 
00113             forAll(b2, facei)
00114             {
00115                 b2[facei] = cos(acos(a12[facei]) - theta[facei]);
00116             }
00117 
00118             scalarField det = 1.0 - a12*a12;
00119 
00120             scalarField a = (b1 - a12*b2)/det;
00121             scalarField b = (b2 - a12*b1)/det;
00122 
00123             nHatp = a*nf + b*nHatp;
00124 
00125             nHatp /= (mag(nHatp) + deltaN_.value());
00126         }
00127     }
00128 }
00129 
00130 
00131 void Foam::threePhaseInterfaceProperties::calculateK()
00132 {
00133     const volScalarField& alpha1 = mixture_.alpha1();
00134 
00135     const fvMesh& mesh = alpha1.mesh();
00136     const surfaceVectorField& Sf = mesh.Sf();
00137 
00138     // Cell gradient of alpha
00139     volVectorField gradAlpha = fvc::grad(alpha1);
00140 
00141     // Interpolated face-gradient of alpha
00142     surfaceVectorField gradAlphaf = fvc::interpolate(gradAlpha);
00143 
00144     // Face unit interface normal
00145     surfaceVectorField nHatfv = gradAlphaf/(mag(gradAlphaf) + deltaN_);
00146     correctContactAngle(nHatfv.boundaryField());
00147 
00148     // Face unit interface normal flux
00149     nHatf_ = nHatfv & Sf;
00150 
00151     // Simple expression for curvature
00152     K_ = -fvc::div(nHatf_);
00153 
00154     // Complex expression for curvature.
00155     // Correction is formally zero but numerically non-zero.
00156     //volVectorField nHat = gradAlpha/(mag(gradAlpha) + deltaN_);
00157     //nHat.boundaryField() = nHatfv.boundaryField();
00158     //K_ = -fvc::div(nHatf_) + (nHat & fvc::grad(nHatfv) & nHat);
00159 }
00160 
00161 
00162 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
00163 
00164 Foam::threePhaseInterfaceProperties::threePhaseInterfaceProperties
00165 (
00166     const threePhaseMixture& mixture
00167 )
00168 :
00169     mixture_(mixture),
00170     cAlpha_
00171     (
00172         readScalar
00173         (
00174             mixture.U().mesh().solutionDict().subDict("PISO").lookup("cAlpha")
00175         )
00176     ),
00177     sigma12_(mixture.lookup("sigma12")),
00178     sigma13_(mixture.lookup("sigma13")),
00179 
00180     deltaN_
00181     (
00182         "deltaN",
00183         1e-8/pow(average(mixture.U().mesh().V()), 1.0/3.0)
00184     ),
00185 
00186     nHatf_
00187     (
00188         (
00189             fvc::interpolate(fvc::grad(mixture.alpha1()))
00190            /(mag(fvc::interpolate(fvc::grad(mixture.alpha1()))) + deltaN_)
00191         ) & mixture.alpha1().mesh().Sf()
00192     ),
00193 
00194     K_
00195     (
00196         IOobject
00197         (
00198             "Kcond",
00199             mixture.alpha1().time().timeName(),
00200             mixture.alpha1().mesh()
00201         ),
00202         -fvc::div(nHatf_)
00203     )
00204 {
00205     calculateK();
00206 }
00207 
00208 
00209 // ************************************************************************* //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines