36 namespace regionModels
38 namespace surfaceFilmModels
47 void contactAngleForce::initialise()
54 if (zeroForcePatches.size())
59 Info<<
" Assigning zero contact force within " << dLim
60 <<
" of patches:" <<
endl;
66 label patchi = iter.key();
102 force(typeName, film, dict),
160 const label cellO = own[facei];
161 const label cellN = nbr[facei];
164 if ((coverage[cellO] > 0.5) && (coverage[cellN] < 0.5))
168 else if ((coverage[cellO] < 0.5) && (coverage[cellN] > 0.5))
173 if (celli != -1 && mask_[celli] > 0.5)
177 gradCoverage[celli]/(
mag(gradCoverage[celli]) + rootVSmall);
178 const scalar cosTheta =
cos(
degToRad(theta[celli]));
179 force[celli] += Ccf_*n*sigma[celli]*(1 - cosTheta)/invDx;
197 if (maskPf[facei] > 0.5)
199 const label cellO = faceCells[facei];
201 if ((coverage[cellO] > 0.5) && (coveragePf[facei] < 0.5))
205 /(
mag(gradCoverage[cellO]) + rootVSmall);
207 const scalar cosTheta =
cos(
degToRad(thetaPf[facei]));
209 Ccf_*n*sigmaPf[facei]*(1 - cosTheta)/invDx[facei];
228 tfvm.
ref() += tForce;
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
#define forAll(list, i)
Loop across all elements in list.
FvWallInfoData< WallInfo, label > label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Base class for surface film models.
Base class for film (stress-based) force models.
A list of keyword definitions, which are a keyword followed by any number of values (e...
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
T & ref() const
Return non-const reference or generate a fatal error.
Unit conversion functions.
static tmp< GeometricField< vector, fvPatchField, volMesh > > New(const word &name, const Internal &, const PtrList< fvPatchField< vector >> &)
Return a temporary field constructed from name,.
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
const dimensionSet dimless
const Time & time() const
Return the reference to the time database.
const Time & time() const
Return the top-level database.
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
const dictionary coeffDict_
Coefficients dictionary.
Macros for easy insertion into run-time selection tables.
scalar degToRad(const scalar deg)
Conversion from degrees to radians.
const dimensionSet dimLength
const labelUList & neighbour() const
Internal face neighbour.
bool writeTime() const
Return true if this is a write time.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
HashSet< label, Hash< label > > labelHashSet
A HashSet with label keys.
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m^2/K^4].
dimensionedScalar cos(const dimensionedScalar &ds)
Calculate the gradient of the given field.
static word timeName(const scalar, const int precision=curPrecision_)
Return time name of given scalar time.
A class for handling words, derived from string.
bool isCoupledPatch(const label regionPatchi) const
Return true if patchi on the local region is a coupled.
const fvMesh & regionMesh() const
Return the region mesh database.
const dimensionSet dimForce
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
const fvPatch & patch() const
Return patch.
dimensionedScalar pos0(const dimensionedScalar &ds)
const labelUList & owner() const
Internal face owner.
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
force(surfaceFilmRegionModel &film)
Construct null.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
fvMatrix< vector > fvVectorMatrix
const dimensionSet dimVolume
surfaceFilmRegionModel & filmModel_
Reference to the film surface film model.
dimensioned< scalar > mag(const dimensioned< Type > &)
virtual const volScalarField & coverage() const =0
Return the film coverage, 1 = covered, 0 = uncovered / [].
static word modelName(Name name, const word &model)
Return the name of the object within the given model.
List< wordRe > wordReList
A List of wordRe (word or regular expression)
A class for managing temporary objects.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
defineTypeNameAndDebug(kinematicSingleLayer, 0)
virtual tmp< volScalarField > sigma() const =0
Return the film surface tension [N/m].
const surfaceScalarField & deltaCoeffs() const
Return reference to cell-centre difference coefficients.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.