contactAngleForce.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) 2020-2023 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "contactAngleForce.H"
30 #include "unitConversion.H"
31 
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace regionModels
37 {
38 namespace areaSurfaceFilmModels
39 {
40 
41 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
42 
44 
45 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
46 
47 contactAngleForce::contactAngleForce
48 (
49  const word& typeName,
50  liquidFilmBase& film,
51  const dictionary& dict
52 )
53 :
54  force(typeName, film, dict),
55  Ccf_(coeffDict_.get<scalar>("Ccf")),
56  hCrit_(coeffDict_.getOrDefault<scalar>("hCrit", GREAT))
57 {}
58 
59 
60 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
61 
63 {
64  auto tforce = tmp<areaVectorField>::New
65  (
66  IOobject
67  (
68  IOobject::scopedName(typeName, "contactForce"),
69  film().regionMesh().time().timeName(),
70  film().regionMesh().thisDb(),
74  ),
75  film().regionMesh(),
77  );
78  vectorField& force = tforce.ref().primitiveFieldRef();
79 
80  const labelUList& own = film().regionMesh().owner();
81  const labelUList& nei = film().regionMesh().neighbour();
82 
84  const scalarField& magSff = magSf.field();
85 
86  const edgeScalarField& invDx = film().regionMesh().deltaCoeffs();
87 
88  const areaScalarField& sigma = film().sigma();
89  const areaScalarField& mu = film().mu();
90  const areaScalarField& rhof = film().rho();
91  const areaVectorField& Uf = film().Uf();
92  const areaScalarField& hf = film().h();
93 
94  tmp<areaScalarField> talpha = film().alpha();
95  const areaScalarField& alpha = talpha();
96 
97  tmp<areaScalarField> ttheta = theta();
98  const areaScalarField& theta = ttheta();
99 
100  tmp<areaVectorField> tgradAlpha = fac::grad(alpha);
101  const areaVectorField& gradAlpha = tgradAlpha();
102 
103  forAll(nei, edgei)
104  {
105  const label faceO = own[edgei];
106  const label faceN = nei[edgei];
107 
108  label facei = -1;
109  if ((alpha[faceO] > 0.5) && (alpha[faceN] < 0.5))
110  {
111  facei = faceO;
112  }
113  else if ((alpha[faceO] < 0.5) && (alpha[faceN] > 0.5))
114  {
115  facei = faceN;
116  }
117 
118  if (facei != -1)
119  {
120  const vector n(normalised(gradAlpha[facei]));
121  const scalar cosTheta = cos(degToRad(theta[facei]));
122 
123  // (MHDX:Eq. 13)
124  force[facei] +=
125  Ccf_*n*sigma[facei]*(1 - cosTheta)
126  /invDx[edgei]/rhof[facei]/magSff[facei];
127 
128  // (NDPC:Eq. 11)
129  if (hf[facei] > hCrit_)
130  {
131  force[facei] -= mu[facei]*Uf[facei]/hCrit_;
132  }
133  }
134  }
135 
136  for (const faPatchScalarField& sigmaBf : sigma.boundaryField())
137  {
138  const faPatch& p = sigmaBf.patch();
139 
140  if (!p.coupled())
141  {
142  const labelUList& faces = p.edgeFaces();
143 
144  forAll(sigmaBf, edgei)
145  {
146  const label face0 = faces[edgei];
147  force[face0] = Zero;
148  }
149  }
150  }
151 
152  if (film().regionMesh().time().writeTime())
153  {
154  tforce().write();
155  gradAlpha.write();
156  }
157 
159 
160  tfvm.ref() += tforce;
161 
162  return tfvm;
163 }
164 
165 
166 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
167 
168 } // End namespace areaSurfaceFilmModels
169 } // End namespace regionModels
170 } // End namespace Foam
171 
172 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
faPatchField< scalar > faPatchScalarField
Base-class for film contact angle force models.
dictionary dict
const edgeScalarField & deltaCoeffs() const
Return reference to difference factors array.
virtual tmp< areaScalarField > theta() const =0
Return the contact angle field.
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
virtual const areaScalarField & sigma() const =0
Access const reference sigma.
Unit conversion functions.
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
const DimensionedField< scalar, areaMesh > & S() const
Return face areas.
Definition: faMesh.C:1133
Base class for film (stress-based) force models.
Definition: force.H:52
Ignore writing from objectRegistry::writeObject()
quaternion normalised(const quaternion &q)
Return the normalised (unit) quaternion of the given quaternion.
Definition: quaternionI.H:674
defineTypeNameAndDebug(kinematicThinFilm, 0)
dimensionedScalar sigma("sigma", dimMass/sqr(dimTime), transportProperties)
Macros for easy insertion into run-time selection tables.
UList< label > labelUList
A UList of labels.
Definition: UList.H:78
virtual bool writeTime() const
Flag to indicate when to write a property.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
const labelUList & neighbour() const
Internal face neighbour.
Definition: faMesh.H:1332
word timeName
Definition: getTimeIndex.H:3
static word scopedName(const std::string &scope, const word &name)
Create scope:name or scope_name string.
Definition: IOobjectI.H:40
dimensionedScalar cos(const dimensionedScalar &ds)
virtual tmp< faVectorMatrix > correct(areaVectorField &U)
Correct.
A class for handling words, derived from Foam::string.
Definition: word.H:63
const areaScalarField & h() const noexcept
Access const reference h.
virtual bool write(const bool writeOnProc=true) const
Write using setting from DB.
const faMesh & regionMesh() const
Return the region mesh database.
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition: tmp.H:206
autoPtr< surfaceVectorField > Uf
const dimensionSet dimForce
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
const dimensionSet dimDensity
const labelUList & owner() const
Internal face owner.
Definition: faMesh.H:1324
const Field< Type > & field() const noexcept
Return const-reference to the field values.
const dimensionedScalar mu
Atomic mass unit.
U
Definition: pEqn.H:72
virtual const areaScalarField & rho() const =0
Access const reference rho.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Nothing to be read.
const areaVectorField & Uf() const noexcept
Access const reference Uf.
label n
virtual const areaScalarField & mu() const =0
Access const reference mu.
volScalarField & p
A class for managing temporary objects.
Definition: HashPtrTable.H:50
tmp< areaScalarField > alpha() const
Wet indicator using h0.
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
constexpr scalar degToRad(const scalar deg) noexcept
Conversion from degrees to radians.
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:180
surfaceScalarField rhof(fvc::interpolate(rho, "div(phi,rho)"))
Do not request registration (bool: false)
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:57
Namespace for OpenFOAM.
tmp< GeometricField< typename outerProduct< vector, Type >::type, faPatchField, areaMesh >> grad(const GeometricField< Type, faePatchField, edgeMesh > &ssf)
Definition: facGrad.C:51
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127
const liquidFilmBase & film() const
Return const access to the film surface film model.