Skip to content
Snippets Groups Projects
Commit 23144121 authored by Evgeny Mortikov's avatar Evgeny Mortikov
Browse files

adding turbulence model def.

parent 48469f17
Branches
No related tags found
No related merge requests found
......@@ -18,7 +18,7 @@ module obl_config
public
!> @brief setup enum def.
!> @brief config enum def.
integer, parameter :: obl_config_kato = 0 !< Kato-Phillips setup
integer, parameter :: obl_config_papa_fluxes = 1 !< Papa-station (fluxes) setup
integer, parameter :: obl_config_papa_meteo = 2 !< Papa-station (meteo) setup
......@@ -27,6 +27,15 @@ module obl_config
character(len = 16), parameter :: obl_config_papa_fluxes_tag = 'papa-fluxes'
character(len = 16), parameter :: obl_config_papa_meteo_tag = 'papa-meteo'
!> @brief model enum def.
integer, parameter :: obl_model_pph = 0 !< pacanowski-philander
integer, parameter :: obl_model_pph_dyn = 1 !< pacanowski-philander (dynamic)
integer, parameter :: obl_model_k_epsilon = 2 !< k-epsilon
character(len = 16), parameter :: obl_model_pph_tag = 'pph'
character(len = 16), parameter :: obl_model_pph_dyn_tag = 'pph-dyn'
character(len = 16), parameter :: obl_model_k_epsilon_tag = 'k-epsilon'
contains
......@@ -47,7 +56,6 @@ contains
end function
! --------------------------------------------------------------------------------
function get_obl_config_tag(id) result(tag)
implicit none
integer :: id
......@@ -64,6 +72,39 @@ contains
end function
! --------------------------------------------------------------------------------
function get_obl_model_id(tag) result(id)
implicit none
character(len=*), intent(in) :: tag
integer :: id
id = - 1
if (trim(tag) == trim(obl_model_pph_tag)) then
id = obl_model_pph
else if (trim(tag) == trim(obl_model_pph_dyn_tag)) then
id = obl_model_pph_dyn
else if (trim(tag) == trim(obl_model_k_epsilon_tag)) then
id = obl_model_k_epsilon
end if
end function
function get_obl_model_tag(id) result(tag)
implicit none
integer :: id
character(len=:), allocatable :: tag
tag = 'undefined'
if (id == obl_model_pph) then
tag = obl_model_pph_tag
else if (id == obl_model_pph_dyn) then
tag = obl_model_pph_dyn_tag
else if (id == obl_model_k_epsilon) then
tag = obl_model_k_epsilon_tag
end if
end function
! --------------------------------------------------------------------------------
subroutine set_grid(grid, config_id, ierr)
!> @brief grid parameters setup
......
......@@ -53,12 +53,10 @@ program obl_main
! = obl_config_papa_fluxes: Papa-station (fluxes)
! = obl_config_papa_meteo: Papa-station (meteo)
integer :: closure_mode ! --- OBL closure def.
! = 1 - pacanowski-philander
! = 2 - pacanowski-philander dynamic
! = 3 - k-epsilon (explicit)
! = 4 - k-epsilon (semiimplicit)
! = 5 - inmom
integer :: obl_model_id ! --- OBL model def.
! = obl_model_pph: pacanowski-philander
! = obl_model_pph_dyn: pacanowski-philander (dynamic)
! = obl_model_k_epsilon: k-epsilon
type(gridDataType) :: grid
......@@ -102,6 +100,7 @@ program obl_main
character(len = 128), parameter :: arg_key_help = '--help'
character(len = 128), parameter :: arg_key_config = "--config"
character(len = 128), parameter :: arg_key_model = "--model"
integer :: i, status
integer :: ierr
......@@ -110,18 +109,17 @@ program obl_main
!< default config = obl_config_kato (Kato-Phiilps)
!< possible values:
!< obl_config_kato
!< obl_config_papa_fluxes
!< obl_config_papa_meteo
!< = obl_config_kato
!< = obl_config_papa_fluxes
!< = obl_config_papa_meteo
obl_config_id = obl_config_kato
!< default closure = 4, k-epsilon (implicit)
!< default closure = 4, k-epsilon
!< poissible values:
!< = 1, Pacanowski-Philander
!< = 2, Pacanowski-Philander (dynamic)
!< = 3, k-epsilon (explicit) *: DEPRECATED
!< = 4, k-epsilon (implicit)
closure_mode = 4
!< = obl_model_pph, pacanowski-philander
!< = obl_model_pph_dyn, pacanowski-philander (dynamic)
!< = obl_model_k_epsilon, k-epsilon
obl_model_id = obl_model_k_epsilon
!< default output = 1, (netcdf)
!< possible values:
......@@ -140,8 +138,10 @@ program obl_main
write(*, *) ' --help'
write(*, *) ' print usage options'
write(*, *) ' --config [key] || [filename]'
write(*, *) ' predefined setup [key] = kato (default) || papa-fluxes || papa'
write(*, *) ' builtin setup [key] = kato (default) || papa-fluxes || papa'
write(*, *) ' use configuration file [filename]'
write(*, *) ' --model [key]'
write(*, *) ' [key] = pph || pph-dyn || k-epsilon (default)'
return
end if
......@@ -172,6 +172,20 @@ program obl_main
return
#endif
end if
else if (trim(arg) == trim(arg_key_model)) then
if (i == num_args) then
write(*, *) ' FAILURE! > missing model [key] argument'
ierr = 1 ! signal ERROR
return
end if
call get_command_argument(i + 1, arg)
obl_model_id = get_obl_model_id(arg)
if (obl_model_id == -1) then
write(*, *) ' FAILURE! > unknown model [key]: ', trim(arg)
ierr = 1 ! signal ERROR
return
end if
endif
enddo
......@@ -243,13 +257,13 @@ program obl_main
call advance_surface_fluxes(bc, time_begin, grid)
call advance_bottom_fluxes(bc, time_begin, grid)
if (closure_mode.eq.1) then
if (obl_model_id.eq.obl_model_pph) then
call define_stability_functions_pph(param_pph, bc, grid)
call init_turbulence_closure_pph(param_pph, bc, grid)
else if (closure_mode.eq.2) then
else if (obl_model_id.eq.obl_model_pph_dyn) then
call define_stability_functions_pph_dyn(param_pph_dyn, bc, grid)
call init_turbulence_closure_pph_dyn(param_pph_dyn, bc, grid)
else if (closure_mode.eq.4) then
else if (obl_model_id.eq.obl_model_k_epsilon) then
call define_stability_functions_k_epsilon(param_k_epsilon, bc, grid)
call init_turbulence_closure_k_epsilon(param_k_epsilon, bc, grid)
endif
......@@ -272,13 +286,13 @@ program obl_main
!< advance turbulence closure
! ----------------------------------------------------------------------------
if (closure_mode.eq.1) then
if (obl_model_id.eq.obl_model_pph) then
call define_stability_functions_pph(param_pph, bc, grid)
call advance_turbulence_closure_pph(param_pph, bc, grid, dt)
else if (closure_mode.eq.2) then
else if (obl_model_id.eq.obl_model_pph_dyn) then
call define_stability_functions_pph_dyn(param_pph_dyn, bc, grid)
call advance_turbulence_closure_pph_dyn(param_pph_dyn, bc, grid, dt)
else if (closure_mode.eq.4) then
else if (obl_model_id.eq.obl_model_k_epsilon) then
call define_stability_functions_k_epsilon(param_k_epsilon, bc, grid)
call advance_turbulence_closure_k_epsilon(param_k_epsilon, bc, grid, dt)
endif
......
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Please register or to comment