state_template_netcdf_s.f90 Source File


Contents


Source Code

submodule(state_template_m) state_template_netcdf_s
    !! NetCDF I/O for state_template 
    use NETCDF
    use parallelisation_setup_m, only : is_rank_info_writer
    use error_handling_grillix_m, only: handle_error, handle_error_netcdf, error_info_t
    use status_codes_grillix_m, only : GRILLIX_ERR_OTHER, GRILLIX_ERR_CMD
    use screen_io_m, only : get_stdout
    implicit none

    integer, parameter :: nsnaps_limit = 99999
    ! Maximum number of allowed snapshots

contains

    module subroutine write_template_state_netcdf(comm_handler, dirname, filename)
        type(comm_handler_t), intent(in) :: comm_handler
        character(len=*), intent(in) , optional :: dirname
        character(len=*), intent(in) , optional :: filename

        character(len=:), allocatable :: dir
        character(len=:), allocatable :: snapsfile
        character(len=5) :: isnaps_c
        integer :: cmd_err, cmd_exit, ierr, nf90_stat, nf90_id
        integer :: id_nplanes_dim, id_np_max_dim, id_frames_dim, id_tau, id_tstep, id_dens, id_parmom, id_pion
        integer, dimension(2) :: starts, counts

        ! Create directory
        if (present(dirname)) then
            dir = dirname // '/'
            if (is_rank_info_writer) then
                call execute_command_line('mkdir -p ' // trim(dir), &
                                        exitstat=cmd_exit, cmdstat=cmd_err)
                if (cmd_err /= 0) then
                    call handle_error("Command execution failed: mkdir -p " // trim(dir), &
                                        GRILLIX_ERR_CMD, __LINE__, __FILE__, &
                                        error_info_t('cmd_err: ',[cmd_err]))
                endif
                if (cmd_exit /= 0) then
                    call handle_error("Command execution returned with warning: mkdir -p" &
                                        // trim(dir), &
                                        GRILLIX_ERR_CMD, __LINE__, __FILE__, &
                                        error_info_t('cmd_exit: ',[cmd_exit]))
                endif
            endif
            call MPI_barrier(comm_handler%get_comm_cart(), ierr)
        else
            dir = ''
        endif

        ! Assemble snapsfile path
        if (present(filename)) then
            snapsfile = dir // filename
        else
            isnaps = isnaps + 1
            if (isnaps > nsnaps_limit) then
                call handle_error("Maximum number of snapshots reached", &
                                        GRILLIX_ERR_OTHER, __LINE__, __FILE__, &
                                        error_info_t('isnaps: ',[isnaps]))
            endif
            write(isnaps_c,'(I5.5)')isnaps
            snapsfile = dir // 'snaps_' // isnaps_c // '.nc'
        endif

        ! Print info to screen
        if (is_rank_info_writer) then
            write(get_stdout(),*)'Writing state to file: ', snapsfile
        endif

        nf90_stat = nf90_create(snapsfile, NF90_NETCDF4 + NF90_NOCLOBBER, nf90_id, &
                                comm=comm_handler%get_comm_cart(), info=MPI_INFO_NULL)
        call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)

        ! Define dimensions
        nf90_stat = nf90_def_dim(nf90_id, 'nframes', 1, id_frames_dim)
        call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)

        nf90_stat = nf90_def_dim(nf90_id, 'nplanes', comm_handler%get_nplanes(), id_nplanes_dim)
        call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)
        
        nf90_stat = nf90_def_dim(nf90_id, 'np_max', np_max, id_np_max_dim)
        call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)

        ! Define variables
        nf90_stat = nf90_def_var(nf90_id, 'tau', NF90_DOUBLE, [id_frames_dim], id_tau)
        call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)

        nf90_stat = nf90_def_var(nf90_id, 'tstep', NF90_INT, [id_frames_dim], id_tstep)
        call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)

        nf90_stat = nf90_def_var(nf90_id, 'dens', NF90_DOUBLE, [id_np_max_dim, id_nplanes_dim], id_dens)
        call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)

        nf90_stat = nf90_def_var(nf90_id, 'parmom', NF90_DOUBLE, [id_np_max_dim, id_nplanes_dim], id_parmom)
        call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)

        nf90_stat = nf90_def_var(nf90_id, 'pion', NF90_DOUBLE, [id_np_max_dim, id_nplanes_dim], id_pion)
        call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)

        nf90_stat = nf90_enddef(nf90_id)
        call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)

        ! Put variables
        if (is_rank_info_writer) then
            nf90_stat = nf90_put_var(nf90_id, id_tau, [tau])
            call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)

            nf90_stat = nf90_put_var(nf90_id, id_tstep, [tstep])
            call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)
        endif

        starts = [1, comm_handler%get_plane()+1]
        counts = [np_max, 1]

        nf90_stat = nf90_put_var(nf90_id, id_dens, dens, start=starts, count=counts)
        call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)

        nf90_stat = nf90_put_var(nf90_id, id_parmom, parmom, start=starts, count=counts)
        call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)

        nf90_stat = nf90_put_var(nf90_id, id_pion, pion, start=starts, count=counts)
        call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)

        ! Close file
        nf90_stat = nf90_close(nf90_id)
        call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)

    end subroutine

    module subroutine read_template_state_netcdf(comm_handler, dirname, filename)
        type(comm_handler_t), intent(in) :: comm_handler
        character(len=*), intent(in) , optional :: dirname
        character(len=*), intent(in) , optional :: filename

        character(len=:), allocatable :: dir
        character(len=:), allocatable :: snapsfile
        character(len=5) :: isnaps_c
        logical :: fexist
        integer :: nf90_stat, nf90_id
        integer :: id_tau, id_tstep, id_dens, id_parmom, id_pion
        integer, dimension(2) :: starts, counts

        ! Create directory
        if (present(dirname)) then
            dir = dirname // '/'
        else
            dir = ''
        endif

        ! Assemble snapsfile path
        if (present(filename)) then
            snapsfile = dir // filename
            isnaps = 0
        else
            do isnaps = nsnaps_limit, 1, -1
                write(isnaps_c,'(I5.5)')isnaps
                snapsfile = dir // 'snaps_' // isnaps_c // '.nc'
                inquire(file=snapsfile, exist=fexist)
                if (fexist) then
                    exit
                endif
            enddo
            ! Check if any snapshotfile was found
            if (isnaps == 0) then
                call handle_error('No snapsfile for continuation found', GRILLIX_ERR_OTHER, __LINE__, __FILE__)
            endif
        endif

        ! Print info to screen
        if (is_rank_info_writer) then
            write(get_stdout(),*)'Reading state from file: ', snapsfile
        endif

        nf90_stat = nf90_open(snapsfile, NF90_NETCDF4, nf90_id, &
                              comm=comm_handler%get_comm_cart(), info=MPI_INFO_NULL)
        call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)

        nf90_stat = nf90_inq_varid(nf90_id, 'tau', id_tau)
        call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)
        nf90_stat = nf90_get_var(nf90_id, id_tau, tau)
        call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)

        nf90_stat = nf90_inq_varid(nf90_id, 'tstep', id_tstep)
        call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)
        nf90_stat = nf90_get_var(nf90_id, id_tstep, tstep)
        call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)

        starts = [1, comm_handler%get_plane()+1]
        counts = [np_max, 1]

        nf90_stat = nf90_inq_varid(nf90_id, 'dens', id_dens)
        call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)
        nf90_stat = nf90_get_var(nf90_id, id_dens, dens, start=starts, count=counts)
        call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)

        nf90_stat = nf90_inq_varid(nf90_id, 'parmom', id_parmom)
        call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)
        nf90_stat = nf90_get_var(nf90_id, id_parmom, parmom, start=starts, count=counts)
        call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)

        nf90_stat = nf90_inq_varid(nf90_id, 'pion', id_pion)
        call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)
        nf90_stat = nf90_get_var(nf90_id, id_pion, pion, start=starts, count=counts)
        call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)

        nf90_stat = nf90_close(nf90_id)
        call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)

    end subroutine

end submodule