interfaceProperties.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2017 OpenFOAM Foundation
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
28#include "interfaceProperties.H"
31#include "surfaceInterpolate.H"
32#include "fvcDiv.H"
33#include "fvcGrad.H"
34#include "fvcSnGrad.H"
35#include "unitConversion.H"
36
37// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
38
39// Correction for the boundary condition on the unit normal nHat on
40// walls to produce the correct contact angle.
41
42// The dynamic contact angle is calculated from the component of the
43// velocity on the direction of the interface, parallel to the wall.
44
45void Foam::interfaceProperties::correctContactAngle
46(
47 surfaceVectorField::Boundary& nHatb,
48 const surfaceVectorField::Boundary& gradAlphaf
49) const
50{
51 const fvMesh& mesh = alpha1_.mesh();
52 const volScalarField::Boundary& abf = alpha1_.boundaryField();
53
54 const fvBoundaryMesh& boundary = mesh.boundary();
55
56 forAll(boundary, patchi)
57 {
58 if (isA<alphaContactAngleTwoPhaseFvPatchScalarField>(abf[patchi]))
59 {
60 alphaContactAngleTwoPhaseFvPatchScalarField& acap =
61 const_cast<alphaContactAngleTwoPhaseFvPatchScalarField&>
62 (
63 refCast<const alphaContactAngleTwoPhaseFvPatchScalarField>
64 (
65 abf[patchi]
66 )
67 );
68
69 fvsPatchVectorField& nHatp = nHatb[patchi];
70 const scalarField theta
71 (
72 degToRad() * acap.theta(U_.boundaryField()[patchi], nHatp)
73 );
74
75 const vectorField nf
76 (
77 boundary[patchi].nf()
78 );
79
80 // Reset nHatp to correspond to the contact angle
81
82 const scalarField a12(nHatp & nf);
83 const scalarField b1(cos(theta));
84
85 scalarField b2(nHatp.size());
86 forAll(b2, facei)
87 {
88 b2[facei] = cos(acos(a12[facei]) - theta[facei]);
89 }
90
91 const scalarField det(1.0 - a12*a12);
92
93 scalarField a((b1 - a12*b2)/det);
94 scalarField b((b2 - a12*b1)/det);
95
96 nHatp = a*nf + b*nHatp;
97 nHatp /= (mag(nHatp) + deltaN_.value());
98
99 acap.gradient() = (nf & nHatp)*mag(gradAlphaf[patchi]);
100 acap.evaluate();
101 }
102 }
103}
104
105
106void Foam::interfaceProperties::calculateK()
107{
108 const fvMesh& mesh = alpha1_.mesh();
109 const surfaceVectorField& Sf = mesh.Sf();
110
111 // Cell gradient of alpha
112 const volVectorField gradAlpha(fvc::grad(alpha1_, "nHat"));
113
114 // Interpolated face-gradient of alpha
115 surfaceVectorField gradAlphaf(fvc::interpolate(gradAlpha));
116
117 //gradAlphaf -=
118 // (mesh.Sf()/mesh.magSf())
119 // *(fvc::snGrad(alpha1_) - (mesh.Sf() & gradAlphaf)/mesh.magSf());
120
121 // Face unit interface normal
122 surfaceVectorField nHatfv(gradAlphaf/(mag(gradAlphaf) + deltaN_));
123 // surfaceVectorField nHatfv
124 // (
125 // (gradAlphaf + deltaN_*vector(0, 0, 1)
126 // *sign(gradAlphaf.component(vector::Z)))/(mag(gradAlphaf) + deltaN_)
127 // );
128 correctContactAngle(nHatfv.boundaryFieldRef(), gradAlphaf.boundaryField());
129
130 // Face unit interface normal flux
131 nHatf_ = nHatfv & Sf;
132
133 // Simple expression for curvature
134 K_ = -fvc::div(nHatf_);
135
136 // Complex expression for curvature.
137 // Correction is formally zero but numerically non-zero.
138 /*
139 volVectorField nHat(gradAlpha/(mag(gradAlpha) + deltaN_));
140 forAll(nHat.boundaryField(), patchi)
141 {
142 nHat.boundaryFieldRef()[patchi] = nHatfv.boundaryField()[patchi];
143 }
144
145 K_ = -fvc::div(nHatf_) + (nHat & fvc::grad(nHatfv) & nHat);
146 */
147}
148
149
150// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
151
153(
154 const volScalarField& alpha1,
155 const volVectorField& U,
156 const IOdictionary& dict
157)
158:
159 transportPropertiesDict_(dict),
160 cAlpha_
161 (
162 alpha1.mesh().solverDict(alpha1.name()).get<scalar>("cAlpha")
163 ),
164
165 sigmaPtr_(surfaceTensionModel::New(dict, alpha1.mesh())),
166
167 deltaN_
168 (
169 "deltaN",
170 1e-8/cbrt(average(alpha1.mesh().V()))
171 ),
172
173 alpha1_(alpha1),
174 U_(U),
175
176 nHatf_
177 (
179 (
180 "nHatf",
181 alpha1_.time().timeName(),
182 alpha1_.mesh()
183 ),
184 alpha1_.mesh(),
186 ),
187
188 K_
189 (
191 (
192 "interfaceProperties:K",
193 alpha1_.time().timeName(),
194 alpha1_.mesh()
195 ),
196 alpha1_.mesh(),
198 )
199{
200 calculateK();
201}
202
203
204// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
205
208{
209 return sigmaPtr_->sigma()*K_;
210}
211
212
215{
216 return fvc::interpolate(sigmaK())*fvc::snGrad(alpha1_);
217}
218
219
222{
223 return pos0(alpha1_ - 0.01)*pos0(0.99 - alpha1_);
224}
225
226
228{
229 calculateK();
230}
231
232
234{
235 alpha1_.mesh().solverDict(alpha1_.name()).readEntry("cAlpha", cAlpha_);
236 sigmaPtr_->readDict(transportPropertiesDict_);
237
238 return true;
239}
240
241
242// ************************************************************************* //
const volScalarField & alpha1
const Mesh & mesh() const
Return mesh.
GeometricBoundaryField< scalar, fvPatchField, volMesh > Boundary
Type of boundary fields.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:57
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
const Type & value() const
Return const reference to value.
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:712
const surfaceVectorField & Sf() const
Return cell face area vectors.
Contains the interface properties.
tmp< surfaceScalarField > surfaceTensionForce() const
tmp< volScalarField > nearInterface() const
Indicator of the proximity of the interface.
tmp< volScalarField > sigmaK() const
bool read()
Read transportProperties dictionary.
A class for managing temporary objects.
Definition: tmp.H:65
U
Definition: pEqn.H:72
faceListList boundary
dynamicFvMesh & mesh
Calculate the divergence of the given field.
Calculate the gradient of the given field.
Calculate the snGrad of the given volField.
word timeName
Definition: getTimeIndex.H:3
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:54
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:47
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:49
GeometricField< vector, fvsPatchField, surfaceMesh > surfaceVectorField
dimensionedScalar det(const dimensionedSphericalTensor &dt)
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:83
const dimensionSet dimless
Dimensionless.
dimensionedScalar pos0(const dimensionedScalar &ds)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:59
dimensioned< Type > average(const DimensionedField< Type, GeoMesh > &df)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
constexpr scalar degToRad() noexcept
Multiplication factor for degrees to radians conversion.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Field< vector > vectorField
Specialisation of Field<T> for vector.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
dimensionedScalar cbrt(const dimensionedScalar &ds)
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
fvsPatchField< vector > fvsPatchVectorField
dimensionedScalar cos(const dimensionedScalar &ds)
dimensionedScalar acos(const dimensionedScalar &ds)
dictionary dict
volScalarField & b
Definition: createFields.H:27
volScalarField & e
Definition: createFields.H:11
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
Unit conversion functions.