Go to the documentation of this file.
39 namespace regionModels
41 namespace surfaceFilmModels
69 c_(
"c",
pow(
dimTime, d_.value() - scalar(1)), coeffDict_),
71 muInf_(
"muInf", mu0_.dimensions(), coeffDict_),
91 mu_.correctBoundaryConditions();
138 "thixotropicViscosity:filmMass",
143 "thixotropicViscosity:deltaMass",
153 a_*
pow((1.0 - lambda_), b_)
164 mu_ = muInf_/(
sqr(1.0 - K_*lambda_) + ROOTVSMALL);
165 mu_.correctBoundaryConditions();
thixotropicViscosity(const thixotropicViscosity &)
Disallow default bitwise copy construct.
const dimensionSet dimPressure
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)
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
dimensionedScalar deltaT() const
Return time step.
const dimensionedScalar mu
Atomic mass unit.
virtual ~thixotropicViscosity()
Destructor.
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
virtual const volVectorField & U() const
Return the film velocity [m/s].
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &vf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
Calculate the divergence of the given field.
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
dimensioned< scalar > mag(const dimensioned< Type > &)
const dimensionedScalar & deltaSmall() const
Return small delta.
Calculate the matrix for the divergence of the given field and flux.
tmp< fvMatrix< Type > > Sp(const DimensionedField< scalar, volMesh > &, const GeometricField< Type, fvPatchField, volMesh > &)
SolverPerformance< Type > solve(const dictionary &)
Solve segregated or coupled returning the solution statistics.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
tmp< volScalarField > deltaMass() const
Return the change in film mass due to sources/sinks.
virtual void correct(const volScalarField &p, const volScalarField &T)
Correct.
void relax(const scalar alpha)
Relax matrix (for steady-state solution).
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
tmp< fvMatrix< Type > > SuSp(const DimensionedField< scalar, volMesh > &, const GeometricField< Type, fvPatchField, volMesh > &)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Thixotropic viscosity model based on the evolution of the structural parameter :
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
virtual const surfaceScalarField & phi() const
Return the film flux [kg.m/s].
A list of keyword definitions, which are a keyword followed by any number of values (e....
Macros for easy insertion into run-time selection tables.
Generic dimensioned Type class.
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Calculate the matrix for implicit and explicit sources.
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
defineTypeNameAndDebug(kinematicSingleLayer, 0)
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
dimensionedScalar sqrt(const dimensionedScalar &ds)
virtual const volVectorField & Uw() const
Return the film wall velocity [m/s].
const volScalarField & delta() const
Return const access to the film thickness / [m].
const volScalarField & alpha() const
Return the film coverage, 1 = covered, 0 = uncovered / [].
addToRunTimeSelectionTable(surfaceFilmModel, kinematicSingleLayer, mesh)
tmp< volScalarField > netMass() const
Return the net film mass available over the next integration.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
const Time & time() const
Return the top-level database.
Base class for surface film models.
Base class for surface film viscosity models.
virtual const volScalarField & deltaRho() const
Return the film thickness*density (helper field) [kg/m3].
Generic GeometricField class.
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Calulate the matrix for the first temporal derivative.
const fvMesh & regionMesh() const
Return the region mesh database.