correctContactAngle.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) 2022 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 "correctContactAngle.H"
27 #include "unitConversion.H"
28 
29 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
30 
32 (
33  const volScalarField& a1,
34  const volScalarField& a2,
35  const volVectorField::Boundary& Ubf,
36  const dimensionedScalar& deltaN,
37  surfaceVectorField::Boundary& nHatbf
38 )
39 {
40  const volScalarField::Boundary& a1bf = a1.boundaryField();
41  const volScalarField::Boundary& a2bf = a2.boundaryField();
42 
43  const fvBoundaryMesh& boundary = a1.mesh().boundary();
44 
45  forAll(boundary, patchi)
46  {
47  const fvPatchScalarField& a1p = a1bf[patchi];
48  const fvPatchScalarField& a2p = a2bf[patchi];
49 
50  const bool a1pIsCa = isA<alphaContactAngleFvPatchScalarField>(a1p);
51  const bool a2pIsCa = isA<alphaContactAngleFvPatchScalarField>(a2p);
52 
53  if (a1pIsCa || a2pIsCa)
54  {
55  if (a1pIsCa != a2pIsCa)
56  {
58  << alphaContactAngleFvPatchScalarField::typeName
59  << " boundary condition specified on patch "
60  << boundary[patchi].name() << " for "
61  << (a1pIsCa ? a1.name() : a2.name()) << " but not for "
62  << (a2pIsCa ? a1.name() : a2.name())
63  << exit(FatalError);
64  }
65 
67  refCast<const alphaContactAngleFvPatchScalarField>(a1p);
69  refCast<const alphaContactAngleFvPatchScalarField>(a2p);
70 
71  const bool a1caHasProps = a1ca.thetaProps().found(a2.group());
72  const bool a2caHasProps = a2ca.thetaProps().found(a1.group());
73 
74  if (!a1caHasProps && !a2caHasProps)
75  {
77  << "Neither "
78  << alphaContactAngleFvPatchScalarField::typeName
79  << " boundary condition specified on patch "
80  << boundary[patchi].name()
81  << " for " << a1.name() << " and " << a2.name()
82  << " contains properties for the other phase"
83  << exit(FatalError);
84  }
85 
86  if
87  (
88  a1caHasProps && a2caHasProps
89  && a1ca.thetaProps()[a2.group()]
90  != a2ca.thetaProps()[a1.group()].reversed()
91  )
92  {
94  << "The "
95  << alphaContactAngleFvPatchScalarField::typeName
96  << " boundary conditions specified on patch "
97  << boundary[patchi].name()
98  << " for " << a1.name() << " and " << a2.name()
99  << " contain inconsistent properties"
100  << exit(FatalError);
101  }
102 
104  tp = a1caHasProps
105  ? a1ca.thetaProps()[a2.group()]
106  : a2ca.thetaProps()[a1.group()].reversed();
107 
108  const vectorField np(a1.mesh().boundary()[patchi].nf());
109 
110  vectorField& nHatp = nHatbf[patchi];
111 
112  // Calculate the contact angle
113  scalarField theta(np.size(), degToRad(tp.theta0()));
114 
115  // Calculate the dynamic contact angle if required
116  if (tp.dynamic())
117  {
118  const scalar uTheta = tp.uTheta();
119  const scalar thetaA = degToRad(tp.thetaA());
120  const scalar thetaR = degToRad(tp.thetaR());
121 
122  // Calculated the component of the velocity parallel to the wall
123  vectorField Uwall
124  (
125  Ubf[patchi].patchInternalField() - Ubf[patchi]
126  );
127  Uwall -= (np & Uwall)*np;
128 
129  // Find the direction of the interface parallel to the wall
130  vectorField nWall(nHatp - (np & nHatp)*np);
131  nWall /= (mag(nWall) + small);
132 
133  // Calculate Uwall resolved normal to the interface parallel to
134  // the interface
135  const scalarField uwall(nWall & Uwall);
136 
137  theta += (thetaA - thetaR)*tanh(uwall/uTheta);
138  }
139 
140  // Reset nHatp to correspond to the contact angle
141  const scalarField a12(nHatp & np);
142  const scalarField b1(cos(theta));
143  scalarField b2(nHatp.size());
144  forAll(b2, facei)
145  {
146  b2[facei] = cos(acos(a12[facei]) - theta[facei]);
147  }
148  const scalarField det(1 - a12*a12);
149  const scalarField a((b1 - a12*b2)/det);
150  const scalarField b((b2 - a12*b1)/det);
151 
152  nHatp = a*np + b*nHatp;
153  nHatp /= (mag(nHatp) + deltaN.value());
154  }
155  }
156 }
157 
158 
159 // ************************************************************************* //
static word group(const word &name)
Return group (extension part of name)
Definition: IOobject.C:134
dimensionedScalar tanh(const dimensionedScalar &ds)
dimensionedScalar acos(const dimensionedScalar &ds)
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
const word & name() const
Return name.
Definition: IOobject.H:315
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
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.
Unit conversion functions.
const dimensionedScalar b
Wien displacement law constant: default SI units: [m K].
Definition: createFields.H:27
scalar thetaR() const
Return the limiting receding contact angle.
scalar theta0() const
Return the equilibrium contact angle theta0.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
dimensionedScalar det(const dimensionedSphericalTensor &dt)
const HashTable< contactAngleProperties > & thetaProps() const
Return the contact angle properties.
scalar degToRad(const scalar deg)
Conversion from degrees to radians.
dimensionedScalar cos(const dimensionedScalar &ds)
const Type & value() const
Return const reference to value.
faceListList boundary(nPatches)
scalar uTheta() const
Return the dynamic contact angle velocity scale.
const Mesh & mesh() const
Return mesh.
label patchi
scalar thetaA() const
Return the limiting advancing contact angle.
Foam::fvBoundaryMesh.
dimensioned< scalar > mag(const dimensioned< Type > &)
Contact-angle boundary condition for multi-phase interface-capturing simulations. ...