static_data_map_s.f90 Source File


Contents

Source Code


Source Code

submodule (static_data_m) static_data_map_s
    !! Handles static data related with map
    implicit none

contains

    module subroutine static_map_init(comm_handler, par_filepath, nf90_filepath)
        type(comm_handler_t), intent(in) :: comm_handler
        character(len=*), intent(in), optional :: par_filepath
        character(len=*), intent(in), optional :: nf90_filepath

        character(len=:), allocatable :: fpath
        logical :: write_to_netcdf = .false.
        integer :: idbg, nf90_stat, nf90_id
        character(len=3) :: plane_c
        character(len=:), allocatable :: nf90_fpath
        real(GP) :: phi_cano, phi_stag

        if (is_rank_info_writer) then
            write(get_stdout(),*)''
            write(get_stdout(),*)repeat('-',80)
            write(get_stdout(),*)'Initializing static data: map'
        endif

        if (.not.equi_is_initialized) then
            call handle_error('equi must be initialized at map initialization', &
                              GRILLIX_ERR_STATIC_DATA, __LINE__, __FILE__)
        endif
        if (.not.multigrid_is_initialized) then
            call handle_error('multigrid must be initialized at map initialization', &
                              GRILLIX_ERR_STATIC_DATA, __LINE__, __FILE__)
        endif

        ! Set debug output level
        if (enable_debug_output) then
            idbg = 1
        else
            idbg = 0
        endif

        if (present(par_filepath)) then
            fpath = par_filepath
        else
            fpath = static_data_parpath
        endif

        ! Select NetCDF file path and enable writing when required
        write(plane_c,'(I3.3)')comm_handler%get_plane()
        if (present(nf90_filepath)) then
            nf90_fpath = nf90_filepath
        else
            if (equi%is_axisymmetric()) then
                nf90_fpath = trim(static_data_dirpath)//'/map.nc'
                if ((cntrl_build_map == BUILD_AND_WRITE) .and. (is_rank_info_writer)) then
                    write_to_netcdf = .true.
                endif
            else
                nf90_fpath = trim(static_data_dirpath)//'/map_plane'//plane_c//'.nc'
                if (cntrl_build_map == BUILD_AND_WRITE) then
                    write_to_netcdf = .true.
                endif
            endif
        endif

        select case(cntrl_build_map)
            case(BUILD_NOWRITE, BUILD_AND_WRITE)
                if (is_rank_info_writer) then
                    write(get_stdout(),*)'Building map from scratch'
                endif
                call map%set_parameters(comm_handler, filename=fpath)
                call map%build(comm_handler, equi, mesh_cano, mesh_stag, dbgout=idbg)

                map_is_initialized = .true.

            case(READ)
                if (is_rank_info_writer) then
                    write(get_stdout(),*)'Reading map from file: ' // nf90_fpath
                endif

                nf90_stat = nf90_open(nf90_fpath, NF90_NETCDF4, nf90_id)
                call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)
                call map%read_netcdf(axisymmetric=equi%is_axisymmetric(), fgid=nf90_id)
                nf90_stat = nf90_close(nf90_id)
                call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)

                if (equi%is_axisymmetric()) then
                    ! Correct toroidal angles for axisymmetric case
                    phi_cano = mesh_cano%get_phi()
                    phi_stag = mesh_stag%get_phi()
                    call map%set_toroidal_location(comm_handler%get_plane(), phi_cano, phi_stag)
                endif

                map_is_initialized = .true.

            case(NONE)
                if (is_rank_info_writer) then
                    write(get_stdout(),*) 'Map not built because it was not requested by the user.'
                endif
            case default
                call handle_error('cntrl_build_map not valid', &
                                  GRILLIX_ERR_STATIC_DATA, __LINE__, __FILE__)
        end select

        if (map_is_initialized) then
#if ENABLE_CUDA
            call map%build_cuda_struct()
#endif
            if (enable_debug_output) then
                call map%display()
            endif
        endif

        if (write_to_netcdf) then
            if (is_rank_info_writer) then
                write(get_stdout(),*)'Writing map to file (series): ' // nf90_fpath
            endif
            nf90_stat = nf90_create(nf90_fpath, NF90_NETCDF4, nf90_id)
            call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)
            call map%write_netcdf(axisymmetric=equi%is_axisymmetric(),fgid=nf90_id)
            nf90_stat = nf90_close(nf90_id)
            call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)
        endif

        if (is_rank_info_writer) then
            write(get_stdout(),*)repeat('-',80)
            write(get_stdout(),*)''
        endif

    end subroutine

end submodule