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