interfaceProperties.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2020 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 "unitConversion.H"
29 #include "surfaceInterpolate.H"
30 #include "fvcDiv.H"
31 #include "fvcGrad.H"
32 #include "fvcSnGrad.H"
33 
34 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
35 
36 // Correction for the boundary condition on the unit normal nHat on
37 // walls to produce the correct contact angle.
38 
39 // The dynamic contact angle is calculated from the component of the
40 // velocity on the direction of the interface, parallel to the wall.
41 
42 void Foam::interfaceProperties::correctContactAngle
43 (
44  surfaceVectorField::Boundary& nHatb,
45  const surfaceVectorField::Boundary& gradAlphaf
46 )
47 {
48  const fvMesh& mesh = alpha1_.mesh();
49  volScalarField::Boundary& a1bf = alpha1_.boundaryFieldRef();
50  volScalarField::Boundary& a2bf = alpha2_.boundaryFieldRef();
51 
52  const fvBoundaryMesh& boundary = mesh.boundary();
53 
54  forAll(boundary, patchi)
55  {
56  if (isA<alphaContactAngleFvPatchScalarField>(a1bf[patchi]))
57  {
58  alphaContactAngleFvPatchScalarField& a1cap =
59  refCast<alphaContactAngleFvPatchScalarField>
60  (
61  a1bf[patchi]
62  );
63 
64  fvsPatchVectorField& nHatp = nHatb[patchi];
65  const scalarField theta
66  (
67  degToRad(a1cap.theta(U_.boundaryField()[patchi], nHatp))
68  );
69 
70  const vectorField nf
71  (
72  boundary[patchi].nf()
73  );
74 
75  // Reset nHatp to correspond to the contact angle
76 
77  const scalarField a12(nHatp & nf);
78  const scalarField b1(cos(theta));
79 
80  scalarField b2(nHatp.size());
81  forAll(b2, facei)
82  {
83  b2[facei] = cos(acos(a12[facei]) - theta[facei]);
84  }
85 
86  const scalarField det(1.0 - a12*a12);
87 
88  scalarField a((b1 - a12*b2)/det);
89  scalarField b((b2 - a12*b1)/det);
90 
91  nHatp = a*nf + b*nHatp;
92  nHatp /= (mag(nHatp) + deltaN_.value());
93 
94  a1cap.gradient() = (nf & nHatp)*mag(gradAlphaf[patchi]);
95  a1cap.evaluate();
96  a2bf[patchi] = 1 - a1cap;
97  }
98  }
99 }
100 
101 
102 void Foam::interfaceProperties::calculateK()
103 {
104  const fvMesh& mesh = alpha1_.mesh();
105  const surfaceVectorField& Sf = mesh.Sf();
106 
107  // Cell gradient of alpha
108  const volVectorField gradAlpha(fvc::grad(alpha1_, "nHat"));
109 
110  // Interpolated face-gradient of alpha
111  surfaceVectorField gradAlphaf(fvc::interpolate(gradAlpha));
112 
113  // gradAlphaf -=
114  // (mesh.Sf()/mesh.magSf())
115  // *(fvc::snGrad(alpha1_) - (mesh.Sf() & gradAlphaf)/mesh.magSf());
116 
117  // Face unit interface normal
118  surfaceVectorField nHatfv(gradAlphaf/(mag(gradAlphaf) + deltaN_));
119  // surfaceVectorField nHatfv
120  // (
121  // (gradAlphaf + deltaN_*vector(0, 0, 1)
122  // *sign(gradAlphaf.component(vector::Z)))/(mag(gradAlphaf) + deltaN_)
123  // );
124  correctContactAngle(nHatfv.boundaryFieldRef(), gradAlphaf.boundaryField());
125 
126  // Face unit interface normal flux
127  nHatf_ = nHatfv & Sf;
128 
129  // Simple expression for curvature
130  K_ = -fvc::div(nHatf_);
131 
132  // Complex expression for curvature.
133  // Correction is formally zero but numerically non-zero.
134  /*
135  volVectorField nHat(gradAlpha/(mag(gradAlpha) + deltaN_));
136  forAll(nHat.boundaryField(), patchi)
137  {
138  nHat.boundaryField()[patchi] = nHatfv.boundaryField()[patchi];
139  }
140 
141  K_ = -fvc::div(nHatf_) + (nHat & fvc::grad(nHatfv) & nHat);
142  */
143 }
144 
145 
146 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
147 
149 (
150  volScalarField& alpha1,
151  volScalarField& alpha2,
152  const volVectorField& U,
153  const IOdictionary& dict
154 )
155 :
156  transportPropertiesDict_(dict),
157 
158  sigmaPtr_(surfaceTensionModel::New(dict, alpha1.mesh())),
159 
160  deltaN_
161  (
162  "deltaN",
163  1e-8/pow(average(alpha1.mesh().V()), 1.0/3.0)
164  ),
165 
166  alpha1_(alpha1),
167  alpha2_(alpha2),
168  U_(U),
169 
170  nHatf_
171  (
172  IOobject
173  (
174  "nHatf",
175  alpha1_.time().timeName(),
176  alpha1_.mesh()
177  ),
178  alpha1_.mesh(),
180  ),
181 
182  K_
183  (
184  IOobject
185  (
186  "interfaceProperties:K",
187  alpha1_.time().timeName(),
188  alpha1_.mesh()
189  ),
190  alpha1_.mesh(),
192  )
193 {
194  calculateK();
195 }
196 
197 
198 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
199 
202 {
203  return sigmaPtr_->sigma()*K_;
204 }
205 
206 
209 {
210  return fvc::interpolate(sigmaK())*fvc::snGrad(alpha1_);
211 }
212 
213 
216 {
217  return pos0(alpha1_ - 0.01)*pos0(0.99 - alpha1_);
218 }
219 
220 
222 {
223  calculateK();
224 }
225 
226 
228 {
229  sigmaPtr_->readDict(transportPropertiesDict_);
230 
231  return true;
232 }
233 
234 
235 // ************************************************************************* //
fvsPatchField< vector > fvsPatchVectorField
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
dimensionedScalar acos(const dimensionedScalar &ds)
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
const dimensionSet dimArea
const Boundary & boundaryField() const
Return const-reference to the boundary field.
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
Unit conversion functions.
GeometricField< vector, fvsPatchField, surfaceMesh > surfaceVectorField
const dimensionedScalar b
Wien displacement law constant: default SI units: [m K].
Definition: createFields.H:27
static autoPtr< surfaceTensionModel > New(const dictionary &dict, const fvMesh &mesh)
bool read()
Read transportProperties dictionary.
Calculate the snGrad of the given volField.
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:636
const dimensionSet dimless
dimensionedScalar det(const dimensionedSphericalTensor &dt)
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:58
tmp< volScalarField > nearInterface() const
Indicator of the proximity of the interface.
tmp< surfaceScalarField > surfaceTensionForce() const
scalar degToRad(const scalar deg)
Conversion from degrees to radians.
const dimensionSet dimLength
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:53
dimensionedScalar cos(const dimensionedScalar &ds)
Calculate the gradient of the given field.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const Type & value() const
Return const reference to value.
dimensioned< Type > average(const DimensionedField< Type, GeoMesh > &df)
Calculate the divergence of the given field.
dimensionedScalar pos0(const dimensionedScalar &ds)
const Mesh & mesh() const
Return mesh.
tmp< volScalarField > sigmaK() const
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.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
label patchi
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const Time & time() const
Return time.
Definition: IOobject.C:318
interfaceProperties(volScalarField &alpha1, volScalarField &alpha2, const volVectorField &U, const IOdictionary &)
Construct from volume fraction field gamma and IOdictionary.
dimensioned< scalar > mag(const dimensioned< Type > &)
Field< vector > vectorField
Specialisation of Field<T> for vector.
const doubleScalar e
Elementary charge.
Definition: doubleScalar.H:105
A class for managing temporary objects.
Definition: PtrList.H:53
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:45