Go to the documentation of this file.
37 namespace regionModels
39 namespace surfaceFilmModels
53 if (zeroForcePatches.
size())
58 Info<<
" Assigning zero contact force within " << dLim
59 <<
" of patches:" <<
endl;
65 label patchI = iter.key();
102 rndGen_(
label(0), -1),
107 coeffDict_.subDict(
"contactAngleDistribution"),
115 typeName +
":contactForceMask",
116 owner_.time().timeName(),
145 typeName +
":contactForce",
170 const label cellO = own[faceI];
171 const label cellN = nbr[faceI];
174 if ((
alpha[cellO] > 0.5) && (
alpha[cellN] < 0.5))
178 else if ((
alpha[cellO] < 0.5) && (
alpha[cellN] > 0.5))
183 if (cellI != -1 &&
mask_[cellI] > 0.5)
187 gradAlpha[cellI]/(
mag(gradAlpha[cellI]) + ROOTVSMALL);
195 if (!
owner().isCoupledPatch(patchI))
204 if (maskf[faceI] > 0.5)
206 label cellO = faceCells[faceI];
208 if ((
alpha[cellO] > 0.5) && (alphaf[faceI] < 0.5))
212 /(
mag(gradAlpha[cellO]) + ROOTVSMALL);
215 Ccf_*
n*
sigma[cellO]*(1.0 - theta)/invDx[faceI];
virtual const volScalarField & sigma() const =0
Return the film surface tension [N/m].
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
bool outputTime() const
Return true if this is an output time (primary or secondary)
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
#define forAll(list, i)
Loop across all elements in list.
A class for managing temporary objects.
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
const dictionary coeffDict_
Coefficients dictionary.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
const Time & time() const
Return the reference to the time database.
Unit conversion functions.
GeometricBoundaryField & boundaryField()
Return reference to GeometricBoundaryField.
static autoPtr< distributionModel > New(const dictionary &dict, cachedRandom &rndGen)
Selector.
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
const surfaceFilmModel & owner() const
Return const access to the owner surface film model.
Ostream & endl(Ostream &os)
Add newline and flush stream.
dimensioned< scalar > mag(const dimensioned< Type > &)
const dimensionSet dimForce
const surfaceScalarField & deltaCoeffs() const
Return reference to cell-centre difference coefficients.
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
fvMatrix< vector > fvVectorMatrix
const dimensionSet dimArea(sqr(dimLength))
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Pre-declare SubField and related Field type.
virtual const volScalarField & magSf() const
Return the face area magnitudes / [m2].
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const word & name() const
Return name.
A list of keyword definitions, which are a keyword followed by any number of values (e....
const labelUList & neighbour() const
Internal face neighbour.
Macros for easy insertion into run-time selection tables.
Generic dimensioned Type class.
const labelUList & owner() const
Internal face owner.
defineTypeNameAndDebug(kinematicSingleLayer, 0)
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if succesful.
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Base class for film (stress-based) force models.
Calculate the gradient of the given field.
const scalarField & deltaCoeffs() const
Return the face - cell distance coeffient.
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
virtual bool correct(volScalarField &y)
Correct the given distance-to-patch field.
surfaceFilmModel & owner_
Reference to the owner surface film model.
addToRunTimeSelectionTable(surfaceFilmModel, kinematicSingleLayer, mesh)
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m2/K4].
virtual const labelUList & faceCells() const
Return faceCells.
const Time & time() const
Return the top-level database.
const fvPatch & patch() const
Return patch.
Base class for surface film models.
void size(const label)
Override size to be inconsistent with allocated storage.
labelHashSet patchSet(const UList< wordRe > &patchNames, const bool warnNotFound=true, const bool usePatchGroups=true) const
Return the set of patch IDs corresponding to the given names.
const dimensionSet dimVolume(pow3(dimLength))
Fast topological mesh-wave method for calculating the distance to nearest patch for all cells and bou...
Generic GeometricField class.
virtual const volScalarField & alpha() const =0
Return the film coverage, 1 = covered, 0 = uncovered / [].
scalar degToRad(const scalar deg)
Conversion from degrees to radians.
dimensionedScalar cos(const dimensionedScalar &ds)
const fvMesh & regionMesh() const
Return the region mesh database.
dimensionedScalar pos(const dimensionedScalar &ds)