Initial commit for version 2.0.x patch release
[OpenFOAM-2.0.x.git] / src / transportModels / twoPhaseInterfaceProperties / alphaContactAngle / dynamicAlphaContactAngle / dynamicAlphaContactAngleFvPatchScalarField.C
blobcd21ec9a5211844abbcd823024a3d59fe493b003
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2004-2010 OpenCFD Ltd.
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
9     This file is part of OpenFOAM.
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.
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.
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/>.
24 \*---------------------------------------------------------------------------*/
26 #include "dynamicAlphaContactAngleFvPatchScalarField.H"
27 #include "addToRunTimeSelectionTable.H"
28 #include "fvPatchFieldMapper.H"
29 #include "volMesh.H"
31 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
33 Foam::dynamicAlphaContactAngleFvPatchScalarField::
34 dynamicAlphaContactAngleFvPatchScalarField
36     const fvPatch& p,
37     const DimensionedField<scalar, volMesh>& iF
40     alphaContactAngleFvPatchScalarField(p, iF),
41     theta0_(0.0),
42     uTheta_(0.0),
43     thetaA_(0.0),
44     thetaR_(0.0)
48 Foam::dynamicAlphaContactAngleFvPatchScalarField::
49 dynamicAlphaContactAngleFvPatchScalarField
51     const dynamicAlphaContactAngleFvPatchScalarField& gcpsf,
52     const fvPatch& p,
53     const DimensionedField<scalar, volMesh>& iF,
54     const fvPatchFieldMapper& mapper
57     alphaContactAngleFvPatchScalarField(gcpsf, p, iF, mapper),
58     theta0_(gcpsf.theta0_),
59     uTheta_(gcpsf.uTheta_),
60     thetaA_(gcpsf.thetaA_),
61     thetaR_(gcpsf.thetaR_)
65 Foam::dynamicAlphaContactAngleFvPatchScalarField::
66 dynamicAlphaContactAngleFvPatchScalarField
68     const fvPatch& p,
69     const DimensionedField<scalar, volMesh>& iF,
70     const dictionary& dict
73     alphaContactAngleFvPatchScalarField(p, iF, dict),
74     theta0_(readScalar(dict.lookup("theta0"))),
75     uTheta_(readScalar(dict.lookup("uTheta"))),
76     thetaA_(readScalar(dict.lookup("thetaA"))),
77     thetaR_(readScalar(dict.lookup("thetaR")))
79     evaluate();
83 Foam::dynamicAlphaContactAngleFvPatchScalarField::
84 dynamicAlphaContactAngleFvPatchScalarField
86     const dynamicAlphaContactAngleFvPatchScalarField& gcpsf
89     alphaContactAngleFvPatchScalarField(gcpsf),
90     theta0_(gcpsf.theta0_),
91     uTheta_(gcpsf.uTheta_),
92     thetaA_(gcpsf.thetaA_),
93     thetaR_(gcpsf.thetaR_)
97 Foam::dynamicAlphaContactAngleFvPatchScalarField::
98 dynamicAlphaContactAngleFvPatchScalarField
100     const dynamicAlphaContactAngleFvPatchScalarField& gcpsf,
101     const DimensionedField<scalar, volMesh>& iF
104     alphaContactAngleFvPatchScalarField(gcpsf, iF),
105     theta0_(gcpsf.theta0_),
106     uTheta_(gcpsf.uTheta_),
107     thetaA_(gcpsf.thetaA_),
108     thetaR_(gcpsf.thetaR_)
112 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
114 Foam::tmp<Foam::scalarField>
115 Foam::dynamicAlphaContactAngleFvPatchScalarField::theta
117     const fvPatchVectorField& Up,
118     const fvsPatchVectorField& nHat
119 ) const
121     if (uTheta_ < SMALL)
122     {
123         return tmp<scalarField>(new scalarField(size(), theta0_));
124     }
126     const vectorField nf(patch().nf());
128     // Calculated the component of the velocity parallel to the wall
129     vectorField Uwall(Up.patchInternalField() - Up);
130     Uwall -= (nf & Uwall)*nf;
132     // Find the direction of the interface parallel to the wall
133     vectorField nWall(nHat - (nf & nHat)*nf);
135     // Normalise nWall
136     nWall /= (mag(nWall) + SMALL);
138     // Calculate Uwall resolved normal to the interface parallel to
139     // the interface
140     scalarField uwall(nWall & Uwall);
142     return theta0_ + (thetaA_ - thetaR_)*tanh(uwall/uTheta_);
146 void Foam::dynamicAlphaContactAngleFvPatchScalarField::write(Ostream& os) const
148     alphaContactAngleFvPatchScalarField::write(os);
149     os.writeKeyword("theta0") << theta0_ << token::END_STATEMENT << nl;
150     os.writeKeyword("uTheta") << uTheta_ << token::END_STATEMENT << nl;
151     os.writeKeyword("thetaA") << thetaA_ << token::END_STATEMENT << nl;
152     os.writeKeyword("thetaR") << thetaR_ << token::END_STATEMENT << nl;
153     writeEntry("value", os);
157 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
159 namespace Foam
161     makePatchTypeField
162     (
163         fvPatchScalarField,
164         dynamicAlphaContactAngleFvPatchScalarField
165     );
169 // ************************************************************************* //