All Classes Namespaces Files Functions Variables Typedefs Friends Macros Groups
qInterfaceProperties.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-2016 OpenFOAM Foundation
9  Copyright (C) 2019 OpenCFD Ltd.
10  Copyright (C) 2016-2019 ISP RAS (www.ispras.ru) UniCFD Group (www.unicfd.ru)
11 -------------------------------------------------------------------------------
12 License
13  This file is part of OpenFOAM.
14 
15  OpenFOAM is free software: you can redistribute it and/or modify it
16  under the terms of the GNU General Public License as published by
17  the Free Software Foundation, either version 3 of the License, or
18  (at your option) any later version.
19 
20  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
21  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
22  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
23  for more details.
24 
25  You should have received a copy of the GNU General Public License
26  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
27 
28 Group
29  grpQInterfaceProp
30 
31 \*---------------------------------------------------------------------------*/
32 
33 #include "qInterfaceProperties.H"
34 #include "alphaContactAngleTwoPhaseFvPatchScalarField.H"
35 #include "mathematicalConstants.H"
36 #include "surfaceInterpolate.H"
37 #include "fvcDiv.H"
38 #include "fvcGrad.H"
39 #include "fvcSnGrad.H"
40 #include "unitConversion.H"
41 
42 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
43 
44 // Correction for the boundary condition on the unit normal nHat on
45 // walls to produce the correct contact angle.
46 
47 // The dynamic contact angle is calculated from the component of the
48 // velocity on the direction of the interface, parallel to the wall.
49 
50 void Foam::qInterfaceProperties::correctContactAngle
51 (
52  surfaceVectorField::Boundary& nHatb,
53  const surfaceVectorField::Boundary& gradAlphaf
54 ) const
55 {
56  const fvMesh& mesh = alpha1_.mesh();
57  const volScalarField::Boundary& abf = alpha1_.boundaryField();
58 
59  const fvBoundaryMesh& boundary = mesh.boundary();
60 
61  forAll(boundary, patchi)
62  {
63  if (isA<alphaContactAngleTwoPhaseFvPatchScalarField>(abf[patchi]))
64  {
65  alphaContactAngleTwoPhaseFvPatchScalarField& acap =
66  const_cast<alphaContactAngleTwoPhaseFvPatchScalarField&>
67  (
68  refCast<const alphaContactAngleTwoPhaseFvPatchScalarField>
69  (
70  abf[patchi]
71  )
72  );
73 
74  fvsPatchVectorField& nHatp = nHatb[patchi];
75  const scalarField theta
76  (
77  degToRad() * acap.theta(U_.boundaryField()[patchi], nHatp)
78  );
79 
80  const vectorField nf
81  (
82  boundary[patchi].nf()
83  );
84 
85  // Reset nHatp to correspond to the contact angle
86 
87  const scalarField a12(nHatp & nf);
88  const scalarField b1(cos(theta));
89 
90  scalarField b2(nHatp.size());
91  forAll(b2, facei)
92  {
93  b2[facei] = cos(acos(a12[facei]) - theta[facei]);
94  }
95 
96  const scalarField det(1.0 - a12*a12);
97 
98  scalarField a((b1 - a12*b2)/det);
99  scalarField b((b2 - a12*b1)/det);
100 
101  nHatp = a*nf + b*nHatp;
102  nHatp /= (mag(nHatp) + deltaN_.value());
103 
104  acap.gradient() = (nf & nHatp)*mag(gradAlphaf[patchi]);
105  acap.evaluate();
106  }
107  }
108 }
109 
110 
111 void Foam::qInterfaceProperties::calculateK()
112 {
113  const fvMesh& mesh = alpha1_.mesh();
114  const surfaceVectorField& Sf = mesh.Sf();
115  autoPtr<surfaceVectorField> tgradAlphaf;
116  {
117  // Cell gradient of alpha
118  const volVectorField gradAlpha(fvc::grad(alpha1_, "nHat"));
119  tgradAlphaf.set
120  (
121  new surfaceVectorField (linearInterpolate(gradAlpha))
122  );
123  }
124  // Interpolated face-gradient of alpha
125  surfaceVectorField& gradAlphaf = tgradAlphaf();
126 
127  //gradAlphaf -=
128  // (mesh.Sf()/mesh.magSf())
129  // *(fvc::snGrad(alpha1_) - (mesh.Sf() & gradAlphaf)/mesh.magSf());
130 
131  // Face unit interface normal
132  surfaceVectorField nHatfv(gradAlphaf/(mag(gradAlphaf) + deltaN_));
133  // surfaceVectorField nHatfv
134  // (
135  // (gradAlphaf + deltaN_*vector(0, 0, 1)
136  // *sign(gradAlphaf.component(vector::Z)))/(mag(gradAlphaf) + deltaN_)
137  // );
138  correctContactAngle(nHatfv.boundaryFieldRef(), gradAlphaf.boundaryField());
139 
140  // Face unit interface normal flux
141  nHatf_ = nHatfv & Sf;
142 
143  // Simple expression for curvature
144  K_ = -fvc::div(nHatf_);
145 
146  // Complex expression for curvature.
147  // Correction is formally zero but numerically non-zero.
148  /*
149  volVectorField nHat(gradAlpha/(mag(gradAlpha) + deltaN_));
150  forAll(nHat.boundaryField(), patchi)
151  {
152  nHat.boundaryFieldRef()[patchi] = nHatfv.boundaryField()[patchi];
153  }
154 
155  K_ = -fvc::div(nHatf_) + (nHat & fvc::grad(nHatfv) & nHat);
156  */
157 }
158 
159 
160 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
161 
162 Foam::qInterfaceProperties::qInterfaceProperties
163 (
164  const volScalarField& alpha1,
165  const volVectorField& U,
166  const IOdictionary& dict
167 )
168 :
169  transportPropertiesDict_(dict),
170  cAlpha_
171  (
172  alpha1.mesh().solverDict(alpha1.name()).get<scalar>("cAlpha")
173  ),
174 
175  sigmaPtr_(surfaceTensionModel::New(dict, alpha1.mesh())),
176 
177  deltaN_
178  (
179  "deltaN",
180  1e-8/cbrt(average(alpha1.mesh().V()))
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(dimArea, Zero)
196  ),
197 
198  K_
199  (
200  IOobject
201  (
202  "qInterfaceProperties:K",
203  alpha1_.time().timeName(),
204  alpha1_.mesh()
205  ),
206  alpha1_.mesh(),
207  dimensionedScalar(dimless/dimLength, Zero)
208  )
209 {
210  calculateK();
211 }
212 
213 
214 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
215 
216 Foam::tmp<Foam::volScalarField>
218 {
219  return sigmaPtr_->sigma()*K_;
220 }
221 
222 
223 Foam::tmp<Foam::surfaceScalarField>
225 {
226  return fvc::interpolate(sigmaK())*fvc::snGrad(alpha1_);
227 }
228 
229 
230 Foam::tmp<Foam::volScalarField>
232 {
233  return pos0(alpha1_ - 0.01)*pos0(0.99 - alpha1_);
234 }
235 
238 {
239  calculateK();
240 }
241 
242 
244 {
245  alpha1_.mesh().solverDict(alpha1_.name()).readEntry("cAlpha", cAlpha_);
246  sigmaPtr_->readDict(transportPropertiesDict_);
247 
248  return true;
249 }
250 
251 
252 // ************************************************************************* //
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &vf, const surfaceScalarField &dir, const word &reconFieldName=word::null)
Interpolate field vf according to direction dir.
volVectorField U(IOobject("U", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh)
bool read()
Read transportProperties dictionary.
tmp< surfaceScalarField > surfaceTensionForce() const
forAll(Y, i)
Definition: QGDYEqn.H:36
volScalarField & e
Definition: createFields.H:50
tmp< surfaceScalarField > div(const volVectorField &vF)
tmp< volScalarField > sigmaK() const
tmp< surfaceVectorField > grad(const volScalarField &vF)
tmp< volScalarField > nearInterface() const
Indicator of the proximity of the interface.