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-2022 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
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
36namespace Foam
37{
38namespace regionModels
39{
40namespace areaSurfaceFilmModels
41{
42
43// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
44
46
47// * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
48
49void contactAngleForce::initialise()
50{
51}
52
53
54// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
55
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 (
68 (
69 typeName + ":fContactForceMask",
70 film.primaryMesh().time().timeName(),
71 film.primaryMesh(),
72 IOobject::NO_READ,
73 IOobject::NO_WRITE
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 (
96 (
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
157 force /= magSf.field();
158
159 if (film().regionMesh().time().writeTime())
160 {
161 tForce().write();
162 gradAlpha.write();
163 }
164
166 (
168 );
169
170 tfvm.ref() += tForce;
171
172 return tfvm;
173}
174
175
176// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
177
178} // End namespace areaSurfaceFilmModels
179} // End namespace regionModels
180} // End namespace Foam
181
182// ************************************************************************* //
label n
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const Field< Type > & field() const
Return field.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
bool writeTime() const noexcept
True if this is a write time.
Definition: TimeStateI.H:67
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
const edgeScalarField & deltaCoeffs() const
Return reference to difference factors array.
const Time & time() const
Return reference to time.
Definition: faMesh.C:673
const DimensionedField< scalar, areaMesh > & S() const
Return face areas.
Definition: faMesh.C:780
const labelUList & owner() const
Internal face owner.
Definition: faMesh.H:719
const labelUList & neighbour() const
Internal face neighbour.
Definition: faMesh.H:725
virtual bool write(const bool valid=true) const
Write using setting from DB.
Base-class for film contact angle force models.
virtual tmp< areaScalarField > theta() const =0
Return the contact angle field.
const liquidFilmBase & film() const
Return const access to the film surface film model.
Base class for film (stress-based) force models.
Definition: force.H:59
virtual const areaScalarField & sigma() const =0
Access const reference sigma.
tmp< areaScalarField > alpha() const
Wet indicator using h0.
virtual const areaScalarField & rho() const =0
Access const reference rho.
const faMesh & regionMesh() const
Return the region mesh database.
A class for managing temporary objects.
Definition: tmp.H:65
T & ref() const
Definition: tmpI.H:227
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
U
Definition: pEqn.H:72
word timeName
Definition: getTimeIndex.H:3
surfaceScalarField rhof(fvc::interpolate(rho, "div(phi,rho)"))
tmp< GeometricField< typename outerProduct< vector, Type >::type, faPatchField, areaMesh > > grad(const GeometricField< Type, faePatchField, edgeMesh > &ssf)
Definition: facGrad.C:56
Namespace for OpenFOAM.
const dimensionSet dimless
Dimensionless.
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:59
constexpr scalar degToRad() noexcept
Multiplication factor for degrees to radians conversion.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
faMatrix< vector > faVectorMatrix
Definition: faMatrices.H:53
const dimensionSet dimForce
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
const dimensionSet dimDensity
dimensionedScalar cos(const dimensionedScalar &ds)
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
Unit conversion functions.