Go to the documentation of this file.
39 Foam::omegaWallFunctionFvPatchScalarField::blendingType
41 Foam::omegaWallFunctionFvPatchScalarField::blendingTypeNames
43 { blendingType::STEPWISE ,
"stepwise" },
44 { blendingType::MAX ,
"max" },
45 { blendingType::BINOMIAL2 ,
"binomial2" },
46 { blendingType::BINOMIAL ,
"binomial" },
47 { blendingType::EXPONENTIAL,
"exponential" },
48 { blendingType::TANH,
"tanh" }
65 const volScalarField::Boundary& bf =
omega.boundaryField();
70 if (isA<omegaWallFunctionFvPatchScalarField>(bf[patchi]))
92 const fvMesh&
mesh = omega.mesh();
114 DynamicList<label> omegaPatches(bf.size());
117 if (isA<omegaWallFunctionFvPatchScalarField>(bf[patchi]))
119 omegaPatches.append(patchi);
121 const labelUList& faceCells = bf[patchi].patch().faceCells();
122 for (
const auto& celli : faceCells)
129 cornerWeights_.setSize(bf.size());
130 for (
const auto& patchi : omegaPatches)
133 cornerWeights_[patchi] = 1.0/wf.patchInternalField();
136 G_.setSize(internalField().size(), 0.0);
137 omega_.setSize(internalField().size(), 0.0);
155 refCast<const omegaWallFunctionFvPatchScalarField>(bf[patchi]);
169 forAll(cornerWeights_, patchi)
171 if (!cornerWeights_[patchi].empty())
182 forAll(cornerWeights_, patchi)
184 if (!cornerWeights_[patchi].empty())
186 omegaWallFunctionFvPatchScalarField& opf = omegaPatch(patchi);
197 const List<scalar>& cornerWeights,
198 const fvPatch&
patch,
203 const label patchi =
patch.index();
205 const nutWallFunctionFvPatchScalarField& nutw =
210 const tmp<scalarField> tnuw = turbModel.nu(patchi);
213 const tmp<volScalarField> tk = turbModel.k();
220 const scalar Cmu25 =
pow025(nutw.Cmu());
225 const label celli =
patch.faceCells()[facei];
226 const scalar
yPlus = Cmu25*
y[facei]*
sqrt(
k[celli])/nuw[facei];
227 const scalar w = cornerWeights[facei];
230 const scalar omegaVis = 6.0*nuw[facei]/(beta1_*
sqr(
y[facei]));
233 const scalar omegaLog =
sqrt(
k[celli])/(Cmu25*nutw.kappa()*
y[facei]);
237 case blendingType::STEPWISE:
239 if (
yPlus > nutw.yPlusLam())
241 omega0[celli] += w*omegaLog;
245 omega0[celli] += w*omegaVis;
250 case blendingType::MAX:
253 omega0[celli] +=
max(omegaVis, omegaLog);
257 case blendingType::BINOMIAL2:
260 omega0[celli] += w*
sqrt(
sqr(omegaVis) +
sqr(omegaLog));
264 case blendingType::BINOMIAL:
269 pow(omegaVis, n_) +
pow(omegaLog, n_),
275 case blendingType::EXPONENTIAL:
279 const scalar invGamma = 1.0/(Gamma + ROOTVSMALL);
282 w*(omegaVis*
exp(-Gamma) + omegaLog*
exp(-invGamma));
286 case blendingType::TANH:
290 const scalar omegab1 = omegaVis + omegaLog;
291 const scalar omegab2 =
292 pow(
pow(omegaVis, 1.2) +
pow(omegaLog, 1.2), 1.0/1.2);
294 omega0[celli] += phiTanh*omegab1 + (1.0 - phiTanh)*omegab2;
299 if (!(blending_ == blendingType::STEPWISE) ||
yPlus > nutw.yPlusLam())
303 *(nutw[facei] + nuw[facei])
305 *Cmu25*
sqrt(
k[celli])
306 /(nutw.kappa()*
y[facei]);
317 const DimensionedField<scalar, volMesh>& iF
320 fixedValueFvPatchField<scalar>(
p, iF),
321 blending_(blendingType::BINOMIAL2),
341 blending_(ptf.blending_),
362 blendingTypeNames.getOrDefault
366 blendingType::BINOMIAL2
371 dict.getCheckOrDefault<scalar>
380 beta1_(
dict.getOrDefault<scalar>(
"beta1", 0.075)),
386 if (
dict.found(
"blended"))
389 <<
"Using deprecated 'blended' keyword"
390 <<
nl <<
" Please use either of the below for the same behaviour:"
391 <<
nl <<
" 'blending binomial2;' for 'blended on;'"
392 <<
nl <<
" 'blending stepwise;' for 'blended off;'"
393 <<
nl <<
" OVERWRITING: 'blended' keyword -> 'blending' enum"
400 blending_ = blendingType::BINOMIAL2;
404 blending_ = blendingType::STEPWISE;
419 blending_(owfpsf.blending_),
423 beta1_(owfpsf.beta1_),
437 blending_(owfpsf.blending_),
441 beta1_(owfpsf.beta1_),
455 if (
patch().index() == master_)
465 return omegaPatch(master_).G();
474 if (
patch().index() == master_)
484 return omegaPatch(master_).omega(
init);
500 internalField().
group()
506 if (
patch().index() == master_)
508 createAveragingWeights();
509 calculateTurbulenceFields(turbModel,
G(
true), omega(
true));
515 typedef DimensionedField<scalar, volMesh> FieldType;
517 FieldType&
G = db().lookupObjectRef<FieldType>(turbModel.GName());
519 FieldType& omega =
const_cast<FieldType&
>(internalField());
523 const label celli =
patch().faceCells()[facei];
525 G[celli] =
G0[celli];
526 omega[celli] = omega0[celli];
548 internalField().
group()
554 if (
patch().index() == master_)
556 createAveragingWeights();
557 calculateTurbulenceFields(turbModel,
G(
true), omega(
true));
563 typedef DimensionedField<scalar, volMesh> FieldType;
565 FieldType&
G = db().lookupObjectRef<FieldType>(turbModel.GName());
567 FieldType& omega =
const_cast<FieldType&
>(internalField());
574 const scalar w = weights[facei];
578 const label celli =
patch().faceCells()[facei];
580 G[celli] = (1.0 - w)*
G[celli] + w*
G0[celli];
581 omega[celli] = (1.0 - w)*omega[celli] + w*omega0[celli];
582 omegaf[facei] = omega[celli];
592 fvMatrix<scalar>& matrix
595 if (manipulatedMatrix())
612 if (manipulatedMatrix())
617 DynamicList<label> constraintCells(weights.size());
618 DynamicList<scalar> constraintValues(weights.size());
621 const DimensionedField<scalar, volMesh>&
fld = internalField();
626 if (tolerance_ < weights[facei])
628 const label celli = faceCells[facei];
630 constraintCells.append(celli);
631 constraintValues.append(
fld[celli]);
638 <<
": number of constrained cells = " << constraintCells.size()
639 <<
" out of " <<
patch().size()
654 os.
writeEntry(
"blending", blendingTypeNames[blending_]);
668 omegaWallFunctionFvPatchScalarField
virtual void createAveragingWeights()
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
fvPatchField< scalar > fvPatchScalarField
omegaWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
const DimensionedField< Type, volMesh > & internalField() const
virtual void operator==(const fvPatchField< Type > &)
virtual tmp< Field< Type > > snGrad() const
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
virtual void manipulateMatrix(fvMatrix< scalar > &matrix)
virtual void calculate(const turbulenceModel &turbulence, const List< scalar > &cornerWeights, const fvPatch &patch, scalarField &G, scalarField &omega)
Defines the attributes of an object for which implicit objectRegistry management is supported,...
virtual tmp< volScalarField > nu() const =0
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
This boundary condition provides a wall constraint on the specific dissipation rate,...
virtual void updateWeightedCoeffs(const scalarField &weights)
constexpr const char *const group
const dimensionedScalar G
A class for managing temporary objects.
static constexpr const zero Zero
bool found(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
static word timeName(const scalar t, const int precision=precision_)
static const word propertiesName
Ostream & endl(Ostream &os)
dimensionedScalar exp(const dimensionedScalar &ds)
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
const nearWallDist & y() const
virtual void write(Ostream &) const
dimensionedScalar pow025(const dimensionedScalar &ds)
scalarField & omega(bool init=false)
dimensionedScalar pow4(const dimensionedScalar &ds)
linear/upwind blended differencing scheme.
Generic templated field type.
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
ThermalDiffusivity< CompressibleTurbulenceModel< fluidThermo > > turbulenceModel
dimensionedScalar tanh(const dimensionedScalar &ds)
const dimensionedScalar G0
DynamicList< T, SizeMin > & append(const T &val)
virtual void write(Ostream &) const
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
virtual tmp< Field< Type > > patchInternalField() const
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
virtual omegaWallFunctionFvPatchScalarField & omegaPatch(const label patchi)
T getCheckOrDefault(const word &keyword, const T &deflt, const Predicate &pred, enum keyType::option matchOpt=keyType::REGEX) const
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
label max(const labelHashSet &set, label maxValue=labelMin)
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< ' ';}gmvFile<< nl;for(const word &name :lagrangianScalarNames){ IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Abstract base class for turbulence models (RAS, LES and laminar).
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
OBJstream os(runTime.globalPath()/outputName)
Macros for easy insertion into run-time selection tables.
Generic dimensioned Type class.
virtual void calculateTurbulenceFields(const turbulenceModel &turbulence, scalarField &G0, scalarField &omega0)
virtual tmp< volScalarField > k() const =0
Mesh data needed to do the Finite Volume discretisation.
static scalar yPlusLam(const scalar kappa, const scalar E)
virtual void manipulateMatrix(fvMatrix< Type > &matrix)
virtual void updateCoeffs()
void setValues(const labelUList &cellLabels, const Type &value)
fvPatchField< vector > fvPatchVectorField
bool changing() const noexcept
dimensionedSymmTensor sqr(const dimensionedVector &dv)
static MinMax< T > ge(const T &minVal)
static const nutWallFunctionFvPatchScalarField & nutw(const turbulenceModel &turbModel, const label patchi)
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
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
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Ostream & writeEntry(const keyType &key, const T &value)
virtual void updateCoeffs()
virtual const labelUList & faceCells() const
const Time & time() const
Foam::fvPatchFieldMapper.
const fvPatch & patch() const
static word groupName(StringType base, const word &group)
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
#define IOWarningInFunction(ios)
UList< label > labelUList
A UList of labels.
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Generic GeometricField class.
The class nutWallFunction is a base class that parents the derived boundary conditions which provide ...
Smooth ATC in cells next to a set of patches supplied by type.
A min/max value pair with additional methods. In addition to conveniently storing values,...
const Boundary & boundaryField() const
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const dimensionSet dimless
const volVectorField & U() const
scalarField & G(bool init=false)