initial commit for version 1.5.x patch release
[OpenFOAM-1.5.x.git] / src / transportModels / interfaceProperties / gammaContactAngle / dynamicGammaContactAngle / dynamicGammaContactAngleFvPatchScalarField.C
blobe96567115077f440755145937710d29d16290ad7
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2008 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 the
13     Free Software Foundation; either version 2 of the License, or (at your
14     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, write to the Free Software Foundation,
23     Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 \*---------------------------------------------------------------------------*/
27 #include "dynamicGammaContactAngleFvPatchScalarField.H"
28 #include "addToRunTimeSelectionTable.H"
29 #include "fvPatchFieldMapper.H"
30 #include "volMesh.H"
32 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
34 Foam::dynamicGammaContactAngleFvPatchScalarField::
35 dynamicGammaContactAngleFvPatchScalarField
37     const fvPatch& p,
38     const DimensionedField<scalar, volMesh>& iF
41     gammaContactAngleFvPatchScalarField(p, iF),
42     theta0_(0.0),
43     uTheta_(0.0),
44     thetaA_(0.0),
45     thetaR_(0.0)
49 Foam::dynamicGammaContactAngleFvPatchScalarField::
50 dynamicGammaContactAngleFvPatchScalarField
52     const dynamicGammaContactAngleFvPatchScalarField& gcpsf,
53     const fvPatch& p,
54     const DimensionedField<scalar, volMesh>& iF,
55     const fvPatchFieldMapper& mapper
58     gammaContactAngleFvPatchScalarField(gcpsf, p, iF, mapper),
59     theta0_(gcpsf.theta0_),
60     uTheta_(gcpsf.uTheta_),
61     thetaA_(gcpsf.thetaA_),
62     thetaR_(gcpsf.thetaR_)
66 Foam::dynamicGammaContactAngleFvPatchScalarField::
67 dynamicGammaContactAngleFvPatchScalarField
69     const fvPatch& p,
70     const DimensionedField<scalar, volMesh>& iF,
71     const dictionary& dict
74     gammaContactAngleFvPatchScalarField(p, iF),
75     theta0_(readScalar(dict.lookup("theta0"))),
76     uTheta_(readScalar(dict.lookup("uTheta"))),
77     thetaA_(readScalar(dict.lookup("thetaA"))),
78     thetaR_(readScalar(dict.lookup("thetaR")))
80     evaluate();
84 Foam::dynamicGammaContactAngleFvPatchScalarField::
85 dynamicGammaContactAngleFvPatchScalarField
87     const dynamicGammaContactAngleFvPatchScalarField& gcpsf
90     gammaContactAngleFvPatchScalarField(gcpsf),
91     theta0_(gcpsf.theta0_),
92     uTheta_(gcpsf.uTheta_),
93     thetaA_(gcpsf.thetaA_),
94     thetaR_(gcpsf.thetaR_)
98 Foam::dynamicGammaContactAngleFvPatchScalarField::
99 dynamicGammaContactAngleFvPatchScalarField
101     const dynamicGammaContactAngleFvPatchScalarField& gcpsf,
102     const DimensionedField<scalar, volMesh>& iF
105     gammaContactAngleFvPatchScalarField(gcpsf, iF),
106     theta0_(gcpsf.theta0_),
107     uTheta_(gcpsf.uTheta_),
108     thetaA_(gcpsf.thetaA_),
109     thetaR_(gcpsf.thetaR_)
113 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
115 Foam::tmp<Foam::scalarField>
116 Foam::dynamicGammaContactAngleFvPatchScalarField::theta
118     const fvPatchVectorField& Up,
119     const fvsPatchVectorField& nHat
120 ) const
122     if (uTheta_ < SMALL)
123     {
124         return tmp<scalarField>(new scalarField(size(), theta0_));
125     }
127     vectorField nf = patch().nf();
129     // Calculated the component of the velocity parallel to the wall
130     vectorField Uwall = Up.patchInternalField() - Up;
131     Uwall -= (nf & Uwall)*nf;
133     // Find the direction of the interface parallel to the wall
134     vectorField nWall = nHat - (nf & nHat)*nf;
136     // Normalise nWall
137     nWall /= (mag(nWall) + SMALL);
139     // Calculate Uwall resolved normal to the interface parallel to
140     // the interface
141     scalarField uwall = nWall & Uwall;
143     return theta0_ + (thetaA_ - thetaR_)*tanh(uwall/uTheta_);
147 void Foam::dynamicGammaContactAngleFvPatchScalarField::write(Ostream& os) const
149     fvPatchScalarField::write(os);
150     os.writeKeyword("theta0") << theta0_ << token::END_STATEMENT << nl;
151     os.writeKeyword("uTheta") << uTheta_ << token::END_STATEMENT << nl;
152     os.writeKeyword("thetaA") << thetaA_ << token::END_STATEMENT << nl;
153     os.writeKeyword("thetaR") << thetaR_ << token::END_STATEMENT << nl;
154     writeEntry("value", os);
158 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
160 namespace Foam
162     makePatchTypeField
163     (
164         fvPatchScalarField,
165         dynamicGammaContactAngleFvPatchScalarField
166     );
170 // ************************************************************************* //