CG  Version 25
Public Types | Public Member Functions | Public Attributes | List of all members
PenaltyWallFunctionBC Class Reference

#include <PenaltyWallFunction.h>

Inheritance diagram for PenaltyWallFunctionBC:
Inheritance graph
[legend]
Collaboration diagram for PenaltyWallFunctionBC:
Collaboration graph
[legend]

Public Types

enum  WallFunctions {
  slipWall, logLaw, simpleLogLaw, wernerWengle,
  fixedUTau, numberOfWallFunctions
}
 

Public Member Functions

 PenaltyWallFunctionBC (const aString &nm)
 
virtual ~PenaltyWallFunctionBC ()
 
virtual bool inputFromGI (GenericGraphicsInterface &gi)
 
virtual bool applyBC (Parameters &parameters, const real &t, const real &dt, realMappedGridFunction &u, const int &grid, int side0, int axis0, realMappedGridFunction *gridVelocity=0)
 
virtual void getShearStresses (const int &nd, Parameters &parameters, const real &distance, const ArraySimpleFixed< real, 3, 1, 1, 1 > &uavg, const ArraySimpleFixed< real, 3, 1, 1, 1 > tng[], real tau[])
 
- Public Member Functions inherited from PenaltySlipWallBC
 PenaltySlipWallBC (const aString &nm)
 
virtual ~PenaltySlipWallBC ()
 
virtual bool setBCCoefficients (Parameters &parameters, const real &t, const real &dt, realMappedGridFunction &u, realMappedGridFunction &coeff, const int &grid, int side0, int axis0, realMappedGridFunction *gridVelocity=0)
 
virtual bool addPenaltyForcing (Parameters &parameters, const real &t, const real &dt, const realMappedGridFunction &u, realMappedGridFunction &dudt, const int &grid, int side0, int axis0, const realMappedGridFunction *gridVelocity=0)
 
virtual const bool isPenaltyBC () const
 
- Public Member Functions inherited from Parameters::BCModifier
 BCModifier (const aString &nm)
 
virtual ~BCModifier ()
 

Public Attributes

real fixedWallDistance
 
real linearLayerTopYPlus
 
WallFunctions wallFunctionType
 
bool includeArtificialDissipationInShearStress
 
bool useFullVelocity
 
- Public Attributes inherited from PenaltySlipWallBC
real normalFlux
 
bool zeroTangentialVelocity
 
DBase::DataBase db
 
- Public Attributes inherited from Parameters::BCModifier
aString name
 

Member Enumeration Documentation

Enumerator
slipWall 
logLaw 
simpleLogLaw 
wernerWengle 
fixedUTau 
numberOfWallFunctions 

Constructor & Destructor Documentation

PenaltyWallFunctionBC::PenaltyWallFunctionBC ( const aString &  nm)
PenaltyWallFunctionBC::~PenaltyWallFunctionBC ( )
virtual

Member Function Documentation

bool PenaltyWallFunctionBC::applyBC ( Parameters parameters,
const real &  t,
const real &  dt,
realMappedGridFunction &  u,
const int &  grid,
int  side0,
int  axis0,
realMappedGridFunction *  gridVelocity = 0 
)
virtual
In 2D, on a curvilinear grid, this looks like:

nx * ( rx*Dor(u.t) + sx*Dos(u.t) ) + ny * (ry*Dor(u.t) + sy*Dos(u.t) ) = tau/nu

  • or - Dor(u.t)*(nx*rx + ny*ry) + Dos(u.t)*(nx*sx + ny*sy) = tau/nu mprod[0] mprod[1] where u is really u-gv. If we are on an axis=0 boundary (i.e. r direction is normal) we can add

(2*side-1)* ( hr^2*(ad21+|ux|*ad22)*Do(u.t) + hr^4*(ad41+|ux|*ad42)*D+D-D+(u.t) )/nu

to the left hand side. NOTE: to be consistent with the D+D- used in the second derivative in the rest of the code we could choose to replace the Do operator above with D-. Similarly, the third derivative in the dissipation would really be D-D+D-, but this would include the second ghost line.

Reimplemented from PenaltySlipWallBC.

References a, A_3D, A_4D, A_5D, ad21, ad22, assert(), axis, InsParameters::BoussinesqModel, cc, d, Parameters::dbase, dx, f, fixedWallDistance, getShearStresses(), Parameters::gridIsMoving(), i1, i2, i3, includeArtificialDissipationInShearStress, isRectangular, mask, maskLocal, maskp, mg, PenaltySlipWallBC::normalFlux, np, nrm, off(), order, OV_APP_TO_PTR_3D, OV_APP_TO_PTR_4D, OV_GET_LOCAL_ARRAY_CONDITIONAL, OV_GET_LOCAL_ARRAY_FROM, OV_RGF_TO_PTR_5D, rx, s, side, solveSmallSystem(), tc, u, U, uc, uLocal, vc, InsParameters::viscoPlasticModel, wc, x, and PenaltySlipWallBC::zeroTangentialVelocity.

void PenaltyWallFunctionBC::getShearStresses ( const int &  nd,
Parameters parameters,
const real &  distance,
const ArraySimpleFixed< real, 3, 1, 1, 1 > &  uavg,
const ArraySimpleFixed< real, 3, 1, 1, 1 >  tng[],
real  tau[] 
)
virtual

Here is the actual implementation of the wall model where we compute the shear stresses on the surface from the "mean" velocity uavg. There are currently three basic models:

1) u+ = (1/kappa) ln(E y+) : which is a traditional log layer model (can also be written (1/kappa) ln (y+) + B 2) u+ = 8.3 (y+)^(1/7) : the Werner-Wengle model 3) u+ = (1/kappa) ln(y/y0) : a simple "rough wall" model sometimes used for atmospheric boundary layers

u+ = uavg/utau, y+=y*utau/nu, y0=roughness height, kappa=Karman's constant, utau = sqrt(tau/rho)

tau = mu * d(uavg .dot. tangent)/dn = mu * normal .dot. grad(uavg.dot.tangent)

There are some options to the above models:

  • If y+<linearLayerTopYPlus then u+=y+

NOTE: this routine actually returns the shear stress divided by the density which is nu*d(u.tangent)/dn

We need to solve a nonlinear equation to get tau: u+ = (1/kappa) ln(E y+) which gives va = (utau/kappa) ln (E distance utau/nu)

  • or - va = (X/kappa) ln( a X ) where X=utau and a = E*distance/nu.

If we write F(X) = -va + X ln (a X)/kappa we want X such that F(X) = 0. Note dF/dX = (ln(a X) + 1/a)/kappa . A simple Newton iteration gives: F(X^k + dX) = F( X^k ) + dX dF(X^k)/dX so that dX = -F(X^k)/(dF(X^k)/dX) as usual.

Guess X^0 = the werner-wengle value

NOTE: I doubt it is really this simple...

References a, PenaltySlipWallBC::db, Parameters::dbase, e, F, fixedUTau, itmax, kappa, linearLayerTopYPlus, logLaw, simpleLogLaw, slipWall, tol, useFullVelocity, wallFunctionType, wernerWengle, and yplus.

Referenced by applyBC().

bool PenaltyWallFunctionBC::inputFromGI ( GenericGraphicsInterface &  gi)
virtual

Member Data Documentation

real PenaltyWallFunctionBC::fixedWallDistance

Referenced by applyBC(), and inputFromGI().

bool PenaltyWallFunctionBC::includeArtificialDissipationInShearStress

Referenced by applyBC(), and inputFromGI().

real PenaltyWallFunctionBC::linearLayerTopYPlus

Referenced by getShearStresses(), and inputFromGI().

bool PenaltyWallFunctionBC::useFullVelocity

Referenced by getShearStresses(), and inputFromGI().

WallFunctions PenaltyWallFunctionBC::wallFunctionType

Referenced by getShearStresses(), and inputFromGI().


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