Write the header and variables in the file "process_xx.dat".
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
integer, | intent(in) | :: | file_handler | |||
real(kind=wp), | intent(in), | dimension(-2:dims%imx+2, -2:dims%jmx+2, -2:dims%kmx+2, 1:dims%n_var), target | :: | state | ||
type(nodetype), | intent(in), | dimension(-2:dims%imx+3,-2:dims%jmx+3,-2:dims%kmx+3) | :: | nodes | ||
type(controltype), | intent(in) | :: | control | |||
type(schemetype), | intent(in) | :: | scheme | |||
type(extent), | intent(in) | :: | dims |
subroutine write_file(file_handler, state, nodes, control, scheme, dims)
!< Write the header and variables in the file "process_xx.dat".
implicit none
integer, intent(in) :: file_handler
type(controltype), intent(in) :: control
type(schemetype), intent(in) :: scheme
type(extent), intent(in) :: dims
type(nodetype), dimension(-2:dims%imx+3,-2:dims%jmx+3,-2:dims%kmx+3), intent(in) :: nodes
real(wp), dimension(-2:dims%imx+2, -2:dims%jmx+2, -2:dims%kmx+2, 1:dims%n_var), intent(in), target :: state
integer :: n
character(len=*), parameter :: err="Write error: Asked to write non-existing variable- "
DebugCall("write_file")
OUT_FILE_UNIT = file_handler
imx = dims%imx
jmx = dims%jmx
kmx = dims%kmx
density(-2:imx+2, -2:jmx+2, -2:kmx+2) => state(:, :, :, 1)
x_speed(-2:imx+2, -2:jmx+2, -2:kmx+2) => state(:, :, :, 2)
y_speed(-2:imx+2, -2:jmx+2, -2:kmx+2) => state(:, :, :, 3)
z_speed(-2:imx+2, -2:jmx+2, -2:kmx+2) => state(:, :, :, 4)
pressure(-2:imx+2, -2:jmx+2, -2:kmx+2) => state(:, :, :, 5)
select case (trim(scheme%turbulence))
case ("none")
!include nothing
continue
case ("sst", "sst2003", "bsl", "des-sst", "kw")
tk(-2:imx+2, -2:jmx+2, -2:kmx+2) => state(:, :, :, 6)
tw(-2:imx+2, -2:jmx+2, -2:kmx+2) => state(:, :, :, 7)
case ("kkl")
tk(-2:imx+2, -2:jmx+2, -2:kmx+2) => state(:, :, :, 6)
tkl(-2:imx+2, -2:jmx+2, -2:kmx+2) => state(:, :, :, 7)
case ("sa", "saBC")
tv(-2:imx+2, -2:jmx+2, -2:kmx+2) => state(:, :, :, 6)
case ("ke")
tk(-2:imx+2, -2:jmx+2, -2:kmx+2) => state(:, :, :, 6)
te(-2:imx+2, -2:jmx+2, -2:kmx+2) => state(:, :, :, 7)
case DEFAULT
Fatal_error
end select
! Transition modeling
select case(trim(scheme%transition))
case('lctm2015')
tgm(-2:imx+2, -2:jmx+2, -2:kmx+2) => state(:, :, :, 8)
case('bc', 'none')
!do nothing
continue
case DEFAULT
Fatal_error
end Select
call write_header(control)
call write_grid(nodes)
do n = 1,control%w_count
select case (trim(control%w_list(n)))
case('Velocity')
call write_scalar(x_speed, "u", -2)
call write_scalar(y_speed, "v", -2)
call write_scalar(z_speed, "w", -2)
case('Density')
call write_scalar(density, "Density", -2)
case('Pressure')
call write_scalar(pressure, "Pressure", -2)
case('Mu')
call write_scalar(mu, "Mu", -2)
case('Mu_t')
call write_scalar(mu_t, "Mu_t", -2)
case('TKE')
call write_scalar(tk, "TKE", -2)
case('Omega')
call write_scalar(tw, "Omega", -2)
case('Kl')
call write_scalar(tkl, "Kl", -2)
case('tv')
call write_scalar(tv, "tv", -2)
case('tgm')
call write_scalar(tgm, "tgm", -2)
case('Wall_distance')
call write_scalar(dist, "Wall_dist", -2)
case('F1')
call write_scalar(sst_F1, "F1", -2)
case('Dudx')
call write_scalar(gradu_x ,"dudx ", 0)
case('Dudy')
call write_scalar(gradu_y ,"dudy ", 0)
case('Dudz')
call write_scalar(gradu_z ,"dudz ", 0)
case('Dvdx')
call write_scalar(gradv_x ,"dvdx ", 0)
case('Dvdy')
call write_scalar(gradv_y ,"dvdy ", 0)
case('Dvdz')
call write_scalar(gradv_z ,"dvdz ", 0)
case('Dwdx')
call write_scalar(gradw_x ,"dwdx ", 0)
case('Dwdy')
call write_scalar(gradw_y ,"dwdy ", 0)
case('Dwdz')
call write_scalar(gradw_z ,"dwdz ", 0)
case('DTdx')
call write_scalar(gradT_x ,"dTdx ", 0)
case('DTdy')
call write_scalar(gradT_y ,"dTdy ", 0)
case('DTdz')
call write_scalar(gradT_z ,"dTdz ", 0)
case('Dtkdx')
call write_scalar(gradtk_x,"dtkdx", 0)
case('Dtkdy')
call write_scalar(gradtk_y,"dtkdy", 0)
case('Dtkdz')
call write_scalar(gradtk_z,"dtkdz", 0)
case('Dtwdx')
call write_scalar(gradtw_x,"dtwdx", 0)
case('Dtwdy')
call write_scalar(gradtw_y,"dtwdy", 0)
case('Dtwdz')
call write_scalar(gradtw_z,"dtwdz", 0)
case('do not write')
! do not write
continue
case Default
print*, err//trim(control%w_list(n))//" to file"
end select
end do
end subroutine write_file