Go to the documentation of this file.
39 namespace regionModels
41 namespace surfaceFilmModels
50 void contactAngleForce::initialise()
52 const wordRes zeroForcePatches
57 if (zeroForcePatches.size())
60 const scalar dLim =
coeffDict_.
get<scalar>(
"zeroForceDistance");
62 Info<<
" Assigning zero contact force within " << dLim
63 <<
" of patches:" <<
endl;
67 for (
const label patchi : patchIDs)
69 Info<<
" " << pbm[patchi].name() <<
endl;
97 contactAngleForce::contactAngleForce
100 surfaceFilmRegionModel& film,
104 force(typeName, film,
dict),
105 Ccf_(coeffDict_.
get<scalar>(
"Ccf")),
110 typeName +
":contactForceMask",
112 filmModel_.regionMesh(),
116 filmModel_.regionMesh(),
140 typeName +
":contactForce",
168 const label cellO = own[facei];
169 const label cellN = nbr[facei];
172 if ((
alpha[cellO] > 0.5) && (
alpha[cellN] < 0.5))
176 else if ((
alpha[cellO] < 0.5) && (
alpha[cellN] > 0.5))
181 if (celli != -1 && mask_[celli] > 0.5)
185 gradAlpha[celli]/(
mag(gradAlpha[celli]) + ROOTVSMALL);
187 force[celli] += Ccf_*
n*
sigma[celli]*(1 - cosTheta)/invDx;
195 const fvPatchField<scalar>& alphaPf =
alpha.boundaryField()[patchi];
196 const fvPatchField<scalar>& maskPf = mask_.
boundaryField()[patchi];
197 const fvPatchField<scalar>& sigmaPf =
sigma.boundaryField()[patchi];
198 const fvPatchField<scalar>& thetaPf =
theta.boundaryField()[patchi];
199 const scalarField& invDx = alphaPf.patch().deltaCoeffs();
200 const labelUList& faceCells = alphaPf.patch().faceCells();
204 if (maskPf[facei] > 0.5)
206 label cellO = faceCells[facei];
208 if ((
alpha[cellO] > 0.5) && (alphaPf[facei] < 0.5))
212 /(
mag(gradAlpha[cellO]) + ROOTVSMALL);
213 const scalar cosTheta =
cos(
degToRad(thetaPf[facei]));
215 Ccf_*
n*sigmaPf[facei]*(1 - cosTheta)/invDx[facei];
229 tmp<fvVectorMatrix> tfvm
234 tfvm.ref() += tForce;
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
A class for handling words, derived from Foam::string.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
A class for managing temporary objects.
static constexpr const zero Zero
const dimensionedScalar alpha
static word timeName(const scalar t, const int precision=precision_)
const dictionary coeffDict_
const Time & time() const
Unit conversion functions.
const polyBoundaryMesh & boundaryMesh() const
Ostream & endl(Ostream &os)
dimensionedScalar pos0(const dimensionedScalar &ds)
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
const dimensionSet dimForce
virtual const surfaceScalarField & deltaCoeffs() const
virtual const volScalarField & sigma() const =0
fvMatrix< vector > fvVectorMatrix
const dimensionSet dimArea(sqr(dimLength))
bool isCoupledPatch(const label regionPatchi) const
Generic templated field type.
Base class for surface film models.
bool writeTime() const noexcept
virtual const volScalarField & magSf() const
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
virtual const volScalarField & alpha() const =0
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
const labelUList & neighbour() const
Macros for easy insertion into run-time selection tables.
Generic dimensioned Type class.
Vector< scalar > vector
A scalar version of the templated Vector.
constexpr scalar degToRad(const scalar deg) noexcept
const labelUList & owner() const
defineTypeNameAndDebug(kinematicSingleLayer, 0)
Base class for film (stress-based) force models.
Calculate the gradient of the given field.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
const scalarField & deltaCoeffs() const
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
virtual const labelUList & faceCells() const
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
const Time & time() const
const fvPatch & patch() const
dimensionedScalar sigma("sigma", dimMass/sqr(dimTime), transportProperties)
UList< label > labelUList
A UList of labels.
HashSet< label, Hash< label > > labelHashSet
A HashSet with label keys and label hasher.
surfaceFilmRegionModel & filmModel_
const dimensionSet dimVolume(pow3(dimLength))
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Generic GeometricField class.
Smooth ATC in cells next to a set of patches supplied by type.
const Boundary & boundaryField() const
const dimensionSet dimless
dimensionedScalar cos(const dimensionedScalar &ds)
const fvMesh & regionMesh() const