Memoery allocation to the gradient variables and setup pointer to the slice to the main gradient variable based on the various models being used.
Setup Pointer to the main array which stores gradient all variables with x, y, z
Setup Pointer to the main array which stores gradient all variables with x, y, z
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
type(controltype), | intent(in) | :: | control | Control parameters |
||
type(schemetype), | intent(in) | :: | scheme | finite-volume Schemes |
||
type(flowtype), | intent(in) | :: | flow | Information about fluid flow: freestream-speed, ref-viscosity,etc. |
||
type(extent), | intent(in) | :: | dims | Extent of the domain:imx,jmx,kmx |
subroutine setup_gradients(control, scheme, flow, dims)
!< Memoery allocation to the gradient variables and
!< setup pointer to the slice to the main gradient variable
!< based on the various models being used.
implicit none
type(controltype), intent(in) :: control
!< Control parameters
type(schemetype) , intent(in) :: scheme
!< finite-volume Schemes
type(flowtype) , intent(in) :: flow
!< Information about fluid flow: freestream-speed, ref-viscosity,etc.
type(extent), intent(in) :: dims
!< Extent of the domain:imx,jmx,kmx
DebugCall("setup_gradients")
imx = dims%imx
jmx = dims%jmx
kmx = dims%kmx
n_var = control%n_var
R_gas = flow%R_gas
if(flow%mu_ref/=0)then
call get_n_grad(scheme)
!call allocate_memory()
call alloc(gradqp_x, 0, imx, 0, jmx, 0, kmx, 1, n_grad, AErrMsg("gradqp_x"))
call alloc(gradqp_y, 0, imx, 0, jmx, 0, kmx, 1, n_grad, AErrMsg("gradqp_y"))
call alloc(gradqp_z, 0, imx, 0, jmx, 0, kmx, 1, n_grad, AErrMsg("gradqp_z"))
! Linking pointer to laminar gradients
!call setup_laminar_grad()
DebugCall('setup_laminar_grad')
!< Setup Pointer to the main array which stores gradient
!< all variables with x, y, z
gradu_x(0:imx, 0:jmx, 0:kmx) => gradqp_x(:, :, :, 1)
gradv_x(0:imx, 0:jmx, 0:kmx) => gradqp_x(:, :, :, 2)
gradw_x(0:imx, 0:jmx, 0:kmx) => gradqp_x(:, :, :, 3)
gradT_x(0:imx, 0:jmx, 0:kmx) => gradqp_x(:, :, :, 4)
gradu_y(0:imx, 0:jmx, 0:kmx) => gradqp_y(:, :, :, 1)
gradv_y(0:imx, 0:jmx, 0:kmx) => gradqp_y(:, :, :, 2)
gradw_y(0:imx, 0:jmx, 0:kmx) => gradqp_y(:, :, :, 3)
gradT_y(0:imx, 0:jmx, 0:kmx) => gradqp_y(:, :, :, 4)
gradu_z(0:imx, 0:jmx, 0:kmx) => gradqp_z(:, :, :, 1)
gradv_z(0:imx, 0:jmx, 0:kmx) => gradqp_z(:, :, :, 2)
gradw_z(0:imx, 0:jmx, 0:kmx) => gradqp_z(:, :, :, 3)
gradT_z(0:imx, 0:jmx, 0:kmx) => gradqp_z(:, :, :, 4)
! Linking pointer to turbulent gradients
select case (trim(scheme%turbulence))
case('none')
!do nothing
continue
case('sa', 'saBC')
!call setup_sa_grad()
DebugCall("setup_sa_grad")
gradtv_x(0:imx, 0:jmx, 0:kmx) => gradqp_x(:, :, :, 5)
gradtv_y(0:imx, 0:jmx, 0:kmx) => gradqp_y(:, :, :, 5)
gradtv_z(0:imx, 0:jmx, 0:kmx) => gradqp_z(:, :, :, 5)
case('sst', 'sst2003')
!< Setup Pointer to the main array which stores gradient
!< all variables with x, y, z
DebugCall('setup_sst_grad')
gradtk_x(0:imx, 0:jmx, 0:kmx) => gradqp_x(:, :, :, 5)
gradtw_x(0:imx, 0:jmx, 0:kmx) => gradqp_x(:, :, :, 6)
gradtk_y(0:imx, 0:jmx, 0:kmx) => gradqp_y(:, :, :, 5)
gradtw_y(0:imx, 0:jmx, 0:kmx) => gradqp_y(:, :, :, 6)
gradtk_z(0:imx, 0:jmx, 0:kmx) => gradqp_z(:, :, :, 5)
gradtw_z(0:imx, 0:jmx, 0:kmx) => gradqp_z(:, :, :, 6)
case('kkl')
!call setup_kkl_grad()
DebugCall('setup_kkl_grad')
gradtk_x(0:imx, 0:jmx, 0:kmx) => gradqp_x(:, :, :, 5)
gradtkl_x(0:imx, 0:jmx, 0:kmx) => gradqp_x(:, :, :, 6)
gradtk_y(0:imx, 0:jmx, 0:kmx) => gradqp_y(:, :, :, 5)
gradtkl_y(0:imx, 0:jmx, 0:kmx) => gradqp_y(:, :, :, 6)
gradtk_z(0:imx, 0:jmx, 0:kmx) => gradqp_z(:, :, :, 5)
gradtkl_z(0:imx, 0:jmx, 0:kmx) => gradqp_z(:, :, :, 6)
case DEFAULT
Fatal_error
end select
!Transition modeling
select case(trim(scheme%transition))
case('lctm2015')
!call setup_lctm2015_grad()
gradtgm_x(0:imx, 0:jmx, 0:kmx) => gradqp_x(:, :, :, n_grad)
gradtgm_y(0:imx, 0:jmx, 0:kmx) => gradqp_y(:, :, :, n_grad)
gradtgm_z(0:imx, 0:jmx, 0:kmx) => gradqp_z(:, :, :, n_grad)
case('none','bc')
!do nothing
continue
case DEFAULT
Fatal_error
end Select
end if
end subroutine setup_gradients