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 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 "faCFD.H"
31 #include "unitConversion.H"
33 
34 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38 namespace regionModels
39 {
40 namespace areaSurfaceFilmModels
41 {
42 
43 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
44 
45 defineTypeNameAndDebug(contactAngleForce, 0);
46 
47 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
48 
49 void contactAngleForce::initialise()
50 {
51 }
52 
53 
54 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
55 
56 contactAngleForce::contactAngleForce
57 (
58  const word& typeName,
59  liquidFilmBase& film,
60  const dictionary& dict
61 )
62 :
63  force(typeName, film, dict),
64  Ccf_(coeffDict_.get<scalar>("Ccf")),
65  mask_
66  (
67  IOobject
68  (
69  typeName + ":fContactForceMask",
70  film.primaryMesh().time().timeName(),
71  film.primaryMesh(),
74  ),
75  film.regionMesh(),
76  dimensionedScalar("mask", dimless, 1.0)
77  )
78 {
79  initialise();
80 }
81 
82 
83 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
84 
86 {}
87 
88 
89 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
90 
92 {
94  (
95  new areaVectorField
96  (
97  IOobject
98  (
99  typeName + ":contactForce",
100  film().primaryMesh().time().timeName(),
101  film().primaryMesh(),
104  ),
105  film().regionMesh(),
107  )
108  );
109 
110  vectorField& force = tForce.ref().primitiveFieldRef();
111 
112  const labelUList& own = film().regionMesh().owner();
113  const labelUList& nbr = film().regionMesh().neighbour();
114 
116 
117  tmp<areaScalarField> talpha = film().alpha();
118  const areaScalarField& sigma = film().sigma();
119 
120  const areaScalarField& rhof = film().rho();
121 
122  tmp<areaScalarField> ttheta = theta();
123  const areaScalarField& theta = ttheta();
124 
125  const areaVectorField gradAlpha(fac::grad(talpha()));
126 
127  forAll(nbr, edgei)
128  {
129  const label faceO = own[edgei];
130  const label faceN = nbr[edgei];
131 
132  label facei = -1;
133  if ((talpha()[faceO] > 0.5) && (talpha()[faceN] < 0.5))
134  {
135  facei = faceO;
136  }
137  else if ((talpha()[faceO] < 0.5) && (talpha()[faceN] > 0.5))
138  {
139  facei = faceN;
140  }
141 
142  if (facei != -1 && mask_[facei] > 0.5)
143  {
144  const scalar invDx = film().regionMesh().deltaCoeffs()[edgei];
145  const vector n
146  (
147  gradAlpha[facei]/(mag(gradAlpha[facei]) + ROOTVSMALL)
148  );
149  const scalar cosTheta = cos(degToRad(theta[facei]));
150 
151  force[facei] +=
152  Ccf_*n*sigma[facei]*(1 - cosTheta)/invDx/rhof[facei];
153  }
154  }
155 
156  forAll(sigma.boundaryField(), patchi)
157  {
158  const faPatchField<scalar>& alphaPf = sigma.boundaryField()[patchi];
159  const faPatchField<scalar>& sigmaPf = sigma.boundaryField()[patchi];
160  const labelUList& faces = alphaPf.patch().edgeFaces();
161 
162  forAll(sigmaPf, edgei)
163  {
164  label face0 = faces[edgei];
165  force[face0] = vector::zero;
166  }
167  }
168 
169 
170  force /= magSf.field();
171 
172  if (film().regionMesh().time().writeTime())
173  {
174  tForce().write();
175  gradAlpha.write();
176  }
177 
179  (
181  );
182 
183  tfvm.ref() += tForce;
184 
185  return tfvm;
186 }
187 
188 
189 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
190 
191 } // End namespace areaSurfaceFilmModels
192 } // End namespace regionModels
193 } // End namespace Foam
194 
195 // ************************************************************************* //
Foam::IOobject::NO_WRITE
Definition: IOobject.H:195
Foam::regionModels::areaSurfaceFilmModels::liquidFilmBase
Definition: liquidFilmBase.H:60
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
Foam::faPatchField
faPatchField<Type> abstract base class. This class gives a fat-interface to all derived classes cover...
Definition: areaFieldsFwd.H:50
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
faCFD.H
Foam::faVectorMatrix
faMatrix< vector > faVectorMatrix
Definition: faMatrices.H:53
Foam::regionModels::regionFaModel::primaryMesh
const fvMesh & primaryMesh() const
Return the reference to the primary mesh database.
Definition: regionFaModelI.H:33
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::DimensionedField::field
const Field< Type > & field() const
Return field.
Definition: DimensionedFieldI.H:90
Foam::dimDensity
const dimensionSet dimDensity
Foam::regionModels::areaSurfaceFilmModels::contactAngleForce::~contactAngleForce
virtual ~contactAngleForce()
Destructor.
Definition: contactAngleForce.C:85
Foam::Time::timeName
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:780
Foam::regionModels::areaSurfaceFilmModels::liquidFilmBase::alpha
tmp< areaScalarField > alpha() const
Wet indicator using h0.
Definition: liquidFilmBase.C:432
Foam::fac::grad
tmp< GeometricField< typename outerProduct< vector, Type >::type, faPatchField, areaMesh >> grad(const GeometricField< Type, faePatchField, edgeMesh > &ssf)
Definition: facGrad.C:56
unitConversion.H
Unit conversion functions.
Foam::regionModels::areaSurfaceFilmModels::force
Base class for film (stress-based) force models.
Definition: force.H:56
Foam::dimForce
const dimensionSet dimForce
Foam::regionModels::areaSurfaceFilmModels::filmSubModelBase::writeTime
virtual bool writeTime() const
Flag to indicate when to write a property.
Definition: filmSubModelBase.C:98
Foam::regionModels::areaSurfaceFilmModels::contactAngleForce::correct
virtual tmp< faVectorMatrix > correct(areaVectorField &U)
Correct.
Definition: contactAngleForce.C:91
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::dimensionedVector
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
Definition: dimensionedVector.H:50
Foam::faPatch::edgeFaces
const labelUList & edgeFaces() const
Return edge-face addressing.
Definition: faPatch.C:358
n
label n
Definition: TABSMDCalcMethod2.H:31
contactAngleForce.H
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:227
Foam::dimArea
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:59
Foam::Field< vector >
Foam::faMesh::S
const DimensionedField< scalar, areaMesh > & S() const
Return face areas.
Definition: faMesh.C:632
Foam::edgeInterpolation::deltaCoeffs
const edgeScalarField & deltaCoeffs() const
Return reference to difference factors array.
Definition: edgeInterpolation.C:101
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:42
timeName
word timeName
Definition: getTimeIndex.H:3
meshWavePatchDistMethod.H
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::regionModels::regionFaModel::regionMesh
const faMesh & regionMesh() const
Return the region mesh database.
Definition: regionFaModelI.H:63
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::faMesh::neighbour
const labelUList & neighbour() const
Internal face neighbour.
Definition: faMesh.H:636
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::degToRad
constexpr scalar degToRad(const scalar deg) noexcept
Conversion from degrees to radians.
Definition: unitConversion.H:48
Foam::faMesh::owner
const labelUList & owner() const
Internal face owner.
Definition: faMesh.H:630
Foam::GeometricField::primitiveFieldRef
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field.
Definition: GeometricField.C:766
Foam::regionModels::areaSurfaceFilmModels::liquidFilmBase::rho
virtual const areaScalarField & rho() const =0
Access const reference rho.
U
U
Definition: pEqn.H:72
Foam::regionModels::areaSurfaceFilmModels::defineTypeNameAndDebug
defineTypeNameAndDebug(kinematicThinFilm, 0)
Foam::Vector< scalar >
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::UList< label >
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:280
Foam::VectorSpace< Vector< scalar >, scalar, 3 >::zero
static const Vector< scalar > zero
Definition: VectorSpace.H:115
Foam::regionModels::areaSurfaceFilmModels::liquidFilmBase::sigma
virtual const areaScalarField & sigma() const =0
Access const reference sigma.
sigma
dimensionedScalar sigma("sigma", dimMass/sqr(dimTime), transportProperties)
rhof
surfaceScalarField rhof(fvc::interpolate(rho, "div(phi,rho)"))
Foam::GeometricField< vector, faPatchField, areaMesh >
Foam::IOobject::NO_READ
Definition: IOobject.H:188
Foam::regionModels::areaSurfaceFilmModels::contactAngleForce::theta
virtual tmp< areaScalarField > theta() const =0
Return the contact angle field.
Foam::regionModels::areaSurfaceFilmModels::filmSubModelBase::film
const liquidFilmBase & film() const
Return const access to the film surface film model.
Definition: filmSubModelBaseI.H:39
Foam::faPatchField::patch
const faPatch & patch() const
Return patch.
Definition: faPatchField.H:279
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54
Foam::dimless
const dimensionSet dimless
Dimensionless.
Definition: dimensionSets.C:189
Foam::cos
dimensionedScalar cos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:265