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 | Copyright (C) 2011-2015 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
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.
15 
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.
20 
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/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "interfaceProperties.H"
27 #include "alphaContactAngleFvPatchScalarField.H"
28 #include "mathematicalConstants.H"
29 #include "surfaceInterpolate.H"
30 #include "fvcDiv.H"
31 #include "fvcGrad.H"
32 #include "fvcSnGrad.H"
33 
34 // * * * * * * * * * * * * * * * Static Member Data * * * * * * * * * * * * //
35 
36 // Symbol to force loading at runtime
37 extern "C"
39 {}
40 
41 
42 const Foam::scalar Foam::interfaceProperties::convertToRad =
44 
45 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
46 
47 // Correction for the boundary condition on the unit normal nHat on
48 // walls to produce the correct contact angle.
49 
50 // The dynamic contact angle is calculated from the component of the
51 // velocity on the direction of the interface, parallel to the wall.
52 
54 (
57 ) const
58 {
59  const fvMesh& mesh = alpha1_.mesh();
60  const volScalarField::GeometricBoundaryField& abf = alpha1_.boundaryField();
61 
62  const fvBoundaryMesh& boundary = mesh.boundary();
63 
65  {
66  if (isA<alphaContactAngleFvPatchScalarField>(abf[patchi]))
67  {
70  (
71  refCast<const alphaContactAngleFvPatchScalarField>
72  (
73  abf[patchi]
74  )
75  );
76 
77  fvsPatchVectorField& nHatp = nHatb[patchi];
78  const scalarField theta
79  (
80  convertToRad*acap.theta(U_.boundaryField()[patchi], nHatp)
81  );
82 
83  const vectorField nf
84  (
85  boundary[patchi].nf()
86  );
87 
88  // Reset nHatp to correspond to the contact angle
89 
90  const scalarField a12(nHatp & nf);
91  const scalarField b1(cos(theta));
92 
93  scalarField b2(nHatp.size());
94  forAll(b2, facei)
95  {
96  b2[facei] = cos(acos(a12[facei]) - theta[facei]);
97  }
98 
99  const scalarField det(1.0 - a12*a12);
100 
101  scalarField a((b1 - a12*b2)/det);
102  scalarField b((b2 - a12*b1)/det);
103 
104  nHatp = a*nf + b*nHatp;
105  nHatp /= (mag(nHatp) + deltaN_.value());
106 
107  acap.gradient() = (nf & nHatp)*mag(gradAlphaf[patchi]);
108  acap.evaluate();
109  }
110  }
111 }
112 
113 
115 {
116  const fvMesh& mesh = alpha1_.mesh();
117  const surfaceVectorField& Sf = mesh.Sf();
118 
119  // Cell gradient of alpha
120  const volVectorField gradAlpha(fvc::grad(alpha1_, "nHat"));
121 
122  // Interpolated face-gradient of alpha
123  surfaceVectorField gradAlphaf(fvc::interpolate(gradAlpha));
124 
125  //gradAlphaf -=
126  // (mesh.Sf()/mesh.magSf())
127  // *(fvc::snGrad(alpha1_) - (mesh.Sf() & gradAlphaf)/mesh.magSf());
128 
129  // Face unit interface normal
130  surfaceVectorField nHatfv(gradAlphaf/(mag(gradAlphaf) + deltaN_));
131  // surfaceVectorField nHatfv
132  // (
133  // (gradAlphaf + deltaN_*vector(0, 0, 1)
134  // *sign(gradAlphaf.component(vector::Z)))/(mag(gradAlphaf) + deltaN_)
135  // );
136  correctContactAngle(nHatfv.boundaryField(), gradAlphaf.boundaryField());
137 
138  // Face unit interface normal flux
139  nHatf_ = nHatfv & Sf;
140 
141  // Simple expression for curvature
142  K_ = -fvc::div(nHatf_);
143 
144  // Complex expression for curvature.
145  // Correction is formally zero but numerically non-zero.
146  /*
147  volVectorField nHat(gradAlpha/(mag(gradAlpha) + deltaN_));
148  forAll(nHat.boundaryField(), patchi)
149  {
150  nHat.boundaryField()[patchi] = nHatfv.boundaryField()[patchi];
151  }
152 
153  K_ = -fvc::div(nHatf_) + (nHat & fvc::grad(nHatfv) & nHat);
154  */
155 }
156 
157 
158 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
159 
161 (
162  const volScalarField& alpha1,
163  const volVectorField& U,
164  const IOdictionary& dict
165 )
166 :
167  transportPropertiesDict_(dict),
168  cAlpha_
169  (
170  readScalar
171  (
172  alpha1.mesh().solverDict(alpha1.name()).lookup("cAlpha")
173  )
174  ),
175  sigma_("sigma", dimensionSet(1, 0, -2, 0, 0), dict),
176 
177  deltaN_
178  (
179  "deltaN",
180  1e-8/pow(average(alpha1.mesh().V()), 1.0/3.0)
181  ),
182 
183  alpha1_(alpha1),
184  U_(U),
185 
186  nHatf_
187  (
188  IOobject
189  (
190  "nHatf",
191  alpha1_.time().timeName(),
192  alpha1_.mesh()
193  ),
194  alpha1_.mesh(),
195  dimensionedScalar("nHatf", dimArea, 0.0)
196  ),
197 
198  K_
199  (
200  IOobject
201  (
202  "interfaceProperties:K",
203  alpha1_.time().timeName(),
204  alpha1_.mesh()
205  ),
206  alpha1_.mesh(),
208  )
209 {
210  calculateK();
211 }
212 
213 
214 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
215 
218 {
219  return fvc::interpolate(sigmaK())*fvc::snGrad(alpha1_);
220 }
221 
222 
225 {
226  return pos(alpha1_ - 0.01)*pos(0.99 - alpha1_);
227 }
228 
229 
231 {
232  alpha1_.mesh().solverDict(alpha1_.name()).lookup("cAlpha") >> cAlpha_;
233  transportPropertiesDict_.lookup("sigma") >> sigma_;
234 
235  return true;
236 }
237 
238 
239 // ************************************************************************* //
Foam::IOdictionary
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:53
Foam::fvc::snGrad
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:45
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
mathematicalConstants.H
Foam::interfaceProperties::nHatf_
surfaceScalarField nHatf_
Definition: interfaceProperties.H:72
Foam::dimless
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
Foam::fvc::grad
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
Foam::interfaceProperties::alpha1_
const volScalarField & alpha1_
Definition: interfaceProperties.H:70
Foam::interfaceProperties::correctContactAngle
void correctContactAngle(surfaceVectorField::GeometricBoundaryField &nHat, surfaceVectorField::GeometricBoundaryField &gradAlphaf) const
Correction for the boundary condition on the unit normal nHat on.
Definition: interfaceProperties.C:54
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::interfaceProperties::interfaceProperties
interfaceProperties(const interfaceProperties &)
Disallow default bitwise copy construct and assignment.
interfacePropertiesLoad
void interfacePropertiesLoad()
Definition: interfaceProperties.C:38
fvcSnGrad.H
Calculate the snGrad of the given volField.
boundary
faceListList boundary(nPatches)
Foam::fvc::interpolate
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &vf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
Definition: surfaceInterpolate.C:110
Foam::fvsPatchField
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:65
fvcDiv.H
Calculate the divergence of the given field.
Foam::alphaContactAngleFvPatchScalarField::theta
virtual tmp< scalarField > theta(const fvPatchVectorField &Up, const fvsPatchVectorField &nHat) const =0
Return the contact angle.
Foam::fvc::div
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
Foam::GeometricField::boundaryField
GeometricBoundaryField & boundaryField()
Return reference to GeometricBoundaryField.
Definition: GeometricField.C:735
Foam::interfaceProperties::calculateK
void calculateK()
Re-calculate the interface curvature.
Definition: interfaceProperties.C:114
Foam::dimensioned::value
const Type & value() const
Return const reference to value.
Definition: dimensionedType.C:261
Foam::dimensionSet
Dimension set for the base types.
Definition: dimensionSet.H:116
Foam::interfaceProperties::K_
volScalarField K_
Definition: interfaceProperties.H:73
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
U
U
Definition: pEqn.H:46
Foam::fvBoundaryMesh
Foam::fvBoundaryMesh.
Definition: fvBoundaryMesh.H:52
Foam::interfaceProperties::deltaN_
const dimensionedScalar deltaN_
Stabilisation for normalisation of the interface normal.
Definition: interfaceProperties.H:68
Foam::interfaceProperties::nearInterface
tmp< volScalarField > nearInterface() const
Indicator of the proximity of the interface.
Definition: interfaceProperties.C:224
Foam::dimArea
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:57
Foam::constant::physicoChemical::b
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:28
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
interfaceProperties.H
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:73
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::interfaceProperties::surfaceTensionForce
tmp< surfaceScalarField > surfaceTensionForce() const
Definition: interfaceProperties.C:217
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:41
surfaceInterpolate.H
Surface Interpolation.
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Foam::e
const double e
Elementary charge.
Definition: doubleFloat.H:94
Foam::alphaContactAngleFvPatchScalarField
Abstract base class for alphaContactAngle boundary conditions.
Definition: alphaContactAngleFvPatchScalarField.H:81
alpha1
volScalarField & alpha1
Definition: createFields.H:15
Foam::readScalar
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if succesful.
Definition: doubleScalar.H:63
Foam::acos
dimensionedScalar acos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:259
fvcGrad.H
Calculate the gradient of the given field.
patchi
label patchi
Definition: getPatchFieldScalar.H:1
Foam::constant::mathematical::pi
const scalar pi(M_PI)
Foam::det
dimensionedScalar det(const dimensionedSphericalTensor &dt)
Definition: dimensionedSphericalTensor.C:60
Foam::alphaContactAngleFvPatchScalarField::evaluate
virtual void evaluate(const Pstream::commsTypes commsType=Pstream::blocking)
Evaluate the patch field.
Definition: alphaContactAngleFvPatchScalarField.C:133
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::interfaceProperties::convertToRad
static const scalar convertToRad
Conversion factor for degrees into radians.
Definition: interfaceProperties.H:98
Foam::GeometricField::GeometricBoundaryField
Definition: GeometricField.H:105
Foam::average
dimensioned< Type > average(const DimensionedField< Type, GeoMesh > &df)
Definition: DimensionedFieldFunctions.C:335
Foam::interfaceProperties::read
bool read()
Read transportProperties dictionary.
Definition: interfaceProperties.C:230
Foam::cos
dimensionedScalar cos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:256
Foam::pos
dimensionedScalar pos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:190