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-2021 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 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38  defineTypeNameAndDebug(interfaceProperties, 0);
39 }
40 
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::interfaceProperties::correctContactAngle
51 (
53  const surfaceVectorField::Boundary& gradAlphaf
54 )
55 {
56  const fvMesh& mesh = alpha1_.mesh();
57  volScalarField::Boundary& a1bf = alpha1_.boundaryFieldRef();
58  volScalarField::Boundary& a2bf = alpha2_.boundaryFieldRef();
59 
60  const fvBoundaryMesh& boundary = mesh.boundary();
61 
62  forAll(boundary, patchi)
63  {
64  if (isA<alphaContactAngleFvPatchScalarField>(a1bf[patchi]))
65  {
66  alphaContactAngleFvPatchScalarField& a1cap =
67  refCast<alphaContactAngleFvPatchScalarField>
68  (
69  a1bf[patchi]
70  );
71 
72  fvsPatchVectorField& nHatp = nHatb[patchi];
73  const scalarField theta
74  (
75  degToRad(a1cap.theta(U_.boundaryField()[patchi], nHatp))
76  );
77 
78  const vectorField nf
79  (
80  boundary[patchi].nf()
81  );
82 
83  // Reset nHatp to correspond to the contact angle
84 
85  const scalarField a12(nHatp & nf);
86  const scalarField b1(cos(theta));
87 
88  scalarField b2(nHatp.size());
89  forAll(b2, facei)
90  {
91  b2[facei] = cos(acos(a12[facei]) - theta[facei]);
92  }
93 
94  const scalarField det(1.0 - a12*a12);
95 
96  scalarField a((b1 - a12*b2)/det);
97  scalarField b((b2 - a12*b1)/det);
98 
99  nHatp = a*nf + b*nHatp;
100  nHatp /= (mag(nHatp) + deltaN_.value());
101 
102  a1cap.gradient() = (nf & nHatp)*mag(gradAlphaf[patchi]);
103  a1cap.evaluate();
104  a2bf[patchi] = 1 - a1cap;
105  }
106  }
107 }
108 
109 
110 void Foam::interfaceProperties::calculateK()
111 {
112  const fvMesh& mesh = alpha1_.mesh();
113  const surfaceVectorField& Sf = mesh.Sf();
114 
115  // Cell gradient of alpha
116  const volVectorField gradAlpha(fvc::grad(alpha1_, "nHat"));
117 
118  // Interpolated face-gradient of alpha
119  surfaceVectorField gradAlphaf(fvc::interpolate(gradAlpha));
120 
121  // gradAlphaf -=
122  // (mesh.Sf()/mesh.magSf())
123  // *(fvc::snGrad(alpha1_) - (mesh.Sf() & gradAlphaf)/mesh.magSf());
124 
125  // Face unit interface normal
126  surfaceVectorField nHatfv(gradAlphaf/(mag(gradAlphaf) + deltaN_));
127  // surfaceVectorField nHatfv
128  // (
129  // (gradAlphaf + deltaN_*vector(0, 0, 1)
130  // *sign(gradAlphaf.component(vector::Z)))/(mag(gradAlphaf) + deltaN_)
131  // );
132  correctContactAngle(nHatfv.boundaryFieldRef(), gradAlphaf.boundaryField());
133 
134  // Face unit interface normal flux
135  nHatf_ = nHatfv & Sf;
136 
137  // Simple expression for curvature
138  K_ = -fvc::div(nHatf_);
139 
140  // Complex expression for curvature.
141  // Correction is formally zero but numerically non-zero.
142  /*
143  volVectorField nHat(gradAlpha/(mag(gradAlpha) + deltaN_));
144  forAll(nHat.boundaryField(), patchi)
145  {
146  nHat.boundaryField()[patchi] = nHatfv.boundaryField()[patchi];
147  }
148 
149  K_ = -fvc::div(nHatf_) + (nHat & fvc::grad(nHatfv) & nHat);
150  */
151 }
152 
153 
154 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
155 
157 (
158  volScalarField& alpha1,
159  volScalarField& alpha2,
160  const volVectorField& U,
161  const IOdictionary& dict
162 )
163 :
164  phasePropertiesDict_(dict),
165 
166  sigmaPtr_(surfaceTensionModel::New(dict, alpha1.mesh())),
167 
168  deltaN_
169  (
170  "deltaN",
171  1e-8/pow(average(alpha1.mesh().V()), 1.0/3.0)
172  ),
173 
174  alpha1_(alpha1),
175  alpha2_(alpha2),
176  U_(U),
177 
178  nHatf_
179  (
180  IOobject
181  (
182  "nHatf",
183  alpha1_.time().timeName(),
184  alpha1_.mesh()
185  ),
186  alpha1_.mesh(),
188  ),
189 
190  K_
191  (
192  IOobject
193  (
194  "interfaceProperties:K",
195  alpha1_.time().timeName(),
196  alpha1_.mesh()
197  ),
198  alpha1_.mesh(),
200  )
201 {
202  calculateK();
203 }
204 
205 
206 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
207 
209 {
210  const volVectorField gradAlpha(fvc::grad(alpha1_));
211 
212  return volVectorField::New("n", gradAlpha/(mag(gradAlpha) + deltaN_));
213 }
214 
215 
218 {
219  return sigmaPtr_->sigma()*K_;
220 }
221 
222 
225 {
226  return fvc::interpolate(sigmaK())*fvc::snGrad(alpha1_);
227 }
228 
229 
232 {
233  return pos0(alpha1_ - 0.01)*pos0(0.99 - alpha1_);
234 }
235 
236 
239 {
240  const fvMesh& mesh = alpha1_.mesh();
241 
243  (
245  (
246  "A",
247  mesh,
249  )
250  );
251  volScalarField::Internal& A = tA.ref();
252 
253  const surfaceVectorField& Sf = mesh.Sf();
254  const labelUList& own = mesh.owner();
255  const labelUList& nei = mesh.neighbour();
256 
257  const surfaceScalarField alphaf(fvc::interpolate(alpha1_));
258  const volVectorField::Internal n(this->n());
259 
260  const scalarField& ialpha = alpha1_;
261  const scalarField& ialphaf = alphaf;
262  scalarField sumnSf(mesh.nCells(), 0);
263 
264  forAll(own, facei)
265  {
266  {
267  const scalar nSf(mag(n[own[facei]] & Sf[facei]));
268  A[own[facei]] += nSf*(ialphaf[facei] - ialpha[own[facei]]);
269  sumnSf[own[facei]] += nSf;
270  }
271  {
272  const scalar nSf(mag(n[nei[facei]] & Sf[facei]));
273  A[nei[facei]] += nSf*(ialphaf[facei] - ialpha[nei[facei]]);
274  sumnSf[nei[facei]] += nSf;
275  }
276  }
277 
278  forAll(mesh.boundary(), patchi)
279  {
280  const labelUList& own = mesh.boundary()[patchi].faceCells();
281  const fvsPatchScalarField& palphaf = alphaf.boundaryField()[patchi];
282 
283  forAll(mesh.boundary()[patchi], facei)
284  {
285  const scalar nSf(mag(n[own[facei]] & Sf[facei]));
286  A[own[facei]] += nSf*(palphaf[facei] - ialpha[own[facei]]);
287  sumnSf[own[facei]] += nSf;
288  }
289  }
290 
291  scalarField& a = A.field();
292  forAll(a, i)
293  {
294  if (sumnSf[i] > small)
295  {
296  a[i] = 2*mag(a[i])/sumnSf[i];
297  }
298  else
299  {
300  a[i] = 0;
301  }
302  }
303 
304  return tA;
305 }
306 
307 
309 {
310  calculateK();
311 }
312 
313 
315 {
316  sigmaPtr_->readDict(phasePropertiesDict_);
317 
318  return true;
319 }
320 
321 
322 // ************************************************************************* //
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 surfaceVectorField & Sf() const
Return cell face area vectors.
const dimensionSet dimArea
static tmp< DimensionedField< Type, GeoMesh > > New(const word &name, const Mesh &mesh, const dimensionSet &)
Return a temporary field constructed from name, mesh.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
void correctContactAngle(const volScalarField &alpha1, const volScalarField &alpha2, const volVectorField::Boundary &Ubf, const dimensionedScalar &deltaN, surfaceVectorField::Boundary &nHatbf)
Correct the contact angle for the two volume fraction fields.
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
Unit conversion functions.
label nCells() const
static tmp< GeometricField< vector, fvPatchField, volMesh > > New(const word &name, const Internal &, const PtrList< fvPatchField< vector >> &)
Return a temporary field constructed from name,.
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 phaseProperties dictionary.
Calculate the snGrad of the given volField.
GeometricBoundaryField< vector, fvsPatchField, surfaceMesh > Boundary
Type of the boundary field.
const dimensionSet dimless
dimensionedScalar det(const dimensionedSphericalTensor &dt)
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:59
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
const labelUList & neighbour() const
Internal face neighbour.
Definition: fvMesh.H:423
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)
tmp< volScalarField::Internal > fraction() const
Interface fraction in a cell.
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: HashTable.H:60
Calculate the divergence of the given field.
dimensionedScalar pos0(const dimensionedScalar &ds)
const Mesh & mesh() const
Return mesh.
defineTypeNameAndDebug(combustionModel, 0)
tmp< volScalarField > sigmaK() const
const Field< Type > & field() 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.
const labelUList & owner() const
Internal face owner.
Definition: fvMesh.H:417
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
label patchi
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:95
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
interfaceProperties(volScalarField &alpha1, volScalarField &alpha2, const volVectorField &U, const IOdictionary &)
Construct from volume fraction field gamma and IOdictionary.
tmp< volVectorField > n() const
dimensioned< scalar > mag(const dimensioned< Type > &)
label n
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:98
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:65
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:45
Namespace for OpenFOAM.
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:800