35#include "surfaceInterpolate.H"
82 scalarField& invR1 = tinvR1.ref().primitiveFieldRef();
85 const scalar rMin = 1
e-6;
96 const scalar rMax = 1e6;
99 if (
mag(invR1[i]) < 1/rMax)
128 label cellO = own[facei];
129 label cellN = nbr[facei];
131 if (
phi[facei] > phiMax[cellO])
133 phiMax[cellO] =
phi[facei];
134 cosAngle[cellO] = -
gHat_ & nf[facei];
136 if (-
phi[facei] > phiMax[cellN])
138 phiMax[cellN] = -
phi[facei];
139 cosAngle[cellN] = -
gHat_ & -nf[facei];
146 const fvPatch&
pp = phip.
patch();
151 label celli = faceCells[i];
152 if (phip[i] > phiMax[celli])
154 phiMax[celli] = phip[i];
155 cosAngle[celli] = -
gHat_ & nf[i];
196 if (debug &&
mesh.time().writeTime())
206 auto& volCosAngle = tvolCosAngle.ref();
208 volCosAngle.primitiveFieldRef() = cosAngle;
209 volCosAngle.correctBoundaryConditions();
219curvatureSeparation::curvatureSeparation
226 gradNHat_(
fvc::grad(film.nHat())),
227 deltaByR1Min_(coeffDict_.getOrDefault<scalar>(
"deltaByR1Min", 0)),
228 definedPatchRadii_(),
229 magG_(
mag(film.
g().value())),
232 if (
magG_ < ROOTVSMALL)
235 <<
"Acceleration due to gravity must be non-zero"
255 if (!uniquePatchIDs.
found(patchi))
257 const scalar radius = prIn[i].second();
260 uniquePatchIDs.
insert(patchi);
299 const scalar Fthreshold = 1
e-10;
306 scalar R1 = 1.0/(invR1[i] + ROOTVSMALL);
307 scalar R2 = R1 +
delta[i];
310 scalar Fi = -
delta[i]*
rho[i]*magSqrU[i]*72.0/60.0*invR1[i];
317 scalar Fs =
sigma[i]/R2;
319 Fnet[i] = Fi +
Fb + Fs;
321 if (Fnet[i] + Fthreshold < 0)
329 massToInject = separated*availableMass;
330 diameterToInject = separated*
delta;
331 availableMass -= separated*availableMass;
335 if (debug &&
mesh.time().writeTime())
345 auto& volFnet = tvolFnet.ref();
347 volFnet.primitiveFieldRef() = Fnet;
348 volFnet.correctBoundaryConditions();
CGAL::indexedFace< K > Fb
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
const polyBoundaryMesh & pbm
const uniformDimensionedVectorField & g
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
void append(const T &val)
Copy append an element to the end of this list.
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
bool insert(const Key &key)
Insert a new entry, not overwriting existing entries.
bool found(const Key &key) const
Same as contains().
@ NO_REGISTER
Do not request registration (bool: false).
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
A 2-tuple for storing two objects of dissimilar types. The container is similar in purpose to std::pa...
A List with indirect addressing. Like IndirectList but does not store addressing.
void size(const label n)
Older name for setAddressableSize.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
ITstream & lookup(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return an entry data stream. FatalIOError if not found, or not a stream.
Smooth ATC in cells next to a set of patches supplied by type.
Mesh data needed to do the Finite Volume discretisation.
static const word & zeroGradientType() noexcept
The type name for zeroGradient patch fields.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
const fvPatch & patch() const noexcept
Return the patch.
A polyBoundaryMesh is a polyPatch list with registered IO, a reference to the associated polyMesh,...
wordList names() const
Return a list of patch names.
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
const fvMesh & regionMesh() const
Return the region mesh database.
Curvature film separation model.
scalar deltaByR1Min_
Minimum gravity driven film thickness (non-dimensionalised delta/R1).
volTensorField gradNHat_
Gradient of surface normals.
tmp< volScalarField > calcInvR1(const volVectorField &U) const
Calculate local (inverse) radius of curvature.
scalar magG_
Magnitude of gravity vector.
tmp< scalarField > calcCosAngle(const surfaceScalarField &phi) const
Calculate the cosine of the angle between gravity vector and.
vector gHat_
Direction of gravity vector.
List< Tuple2< label, scalar > > definedPatchRadii_
List of radii for patches - if patch not defined, radius.
virtual ~curvatureSeparation()
Destructor.
const surfaceFilmRegionModel & film() const
Return const access to the film surface film model.
Base class for film injection models, handling mass transfer from the film.
void addToInjectedMass(const scalar dMass)
Add to injected mass.
Kinematic form of single-cell layer surface film model.
Base class for surface film models.
virtual const volScalarField & sigma() const =0
Return the film surface tension [N/m].
virtual const volVectorField & U() const =0
Return the film velocity [m/s].
virtual const volScalarField & rho() const =0
Return the film density [kg/m3].
const dimensionedVector & g() const
Return the acceleration due to gravity.
virtual const volScalarField & delta() const =0
Return the film thickness [m].
const dictionary coeffDict_
Coefficients dictionary.
const dictionary & dict() const
Return const access to the cloud dictionary.
A class for managing temporary objects.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Calculate the divergence of the given field.
Calculate the gradient of the given field.
Namespace for handling debugging switches.
Namespace of functions to calculate explicit derivatives.
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
List< word > wordList
List of word.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
GeometricField< vector, fvPatchField, volMesh > volVectorField
const dimensionSet dimless
Dimensionless.
List< label > labelList
A List of labels.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
const dimensionSet dimVelocity
dimensionSet clamp(const dimensionSet &a, const dimensionSet &range)
MinMax< scalar > scalarMinMax
A scalar min/max range.
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
const dimensionSet dimForce
Field< vector > vectorField
Specialisation of Field<T> for vector.
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &f1, const label comm)
static constexpr const zero Zero
Global zero (0).
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
labelList findIndices(const ListType &input, typename ListType::const_reference val, label start=0)
Linear search to find all occurrences of given element.
errorManipArg< error, int > exit(error &err, const int errNo=1)
UList< label > labelUList
A UList of labels.
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
fvsPatchField< scalar > fvsPatchScalarField
dimensionedScalar sigma("sigma", dimMass/sqr(dimTime), transportProperties)
#define forAll(list, i)
Loop across all elements in list.
#define forAllReverse(list, i)
Reverse loop across all elements in list.