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

interfaceProperties.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-2011 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 
00024 \*---------------------------------------------------------------------------*/
00025 
00026 #include "interfaceProperties.H"
00027 #include <twoPhaseInterfaceProperties/alphaContactAngleFvPatchScalarField.H>
00028 #include <OpenFOAM/mathematicalConstants.H>
00029 #include <finiteVolume/surfaceInterpolate.H>
00030 #include <finiteVolume/fvcDiv.H>
00031 #include <finiteVolume/fvcGrad.H>
00032 #include <finiteVolume/fvcSnGrad.H>
00033 
00034 // * * * * * * * * * * * * * * * Static Member Data  * * * * * * * * * * * * //
00035 
00036 const Foam::scalar Foam::interfaceProperties::convertToRad =
00037     Foam::mathematicalConstant::pi/180.0;
00038 
00039 
00040 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
00041 
00042 // Correction for the boundary condition on the unit normal nHat on
00043 // walls to produce the correct contact angle.
00044 
00045 // The dynamic contact angle is calculated from the component of the
00046 // velocity on the direction of the interface, parallel to the wall.
00047 
00048 void Foam::interfaceProperties::correctContactAngle
00049 (
00050     surfaceVectorField::GeometricBoundaryField& nHatb,
00051     surfaceVectorField::GeometricBoundaryField& gradAlphaf
00052 ) const
00053 {
00054     const fvMesh& mesh = alpha1_.mesh();
00055     const volScalarField::GeometricBoundaryField& abf = alpha1_.boundaryField();
00056 
00057     const fvBoundaryMesh& boundary = mesh.boundary();
00058 
00059     forAll(boundary, patchi)
00060     {
00061         if (isA<alphaContactAngleFvPatchScalarField>(abf[patchi]))
00062         {
00063             alphaContactAngleFvPatchScalarField& acap =
00064                 const_cast<alphaContactAngleFvPatchScalarField&>
00065                 (
00066                     refCast<const alphaContactAngleFvPatchScalarField>
00067                     (
00068                         abf[patchi]
00069                     )
00070                 );
00071 
00072             fvsPatchVectorField& nHatp = nHatb[patchi];
00073             scalarField theta =
00074                 convertToRad*acap.theta(U_.boundaryField()[patchi], nHatp);
00075 
00076             vectorField nf = boundary[patchi].nf();
00077 
00078             // Reset nHatp to correspond to the contact angle
00079 
00080             scalarField a12 = nHatp & nf;
00081 
00082             scalarField b1 = cos(theta);
00083 
00084             scalarField b2(nHatp.size());
00085 
00086             forAll(b2, facei)
00087             {
00088                 b2[facei] = cos(acos(a12[facei]) - theta[facei]);
00089             }
00090 
00091             scalarField det = 1.0 - a12*a12;
00092 
00093             scalarField a = (b1 - a12*b2)/det;
00094             scalarField b = (b2 - a12*b1)/det;
00095 
00096             nHatp = a*nf + b*nHatp;
00097 
00098             nHatp /= (mag(nHatp) + deltaN_.value());
00099 
00100             acap.gradient() = (nf & nHatp)*mag(gradAlphaf[patchi]);
00101             acap.evaluate();
00102         }
00103     }
00104 }
00105 
00106 
00107 void Foam::interfaceProperties::calculateK()
00108 {
00109     const fvMesh& mesh = alpha1_.mesh();
00110     const surfaceVectorField& Sf = mesh.Sf();
00111 
00112     // Cell gradient of alpha
00113     volVectorField gradAlpha = fvc::grad(alpha1_);
00114 
00115     // Interpolated face-gradient of alpha
00116     surfaceVectorField gradAlphaf = fvc::interpolate(gradAlpha);
00117     //gradAlphaf -=
00118     //    (mesh.Sf()/mesh.magSf())
00119     //   *(fvc::snGrad(alpha1_) - (mesh.Sf() & gradAlphaf)/mesh.magSf());
00120 
00121     // Face unit interface normal
00122     surfaceVectorField nHatfv = gradAlphaf/(mag(gradAlphaf) + deltaN_);
00123     correctContactAngle(nHatfv.boundaryField(), gradAlphaf.boundaryField());
00124 
00125     // Face unit interface normal flux
00126     nHatf_ = nHatfv & Sf;
00127 
00128     // Simple expression for curvature
00129     K_ = -fvc::div(nHatf_);
00130 
00131     // Complex expression for curvature.
00132     // Correction is formally zero but numerically non-zero.
00133     /*
00134     volVectorField nHat = gradAlpha/(mag(gradAlpha) + deltaN_);
00135     forAll(nHat.boundaryField(), patchi)
00136     {
00137         nHat.boundaryField()[patchi] = nHatfv.boundaryField()[patchi];
00138     }
00139 
00140     K_ = -fvc::div(nHatf_) + (nHat & fvc::grad(nHatfv) & nHat);
00141     */
00142 }
00143 
00144 
00145 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
00146 
00147 Foam::interfaceProperties::interfaceProperties
00148 (
00149     const volScalarField& alpha1,
00150     const volVectorField& U,
00151     const IOdictionary& dict
00152 )
00153 :
00154     transportPropertiesDict_(dict),
00155     cAlpha_
00156     (
00157         readScalar
00158         (
00159             alpha1.mesh().solutionDict().subDict("PISO").lookup("cAlpha")
00160         )
00161     ),
00162     sigma_(dict.lookup("sigma")),
00163 
00164     deltaN_
00165     (
00166         "deltaN",
00167         1e-8/pow(average(alpha1.mesh().V()), 1.0/3.0)
00168     ),
00169 
00170     alpha1_(alpha1),
00171     U_(U),
00172 
00173     nHatf_
00174     (
00175         IOobject
00176         (
00177             "nHatf",
00178             alpha1_.time().timeName(),
00179             alpha1_.mesh()
00180         ),
00181         alpha1_.mesh(),
00182         dimensionedScalar("nHatf", dimArea, 0.0)
00183     ),
00184 
00185     K_
00186     (
00187         IOobject
00188         (
00189             "Kcond",
00190             alpha1_.time().timeName(),
00191             alpha1_.mesh()
00192         ),
00193         alpha1_.mesh(),
00194         dimensionedScalar("Kcond", dimless/dimLength, 0.0)
00195     )
00196 {
00197     calculateK();
00198 }
00199 
00200 
00201 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines