Go to the documentation of this file.
38 namespace regionModels
40 namespace areaSurfaceFilmModels
49 void contactAngleForce::initialise()
56 contactAngleForce::contactAngleForce
64 Ccf_(coeffDict_.
get<scalar>(
"Ccf")),
69 typeName +
":fContactForceMask",
70 film.primaryMesh().time().
timeName(),
99 typeName +
":contactForce",
101 film().primaryMesh(),
129 const label faceO = own[edgei];
130 const label faceN = nbr[edgei];
133 if ((talpha()[faceO] > 0.5) && (talpha()[faceN] < 0.5))
137 else if ((talpha()[faceO] < 0.5) && (talpha()[faceN] > 0.5))
142 if (facei != -1 && mask_[facei] > 0.5)
147 gradAlpha[facei]/(
mag(gradAlpha[facei]) + ROOTVSMALL)
152 Ccf_*
n*
sigma[facei]*(1 - cosTheta)/invDx/
rhof[facei];
164 label face0 = faces[edgei];
170 force /= magSf.
field();
183 tfvm.ref() += tForce;
Defines the attributes of an object for which implicit objectRegistry management is supported,...
faPatchField<Type> abstract base class. This class gives a fat-interface to all derived classes cover...
A class for handling words, derived from Foam::string.
faMatrix< vector > faVectorMatrix
const fvMesh & primaryMesh() const
A class for managing temporary objects.
static constexpr const zero Zero
const Field< Type > & field() const
const dimensionSet dimDensity
static word timeName(const scalar t, const int precision=precision_)
tmp< areaScalarField > alpha() const
tmp< GeometricField< typename outerProduct< vector, Type >::type, faPatchField, areaMesh >> grad(const GeometricField< Type, faePatchField, edgeMesh > &ssf)
Unit conversion functions.
Base class for film (stress-based) force models.
const dimensionSet dimForce
virtual bool writeTime() const
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
const labelUList & edgeFaces() const
const dimensionSet dimArea(sqr(dimLength))
Generic templated field type.
const DimensionedField< scalar, areaMesh > & S() const
const edgeScalarField & deltaCoeffs() const
const liquidFilmBase & film() const
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Macros for easy insertion into run-time selection tables.
const faMesh & regionMesh() const
Generic dimensioned Type class.
const labelUList & neighbour() const
constexpr scalar degToRad(const scalar deg) noexcept
const labelUList & owner() const
virtual const areaScalarField & rho() const =0
defineTypeNameAndDebug(kinematicThinFilm, 0)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
const Time & time() const
virtual const areaScalarField & sigma() const =0
dimensionedScalar sigma("sigma", dimMass/sqr(dimTime), transportProperties)
surfaceScalarField rhof(fvc::interpolate(rho, "div(phi,rho)"))
Generic GeometricField class.
const faPatch & patch() const
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const dimensionSet dimless
dimensionedScalar cos(const dimensionedScalar &ds)