Skip to content
Snippets Groups Projects
obl_main.f90 11.2 KiB
Newer Older
  • Learn to ignore specific revisions
  • #include "obl_def.fi"
    
    
    Daria Gladskikh's avatar
    Daria Gladskikh committed
        !< @brief main program for calculations for ocean boundary layer
    
    #ifdef USE_CONFIG_PARSER
    
        use iso_c_binding, only: C_NULL_CHAR
        use config_parser
    
    Daria Gladskikh's avatar
    Daria Gladskikh committed
        ! modules used
    
    Evgeny Mortikov's avatar
    Evgeny Mortikov committed
        use obl_grid
    
        use obl_state
    
        use obl_state_eq
    
        use obl_common_phys
    
        use obl_tforcing
    
        use obl_tslice
        use obl_tseries
    
    Evgeny Mortikov's avatar
    Evgeny Mortikov committed
        use obl_output
    
    Evgeny Mortikov's avatar
    Evgeny Mortikov committed
        use obl_fluxes
    
        use obl_diag
    
        use obl_pph, only: pphParamType, &
            define_stability_functions_pph => define_stability_functions, & 
            advance_turbulence_closure_pph => advance_turbulence_closure
        use obl_pph_dyn, only: pphDynParamType, &
            define_stability_functions_pph_dyn => define_stability_functions, & 
            advance_turbulence_closure_pph_dyn => advance_turbulence_closure
        use obl_k_epsilon, only: kepsilonParamType, &
            define_stability_functions_k_epsilon => define_stability_functions, & 
            advance_turbulence_closure_k_epsilon => advance_turbulence_closure, &
            TKE_init, EPS_init, TKE_bc, EPS_bc
    
    Daria Gladskikh's avatar
    Daria Gladskikh committed
        !use vertical_mixing, default = off
        !use vermix
    
    Daria Gladskikh's avatar
    Daria Gladskikh committed
        ! directives list
    
        type(gridDataType) :: grid
    
    
        ! turbulence closure parameters
    
        type(pphParamType) :: param_pph
        type(pphDynParamType) :: param_pph_dyn
    
    Evgeny Mortikov's avatar
    Evgeny Mortikov committed
        type(kepsilonParamType) :: param_k_epsilon
    
        !< boundary conditions data
        type(oblBcType) :: bc
    
    
        real :: domain_height
        integer :: grid_cz
    
    
        !< output 
        type(oblOutputStx) :: scm_output
    
    
    Daria Gladskikh's avatar
    Daria Gladskikh committed
        real :: dt !< time step [s] 
        integer :: i, k  !< counters          
        integer :: status, num !< for file input/output
        real :: time_begin, time_end, time_current  !< time for output
    
    Daria Gladskikh's avatar
    Daria Gladskikh committed
        integer :: closure_mode  !< closure type:
                                 !1 - pacanowski-philander, 2 - pacanowski-philander+, 
                                 !3 - k-epsilon explicit, 4 - k-epsilon semiimplicit, 5 - inmom
    
    Evgeny Mortikov's avatar
    Evgeny Mortikov committed
        integer, parameter :: output_mode = 3   ! 1 -- netcdf, 2 -- ascii, 3 -- tecplot
    
    Evgeny Mortikov's avatar
    Evgeny Mortikov committed
    
    
        integer :: obl_setup     ! 1 - Kato-Phillips, 2 - Papa station (fluxes), 3 - Papa station (meteo)
    
    
        real, parameter :: Cd0 = 0.001           ! bottom drag coefficient
    
    Daria Gladskikh's avatar
    Daria Gladskikh committed
        real :: mld !< mixed layer depth, [m]
        real :: lab_mld !< theoretical mixed layer depth, [m]
    
    Daria Gladskikh's avatar
    Daria Gladskikh committed
        real :: N2_0
    
        
        ! command line arguments
        ! --------------------------------------------------------------------------------
        integer :: num_args
        character(len = 128) :: arg
    
        character(len = 128), parameter :: arg_key_help = '--help'
    
        character(len = 128), parameter :: arg_key_setup = '--setup'
    
        character(len = 128), parameter :: arg_key_config = "--config"
    
        integer :: ierr
        ! --------------------------------------------------------------------------------
    
    
    Evgeny Mortikov's avatar
    Evgeny Mortikov committed
        ! screen output parameters
        integer, parameter :: nscreen = 1000
    
        ! file output parameters
        integer, parameter :: noutput = 60
    
    
        obl_setup = setup_kato !< ! 0 - Kato-Phillips (default), 1 - Papa station (fluxes), 2 - Papa station (meteo)
    
    Daria Gladskikh's avatar
    Daria Gladskikh committed
        closure_mode = 4 !< 1 - pacanowski-philander, 2 - pacanowski-philander+, 3 - k-epsilon explicit, 4 - k-epsilon semiimplicit, 5 - inmom
    
    
    
        ! --- command line arguments processing
        num_args = command_argument_count()
        do i = 1, num_args
    
            call get_command_argument(i, arg)
            if (trim(arg) == trim(arg_key_help)) then
                write(*, *) ' obl model, usage:'
                write(*, *) ' --help'
                write(*, *) '    print usage options'
    
                write(*, *) ' --setup [key]'
                write(*, *) '    key = kato (default) || papa-fluxes || papa'
    
                write(*, *) ' --config [filename]'
                write(*, *) '    use configuration file'
                return
            end if
    
    
            if (trim(arg) == trim(arg_key_setup)) then
                if (i == num_args) then
                    write(*, *) ' FAILURE! > missing setup [key] argument'
                    ierr = 1        ! signal ERROR
                    return
                end if
                
                call get_command_argument(i + 1, arg)
                obl_setup = get_setup_id(arg)
                if (obl_setup == -1) then
                    write(*, *) ' FAILURE! > unknown setup [key]: ', trim(arg)
                    ierr = 1        ! signal ERROR
                    return
                end if
            else if (trim(arg) == trim(arg_key_config)) then
    
                if (i == num_args) then
                    write(*, *) ' FAILURE! > missing configuration file [key] argument'
                    ierr = 1        ! signal ERROR
                    return
                end if
            
                call get_command_argument(i + 1, arg)
    #ifdef USE_CONFIG_PARSER
                call c_config_run(trim(arg)//C_NULL_CHAR, status)
                if (status == 0) then
                    write(*, *) ' FAILURE! > unable to parse configuration file: ', trim(arg)
                    ierr = 1        ! signal ERROR
                    return
                end if
    
    
                !< forcing configuration file setup
                obl_setup = -1
    
    Daria Gladskikh's avatar
    Daria Gladskikh committed
    
    
        ! ---------------------------------------------------------------------------- 
    
        call set_grid(grid, obl_setup, ierr)
        if (ierr /= 0) then
            return
        endif 
    
        !< debug grid print
    
    Evgeny Mortikov's avatar
    Evgeny Mortikov committed
        call print_grid(grid)
    
        ! ---------------------------------------------------------------------------- 
    
        !< setting model time
        ! ----------------------------------------------------------------------------  
        call set_time(time_begin, time_end, dt, obl_setup, ierr)
        if (ierr /= 0) then
            return
        endif 
        time_current = time_begin
        ! ----------------------------------------------------------------------------  
    
    Daria Gladskikh's avatar
    Daria Gladskikh committed
        ! allocation
    
        call allocate_state_vec(grid%cz)
    
        ! initialize scm
        call init_scm_vec(grid%cz)
    
    
        !< setting phys
        ! ----------------------------------------------------------------------------
    
        call set_phys(obl_setup, ierr)
        if (ierr /= 0) then
            return
    
        endif 
        ! ----------------------------------------------------------------------------
    
    Evgeny Mortikov's avatar
    Evgeny Mortikov committed
    
    
        !< setting forcing
        ! ---------------------------------------------------------------------------- 
        call set_forcing(obl_setup, ierr)
        if (ierr /= 0) then
            return
        endif 
        ! ----------------------------------------------------------------------------
    
        !< initialization of main fields
        ! ---------------------------------------------------------------------------- 
    
        call set_initial_conditions(grid, obl_setup, ierr)
        if (ierr /= 0) then
            return
    
        endif
    
        Theta_dev = Theta - Theta_ref
        Salin_dev = Salin - Salin_ref
    
        call get_density(Rho, Theta_dev, Salin_dev, grid%cz)
    
        ! ----------------------------------------------------------------------------
    
        !< initialization of turbulence closure
    
        if (closure_mode.eq.3 .or. closure_mode.eq.4) then
    
    Evgeny Mortikov's avatar
    Evgeny Mortikov committed
            call TKE_init(TKE, param_k_epsilon, grid%cz)
            call eps_init(EPS, param_k_epsilon, grid%cz, grid%height)
    
        ! ---------------------------------------------------------------------------- 
    
    Evgeny Mortikov's avatar
    Evgeny Mortikov committed
        !< setting reference data
        ! ----------------------------------------------------------------------------
        N2_0 = 0.00044
        ! ----------------------------------------------------------------------------
    
        do while (time_current < time_end )
            ! ----------------------------------------------------------------------------
    
            
            !< define fluxes & dynamic scales [surface]
            ! ----------------------------------------------------------------------------
    
    Evgeny Mortikov's avatar
    Evgeny Mortikov committed
            call advance_surface_fluxes(bc, time_current, grid)
    
            ! ----------------------------------------------------------------------------
    
            !< define fluxes & dynamic scales [bottom]
            ! ----------------------------------------------------------------------------
    
    Evgeny Mortikov's avatar
    Evgeny Mortikov committed
            call advance_bottom_fluxes(bc, time_current, grid)
    
            ! ----------------------------------------------------------------------------
    
    
            !< advance turbulence closure
            ! ----------------------------------------------------------------------------
    
            if (closure_mode.eq.1) 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
    
                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
    
                call define_stability_functions_k_epsilon(param_k_epsilon, bc, grid)
                call advance_turbulence_closure_k_epsilon(param_k_epsilon, bc, grid, dt)
    
            ! ----------------------------------------------------------------------------
    
            !< advance single-column model eq.
    
            ! ----------------------------------------------------------------------------
    
            call advance_scm_eq(bc, grid, dt)
    
            ! ----------------------------------------------------------------------------
    
    Evgeny Mortikov's avatar
    Evgeny Mortikov committed
            !> advance time
    
    Evgeny Mortikov's avatar
    Evgeny Mortikov committed
            i = i + 1
            time_current = time_current + dt
    
    
    Evgeny Mortikov's avatar
    Evgeny Mortikov committed
            !> advance screen output
    
    Evgeny Mortikov's avatar
    Evgeny Mortikov committed
            if (mod(i, nscreen) == 0) then
    
                call get_mld(mld, N2, grid%dz, grid%cz)
                call get_mld_ref(lab_mld, bc%U_dynH, N2_0, time_current, grid%height)
    
    Evgeny Mortikov's avatar
    Evgeny Mortikov committed
    
    
                ! *: add finite check
    
    
    Daria Gladskikh's avatar
    Daria Gladskikh committed
                write(*, '(a,g0)') ' mld = ', mld
    
    Evgeny Mortikov's avatar
    Evgeny Mortikov committed
                write(*, '(a,g0)') ' Theta(surface) = ', Theta_dev(grid%cz) + Theta_ref
                write(*, '(a,g0,a,g0,a)') ' current time = ', time_current / 3600.0, ' HRS [ ', & 
                    (time_current / time_end) * 100.0, '% ]' 
                write(*, '(a)') '-------------------------------------------------'
    
    Evgeny Mortikov's avatar
    Evgeny Mortikov committed
            !> advance file output
    
            if (mod(i, noutput) == 0) then
    
                call push_state_vec(scm_output, grid%cz)
    
                call get_mld(mld, N2, grid%dz, grid%cz)
                call get_mld_ref(lab_mld, bc%U_dynH, N2_0, time_current, grid%height)
    
                call push_value_to_tseries(scm_output%mld, mld)
                call push_value_to_tseries(scm_output%mld_ref, lab_mld)
    
                call push_value_to_tseries(scm_output%tau_x, bc%u_momentum_fluxH * Rho_ref)
                call push_value_to_tseries(scm_output%tau_y, bc%v_momentum_fluxH * Rho_ref)
    
    Evgeny Mortikov's avatar
    Evgeny Mortikov committed
    
    
                call push_value_to_tseries(scm_output%Theta_surf, Theta_dev(grid%cz) + Theta_ref)
    
                call push_value_to_tseries(scm_output%time, time_current / 3600.0)
    
    Evgeny Mortikov's avatar
    Evgeny Mortikov committed
        if (output_mode.eq.1) then
            write(*, *) ' >> writing netcdf output ...'
    
            call write_netcdf(scm_output, grid%z)
    
    Evgeny Mortikov's avatar
    Evgeny Mortikov committed
        endif
    
    Evgeny Mortikov's avatar
    Evgeny Mortikov committed
        if (output_mode.eq.2) then
            write(*, *) ' >> writing ascii output ...'
    
            call write_ascii(scm_output, grid%z)
    
    Evgeny Mortikov's avatar
    Evgeny Mortikov committed
        endif
    
    
        if (output_mode.eq.3) then
            write(*, *) ' >> writing tecplot output ...'
    
            call write_tecplot(scm_output, grid%z)
    
    Evgeny Mortikov's avatar
    Evgeny Mortikov committed
        endif
    
        !> deallocate state
    
        call deallocate_state_vec
    
    
        !> deallocate scm
        call deallocate_scm_vec
    
    Evgeny Mortikov's avatar
    Evgeny Mortikov committed
        !> deallocate time-dependent forcing
        call deallocate_fluxes_data
    
    Evgeny Mortikov's avatar
    Evgeny Mortikov committed
        !> removing time slice data
    
        call output_cleanup(scm_output)
    
    Evgeny Mortikov's avatar
    Evgeny Mortikov committed
    
    
    Evgeny Mortikov's avatar
    Evgeny Mortikov committed
        ! > removing grid data
        call deallocate_grid(grid)