Go to the documentation of this file.
35 #include "surfaceInterpolate.H"
45 namespace regionModels
47 namespace surfaceFilmModels
84 scalarField& invR1 = tinvR1.ref().primitiveFieldRef();
87 const scalar rMin = 1
e-6;
98 const scalar rMax = 1e6;
101 if (
mag(invR1[i]) < 1/rMax)
130 label cellO = own[facei];
131 label cellN = nbr[facei];
133 if (
phi[facei] > phiMax[cellO])
135 phiMax[cellO] =
phi[facei];
136 cosAngle[cellO] = -
gHat_ & nf[facei];
138 if (-
phi[facei] > phiMax[cellN])
140 phiMax[cellN] = -
phi[facei];
141 cosAngle[cellN] = -
gHat_ & -nf[facei];
148 const fvPatch& pp = phip.patch();
149 const labelList& faceCells = pp.faceCells();
153 label celli = faceCells[i];
154 if (phip[i] > phiMax[celli])
156 phiMax[celli] = phip[i];
157 cosAngle[celli] = -
gHat_ & nf[i];
205 mesh.time().timeName(),
211 zeroGradientFvPatchScalarField::typeName
213 volCosAngle.primitiveFieldRef() = cosAngle;
214 volCosAngle.correctBoundaryConditions();
218 return max(
min(cosAngle, scalar(1)), scalar(-1));
224 curvatureSeparation::curvatureSeparation
231 gradNHat_(fvc::
grad(film.nHat())),
232 deltaByR1Min_(coeffDict_.getOrDefault<scalar>(
"deltaByR1Min", 0)),
233 definedPatchRadii_(),
234 magG_(
mag(film.
g().value())),
237 if (
magG_ < ROOTVSMALL)
240 <<
"Acceleration due to gravity must be non-zero"
258 const label patchi = patchIDs[j];
260 if (!uniquePatchIDs.found(patchi))
262 const scalar radius = prIn[i].second();
265 uniquePatchIDs.
insert(patchi);
290 refCast<const kinematicSingleLayer>(this->
film());
304 const scalar Fthreshold = 1
e-10;
311 scalar R1 = 1.0/(invR1[i] + ROOTVSMALL);
312 scalar R2 = R1 +
delta[i];
315 scalar Fi = -
delta[i]*
rho[i]*magSqrU[i]*72.0/60.0*invR1[i];
322 scalar Fs =
sigma[i]/R2;
324 Fnet[i] = Fi +
Fb + Fs;
326 if (Fnet[i] + Fthreshold < 0)
334 massToInject = separated*availableMass;
335 diameterToInject = separated*
delta;
336 availableMass -= separated*availableMass;
347 mesh.time().timeName(),
353 zeroGradientFvPatchScalarField::typeName
355 volFnet.primitiveFieldRef() = Fnet;
356 volFnet.correctBoundaryConditions();
List< label > labelList
A List of labels.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
List< Tuple2< label, scalar > > definedPatchRadii_
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
virtual ~curvatureSeparation()
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
A class for managing temporary objects.
static constexpr const zero Zero
const dimensionSet dimVelocity
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
fvsPatchField< scalar > fvsPatchScalarField
virtual const volVectorField & U() const
tmp< volScalarField > calcInvR1(const volVectorField &U) const
virtual const volScalarField & delta() const =0
const dictionary coeffDict_
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
virtual const volVectorField & U() const =0
Calculate the divergence of the given field.
tmp< GeometricField< typename outerProduct< vector, Type >::type, faPatchField, areaMesh >> grad(const GeometricField< Type, faePatchField, edgeMesh > &ssf)
Kinematic form of single-cell layer surface film model.
void append(const T &val)
const polyBoundaryMesh & boundaryMesh() const
const Type & value() const
const dimensionSet dimForce
addToRunTimeSelectionTable(surfaceFilmRegionModel, kinematicSingleLayer, mesh)
A HashTable with keys but without contents that is similar to std::unordered_set.
label min(const labelHashSet &set, label minValue=labelMax)
const dimensionedVector & g() const
Field< vector > vectorField
Specialisation of Field<T> for vector.
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
virtual const volScalarField & sigma() const =0
An indexed form of CGAL::Triangulation_face_base_2<K> used to keep track of the vertices in the trian...
tmp< vectorField > nf() const
Generic templated field type.
Base class for surface film models.
DynamicList< T, SizeMin > & append(const T &val)
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
void addToInjectedMass(const scalar dMass)
label max(const labelHashSet &set, label maxValue=labelMin)
virtual const surfaceScalarField & phi() const
ITstream & lookup(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Base class for film injection models, handling mass transfer from the film.
Macros for easy insertion into run-time selection tables.
virtual const volScalarField & rho() const =0
Generic dimensioned Type class.
Mesh data needed to do the Finite Volume discretisation.
const uniformDimensionedVectorField & g
errorManipArg< error, int > exit(error &err, const int errNo=1)
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
void correctBoundaryConditions()
const volScalarField & sigma() const
Curvature film separation model.
defineTypeNameAndDebug(kinematicSingleLayer, 0)
#define FatalErrorInFunction
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const fvPatch & patch() const
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Calculate the gradient of the given field.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
fileName::Type type(const fileName &name, const bool followLink=true)
virtual const volVectorField & nHat() const
const volScalarField & delta() const
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
const dimensionedScalar e
const surfaceFilmRegionModel & film() const
bool insert(const Key &key)
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
virtual const labelUList & faceCells() const
#define forAllReverse(list, i)
A List with indirect addressing.
Operations on lists of strings.
labelList findIndices(const ListType &input, typename ListType::const_reference val, label start=0)
dimensionedScalar sigma("sigma", dimMass/sqr(dimTime), transportProperties)
A 2-tuple for storing two objects of dissimilar types. The container is similar in purpose to std::pa...
virtual const volScalarField & rho() const
Generic GeometricField class.
Smooth ATC in cells next to a set of patches supplied by type.
const dimensionSet dimless
tmp< scalarField > calcCosAngle(const surfaceScalarField &phi) const
const fvMesh & regionMesh() const