calcVisc_incomp_PL Subroutine

private subroutine calcVisc_incomp_PL(nNwtn, viscKine, omega, state, neigh, auxField, densPos, velPos, nSize, nSolve, nScalars, nAuxScalars, layout, convFac)

Calculate kinematic viscosity from nonNewtonian power-law model for incompressible model $\mu = K shearRate^(n-1)$

Shear rate is computed from strain rate which is computed from nonEquilibrium PDF which in turn computed from pre-collision PDF

Arguments

TypeIntentOptionalAttributesName
class(mus_nNwtn_type), intent(in) :: nNwtn

contains non-Newtonian model parameters loaded from config file

real(kind=rk), intent(inout) :: viscKine(:)

non-Netonian model

real(kind=rk), intent(in) :: omega(:)

Kinematic viscosity omega from last timestep

real(kind=rk), intent(in) :: state(:)

state array

integer, intent(in) :: neigh(:)

neigh array to obtain precollision pdf

real(kind=rk), intent(in) :: auxField(:)

Auxiliary field variable array

integer, intent(in) :: densPos

position of density in auxField

integer, intent(in) :: velPos(3)

position of velocity components in auxField

integer, intent(in) :: nSize

number of elements in state array

integer, intent(in) :: nSolve

Number of element to solve in this level

integer, intent(in) :: nScalars

number of scalars in state array

integer, intent(in) :: nAuxScalars

number of scalars in auxField array

type(mus_scheme_layout_type), intent(in) :: layout

scheme layout

type(mus_convertFac_type), intent(in) :: convFac

conversion factor to convert lattice to physical units


Contents


Variables

TypeVisibilityAttributesNameInitial
integer, private :: iElem
integer, private :: iDir
integer, private :: QQ
integer, private :: elemOff
real(kind=rk), private :: rho
real(kind=rk), private :: vel(3)
real(kind=rk), private :: f_preCol(layout%fStencil%QQ)

precollision PDF

real(kind=rk), private :: fEq(layout%fStencil%QQ)
real(kind=rk), private :: nEq(layout%fStencil%QQ)
real(kind=rk), private :: nEqTens(6)
real(kind=rk), private :: nEqTensMag
real(kind=rk), private :: shearRate
real(kind=rk), private :: strainRate
real(kind=rk), private :: viscDynaPhy
real(kind=rk), private :: coeffSR