Loading...
Searching...
No Matches
atmAlphatkWallFunctionFvPatchScalarField Class Reference

This boundary condition provides a wall constraint on the kinematic turbulent thermal conductivity (i.e. alphat) for atmospheric boundary layer modelling. It assumes a logarithmic distribution of the potential temperature within the first cell. More...

#include <atmAlphatkWallFunctionFvPatchScalarField.H>

Inheritance diagram for atmAlphatkWallFunctionFvPatchScalarField:
Collaboration diagram for atmAlphatkWallFunctionFvPatchScalarField:

Public Member Functions

 TypeName ("atmAlphatkWallFunction")
 Runtime type information.
 atmAlphatkWallFunctionFvPatchScalarField (const fvPatch &, const DimensionedField< scalar, volMesh > &)
 Construct from patch and internal field.
 atmAlphatkWallFunctionFvPatchScalarField (const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &)
 Construct from patch, internal field and dictionary.
 atmAlphatkWallFunctionFvPatchScalarField (const atmAlphatkWallFunctionFvPatchScalarField &, const fvPatch &, const DimensionedField< scalar, volMesh > &, const fvPatchFieldMapper &)
 Construct by mapping given atmAlphatkWallFunctionFvPatchScalarField onto a new patch.
 atmAlphatkWallFunctionFvPatchScalarField (const atmAlphatkWallFunctionFvPatchScalarField &)
 Construct as copy.
 atmAlphatkWallFunctionFvPatchScalarField (const atmAlphatkWallFunctionFvPatchScalarField &, const DimensionedField< scalar, volMesh > &)
 Construct as copy setting internal field reference.
virtual tmp< fvPatchField< scalar > > clone () const
 Return a clone.
virtual tmp< fvPatchField< scalar > > clone (const DimensionedField< scalar, volMesh > &iF) const
 Clone with an internal field reference.
virtual void updateCoeffs ()
 Update the coefficients associated with the patch field.
virtual void autoMap (const fvPatchFieldMapper &)
 Map (and resize as needed) from self given a mapping object.
virtual void rmap (const fvPatchScalarField &, const labelList &)
 Reverse map the given fvPatchField onto this fvPatchField.
virtual void write (Ostream &) const
 Write.

Protected Member Functions

virtual void checkType ()
 Check the type of the patch.
void writeLocalEntries (Ostream &) const
 Write local wall function variables.

Protected Attributes

const scalar Cmu_
 Empirical model constant.
const scalar kappa_
 von Kármán constant
autoPtr< Function1< scalar > > Pr_
 Molecular Prandtl number.
autoPtr< PatchFunction1< scalar > > Prt_
 Turbulent Prandtl number field.
autoPtr< PatchFunction1< scalar > > z0_
 Surface roughness length [m].

Static Protected Attributes

static scalar tolerance_ = 0.01
 Solution parameters.
static label maxIters_ = 10

Detailed Description

This boundary condition provides a wall constraint on the kinematic turbulent thermal conductivity (i.e. alphat) for atmospheric boundary layer modelling. It assumes a logarithmic distribution of the potential temperature within the first cell.

Required fields:

  alphat    | Kinematic turbulent thermal conductivity    [m2/s]
Usage
Example of the boundary condition specification:
<patchName>
{
    // Mandatory entries
    type            atmAlphatkWallFunction;
    Pr              <Function1<scalar>>;
    Prt             <PatchFunction1<scalar>>;
    z0              <PatchFunction1<scalar>>;

    // Optional entries
    Cmu             <scalar>;
    kappa           <scalar>;

    // Inherited entries
    ...
}

where the entries mean:

Property Description Type Reqd Deflt
type Type name: atmAlphatkWallFunction word yes -
Pr Molecular Prandtl number Function1<scalar> yes -
Prt Turbulent Prandtl number PatchFunction1<scalar> yes -
z0 Surface roughness length [m] PatchFunction1<scalar> yes -
Cmu Empirical model constant scalar no 0.09
kappa von Kármán constant scalar no 0.41

The inherited entries are elaborated in:

Source files

Definition at line 140 of file atmAlphatkWallFunctionFvPatchScalarField.H.

Constructor & Destructor Documentation

◆ atmAlphatkWallFunctionFvPatchScalarField() [1/5]

atmAlphatkWallFunctionFvPatchScalarField ( const fvPatch & p,
const DimensionedField< scalar, volMesh > & iF )

Construct from patch and internal field.

Definition at line 80 of file atmAlphatkWallFunctionFvPatchScalarField.C.

References checkType(), Cmu_, kappa_, p, Pr_, Prt_, and z0_.

Referenced by atmAlphatkWallFunctionFvPatchScalarField(), atmAlphatkWallFunctionFvPatchScalarField(), atmAlphatkWallFunctionFvPatchScalarField(), and TypeName().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ atmAlphatkWallFunctionFvPatchScalarField() [2/5]

atmAlphatkWallFunctionFvPatchScalarField ( const fvPatch & p,
const DimensionedField< scalar, volMesh > & iF,
const dictionary & dict )

Construct from patch, internal field and dictionary.

Definition at line 118 of file atmAlphatkWallFunctionFvPatchScalarField.C.

References checkType(), Cmu_, dict, kappa_, Foam::New(), p, Pr_, Prt_, and z0_.

Here is the call graph for this function:

◆ atmAlphatkWallFunctionFvPatchScalarField() [3/5]

atmAlphatkWallFunctionFvPatchScalarField ( const atmAlphatkWallFunctionFvPatchScalarField & ptf,
const fvPatch & p,
const DimensionedField< scalar, volMesh > & iF,
const fvPatchFieldMapper & mapper )

Construct by mapping given atmAlphatkWallFunctionFvPatchScalarField onto a new patch.

Definition at line 98 of file atmAlphatkWallFunctionFvPatchScalarField.C.

References atmAlphatkWallFunctionFvPatchScalarField(), checkType(), clone(), Cmu_, kappa_, p, Pr_, Prt_, and z0_.

Here is the call graph for this function:

◆ atmAlphatkWallFunctionFvPatchScalarField() [4/5]

atmAlphatkWallFunctionFvPatchScalarField ( const atmAlphatkWallFunctionFvPatchScalarField & wfpsf)

Construct as copy.

Definition at line 153 of file atmAlphatkWallFunctionFvPatchScalarField.C.

References atmAlphatkWallFunctionFvPatchScalarField(), checkType(), clone(), Cmu_, kappa_, Pr_, Prt_, and z0_.

Here is the call graph for this function:

◆ atmAlphatkWallFunctionFvPatchScalarField() [5/5]

atmAlphatkWallFunctionFvPatchScalarField ( const atmAlphatkWallFunctionFvPatchScalarField & wfpsf,
const DimensionedField< scalar, volMesh > & iF )

Construct as copy setting internal field reference.

Definition at line 170 of file atmAlphatkWallFunctionFvPatchScalarField.C.

References atmAlphatkWallFunctionFvPatchScalarField(), checkType(), clone(), Cmu_, kappa_, Pr_, Prt_, and z0_.

Here is the call graph for this function:

Member Function Documentation

◆ checkType()

void checkType ( )
protectedvirtual

◆ writeLocalEntries()

void writeLocalEntries ( Ostream & os) const
protected

Write local wall function variables.

Definition at line 55 of file atmAlphatkWallFunctionFvPatchScalarField.C.

References Cmu_, kappa_, os(), Pr_, Prt_, and z0_.

Referenced by write().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ TypeName()

TypeName ( "atmAlphatkWallFunction" )

Runtime type information.

References atmAlphatkWallFunctionFvPatchScalarField().

Here is the call graph for this function:

◆ clone() [1/2]

virtual tmp< fvPatchField< scalar > > clone ( ) const
inlinevirtual

Return a clone.

Definition at line 257 of file atmAlphatkWallFunctionFvPatchScalarField.H.

References fvPatchField< Type >::Clone().

Referenced by atmAlphatkWallFunctionFvPatchScalarField(), atmAlphatkWallFunctionFvPatchScalarField(), and atmAlphatkWallFunctionFvPatchScalarField().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ clone() [2/2]

virtual tmp< fvPatchField< scalar > > clone ( const DimensionedField< scalar, volMesh > & iF) const
inlinevirtual

Clone with an internal field reference.

Definition at line 265 of file atmAlphatkWallFunctionFvPatchScalarField.H.

References fvPatchField< Type >::Clone().

Here is the call graph for this function:

◆ updateCoeffs()

void updateCoeffs ( )
virtual

Update the coefficients associated with the patch field.

Definition at line 190 of file atmAlphatkWallFunctionFvPatchScalarField.C.

References Cmu_, e, Foam::exit(), Foam::FatalError, FatalErrorInFunction, forAll, IOobject::groupName(), k, kappa_, Foam::log(), Foam::max(), Foam::pow025(), Pr(), Pr_, turbulenceModel::propertiesName, Prt(), Prt_, Foam::sqrt(), fvPatchField< Type >::updateCoeffs(), y, and z0_.

Here is the call graph for this function:

◆ autoMap()

void autoMap ( const fvPatchFieldMapper & m)
virtual

Map (and resize as needed) from self given a mapping object.

Definition at line 273 of file atmAlphatkWallFunctionFvPatchScalarField.C.

References Prt_, and z0_.

◆ rmap()

void rmap ( const fvPatchScalarField & ptf,
const labelList & addr )
virtual

Reverse map the given fvPatchField onto this fvPatchField.

Definition at line 291 of file atmAlphatkWallFunctionFvPatchScalarField.C.

References Prt_, Foam::refCast(), and z0_.

Here is the call graph for this function:

◆ write()

void write ( Ostream & os) const
virtual

Write.

Definition at line 313 of file atmAlphatkWallFunctionFvPatchScalarField.C.

References os(), fvPatchField< Type >::write(), writeLocalEntries(), and fvPatchField< Type >::writeValueEntry().

Here is the call graph for this function:

Member Data Documentation

◆ Cmu_

◆ kappa_

◆ Pr_

◆ Prt_

◆ z0_

◆ tolerance_

scalar tolerance_ = 0.01
staticprotected

Solution parameters.

Definition at line 178 of file atmAlphatkWallFunctionFvPatchScalarField.H.

◆ maxIters_

label maxIters_ = 10
staticprotected

Definition at line 179 of file atmAlphatkWallFunctionFvPatchScalarField.H.


The documentation for this class was generated from the following files: