27 #include "surfaceInterpolate.H"
52 Foam::fv::interfaceTurbulenceDamping::interfaceFraction
59 tmp<volScalarField::Internal> tA
89 const scalar nSf(
mag(
n[own[facei]] & Sf[facei]));
90 A[own[facei]] += nSf*(ialphaf[facei] - ialpha[own[facei]]);
91 sumnSf[own[facei]] += nSf;
94 const scalar nSf(
mag(
n[nei[facei]] & Sf[facei]));
95 A[nei[facei]] += nSf*(ialphaf[facei] - ialpha[nei[facei]]);
96 sumnSf[nei[facei]] += nSf;
107 const scalar nSf(
mag(
n[own[facei]] & Sf[facei]));
108 A[own[facei]] += nSf*(palphaf[facei] - ialpha[own[facei]]);
109 sumnSf[own[facei]] += nSf;
116 if (sumnSf[i] > small)
118 a[i] = 2*
mag(a[i])/sumnSf[i];
130 template<
class RhoType>
131 void Foam::fv::interfaceTurbulenceDamping::addRhoSup
135 fvMatrix<scalar>& eqn
140 Info<<
type() <<
": applying source to " << eqn.psi().name() <<
endl;
144 phase_.fluid().movingPhases();
148 movingPhases[0]*
sqr(movingPhases[0].fluidThermo().nu()()())
151 for (
label phasei=1; phasei<movingPhases.size(); phasei++)
155 *
sqr(movingPhases[phasei].fluidThermo().nu()()());
158 if (field.name() ==
"epsilon")
160 eqn +=
rho*interfaceFraction(phase_)*C2_*aSqrnu*turbulence_.k()()
163 else if (field.name() ==
"omega")
165 eqn +=
rho*interfaceFraction(phase_)*beta_*aSqrnu
166 /(
sqr(betaStar_)*
pow4(delta_));
171 <<
"Support for field " << field.name() <<
" is not implemented"
181 const word& sourceName,
182 const word& modelType,
188 phaseName_(
dict.lookup(
"phase")),
207 fieldName_ = epsilonName;
211 fieldName_ = omegaName;
216 <<
"Cannot find either " << epsilonName <<
" or " << omegaName
236 addRhoSup(
one(), field, eqn);
247 addRhoSup(
rho(), field, eqn);
266 alpha*
sqr(phase_.fluidThermo().nu()()())
271 eqn +=
rho()*interfaceFraction(
alpha)
272 *C2_*aSqrnu*turbulence_.k()()/
pow4(delta_);
276 eqn +=
rho()*interfaceFraction(
alpha)
277 *beta_*aSqrnu/(
sqr(betaStar_)*
pow4(delta_));
282 <<
"Support for field " << field.
name() <<
" is not implemented"
#define forAll(list, i)
Loop across all elements in list.
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
static tmp< DimensionedField< Type, GeoMesh, PrimitiveField > > New(const word &name, const Mesh &mesh, const dimensionSet &, const PrimitiveField< Type > &)
Return a temporary field constructed from name, mesh,.
Generic GeometricField class.
DimensionedField< Type, GeoMesh, PrimitiveField > Internal
Type of the internal field from which this GeometricField is derived.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
const word & name() const
Return name.
static word groupName(Name name, const word &group)
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Mesh data needed to do the Finite Volume discretisation.
const labelUList & owner() const
Internal face owner.
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
const surfaceVectorField & Sf() const
Return cell face area vectors.
const labelUList & neighbour() const
Internal face neighbour.
Finite volume model abstract base class.
const fvMesh & mesh() const
Return const access to the mesh database.
Free-surface phase turbulence damping function.
virtual bool movePoints()
Update for mesh motion.
virtual wordList addSupFields() const
Return the list of fields for which the option adds source term.
virtual void addSup(const volScalarField &field, fvMatrix< scalar > &eqn) const
Add source to mixture epsilon or omega equation.
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
interfaceTurbulenceDamping(const word &sourceName, const word &modelType, const fvMesh &mesh, const dictionary &dict)
Construct from explicit source name and mesh.
Abstract base class for turbulence models (RAS, LES and laminar).
bool foundObject(const word &name) const
Is the named Type in registry.
A class representing the concept of 1 (scalar(1)) used to avoid unnecessary manipulations for objects...
const phaseSystem & fluid() const
Return the system to which this phase belongs.
const dimensionedScalar & deltaN() const
Stabilisation for normalisation of the interface normal.
UPtrList< phaseModel > phaseModelPartialList
Partial list of phase models.
Class containing mesh-to-mesh mapping information after a mesh distribution where we send parts of me...
Class containing mesh-to-mesh mapping information.
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
A class for managing temporary objects.
A class for handling words, derived from string.
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Calculate the gradient of the given field.
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
addToRunTimeSelectionTable(fvConstraint, bound, dictionary)
defineTypeNameAndDebug(bound, 0)
static tmp< SurfaceField< Type > > interpolate(const VolField< Type > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
tmp< VolField< typename outerProduct< vector, Type >::type > > grad(const SurfaceField< Type > &ssf)
static const coefficient A("A", dimPressure, 611.21)
errorManipArg< error, int > exit(error &err, const int errNo=1)
List< word > wordList
A List of words.
VolField< vector > volVectorField
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Ostream & endl(Ostream &os)
Add newline and flush stream.
const dimensionSet dimless
SurfaceField< scalar > surfaceScalarField
const dimensionSet dimLength
void pow4(LagrangianPatchField< scalar > &f, const LagrangianPatchField< scalar > &f1)
void mag(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
VolField< scalar > volScalarField
void sqr(LagrangianPatchField< typename outerProduct< Type, Type >::type > &f, const LagrangianPatchField< Type > &f1)
SurfaceField< vector > surfaceVectorField
UList< label > labelUList
fvsPatchField< scalar > fvsPatchScalarField
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.