static_data_multigrid_s.f90 Source File


Contents


Source Code

submodule (static_data_m) static_data_multigrid_s
    !! Handles static data related with multigrid
    implicit none

    ! Parameters for multigrid meshes
    real(GP) :: spacing_f = GP_NAN
    !! Mesh spacing on finest level
    integer :: nlvls = 3
    !! Number of multigrid levels
    integer :: size_neighbor = 2
    !! Number of neighbor points stored
    integer :: size_ghost_layer = 2
    !! Depth of ghost layer
    integer :: reorder_size = 8
    !! Block size for reordering
    logical :: extend_beyond_wall = .false.
    !! Switch if mesh is extended beyond wall

    namelist / params_multigrid / spacing_f, nlvls, size_neighbor, &
        size_ghost_layer, reorder_size, extend_beyond_wall

contains

    subroutine read_params_multigrid(par_filepath)
        character(len=*), intent(in), optional :: par_filepath

        character(len=:), allocatable :: fpath

        integer :: io_unit, io_error
        character(len=256) :: io_errmsg

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

        ! Read multigrid parameters
        open(newunit=io_unit, file=fpath, status='old', action='read', iostat=io_error, iomsg=io_errmsg)
        if (io_error /= 0) then
            call handle_error(io_errmsg, GRILLIX_ERR_NAMELIST, __LINE__, __FILE__)
        endif
        read(io_unit, nml=params_multigrid, iostat=io_error, iomsg=io_errmsg)
        if (io_error /= 0) then
            call handle_error(io_errmsg, GRILLIX_ERR_NAMELIST, __LINE__, __FILE__)
        endif
        close(io_unit)

    end subroutine

    subroutine display_params_multigrid()

        if (is_rank_info_writer) then
            write(get_stdout(), nml=params_multigrid)
        endif

    end subroutine

    module subroutine static_multigrid_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=3) :: plane_c
        real(GP) :: dphi, phi_cano, phi_stag
        integer :: idbg
        integer, allocatable, dimension(:) :: mgr_reorder_size
        logical :: write_to_netcdf = .false.

        integer :: nf90_id, nf90_stat, nf90_grpid
        character(len=:), allocatable :: nf90_fpath

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

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

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

        ! Read multigrid parameters
        if (present(par_filepath)) then
            call read_params_multigrid(par_filepath)
        else
            call read_params_multigrid()
        endif
        call display_params_multigrid()

        ! Toroidal location of mesh
        dphi = TWO_PI / comm_handler%get_nplanes()
        phi_cano = dphi * comm_handler%get_plane()
        phi_stag = phi_cano + dphi / 2.0_GP
        write(plane_c,'(I3.3)')comm_handler%get_plane()

        ! Select NetCDF file path and enable writing when required
        if (present(nf90_filepath)) then
            nf90_fpath = nf90_filepath
        else
            if (equi%is_axisymmetric()) then
                nf90_fpath = trim(static_data_dirpath)//'/multigrid.nc'
                if ((cntrl_build_multigrid == BUILD_AND_WRITE) .and. (is_rank_info_writer)) then
                    write_to_netcdf = .true.
                endif
            else
                nf90_fpath = trim(static_data_dirpath)//'/multigrids_plane'//plane_c//'.nc'
                if (cntrl_build_multigrid == BUILD_AND_WRITE) then
                    write_to_netcdf = .true.
                endif
            endif
        endif

        select case(cntrl_build_multigrid)
            case(BUILD_NOWRITE, BUILD_AND_WRITE)
                if (is_rank_info_writer) then
                    write(get_stdout(),*)'Building multigrid data from scratch'
                endif
                allocate(mgr_reorder_size(nlvls), source=reorder_size)
                call multigrid_cano%init(equi, phi_cano, &
                                        nlvls, spacing_f, size_neighbor, size_ghost_layer, &
                                        mesh_cano, mgr_reorder_size, extend_beyond_wall, idbg)

                call multigrid_stag%init(equi, phi_stag, &
                                        nlvls, spacing_f, size_neighbor, size_ghost_layer, &
                                        mesh_stag, mgr_reorder_size, extend_beyond_wall, idbg)

                deallocate(mgr_reorder_size)

                multigrid_is_initialized = .true.

            case(READ)
                if (is_rank_info_writer) then
                    write(get_stdout(),*)'Reading multigrid data from file (series)' // nf90_fpath
                endif

                nf90_stat = nf90_open(nf90_fpath, NF90_NETCDF4, nf90_id)
                call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)

                if (equi%is_axisymmetric()) then
                    ! Read multigrid and mesh from single file and correct phi angle
                    call multigrid_cano%read_netcdf(nf90_id, mesh_cano)
                    call multigrid_stag%read_netcdf(nf90_id, mesh_stag)

                    call multigrid_cano%set_phi_coarse(phi_cano)
                    call mesh_cano%set_phi(phi_cano)
                    call multigrid_stag%set_phi_coarse(phi_stag)
                    call mesh_stag%set_phi(phi_stag)
                    call mesh_stag%set_phi(phi_stag)
                else
                    ! Read Multigrid and mesh, canonical and staggered from different NetCDF groups
                    nf90_stat = nf90_inq_grp_ncid(nf90_id, 'canonical', nf90_grpid)
                    call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)
                    call multigrid_cano%read_netcdf(nf90_grpid, mesh_cano)
                    nf90_stat = nf90_inq_grp_ncid(nf90_id, 'staggered', nf90_grpid)
                    call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)
                    call multigrid_stag%read_netcdf(nf90_grpid, mesh_stag)
                endif

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

                ! Create Cartesian index arrays used to accelerate cart_to_mesh_index for further build stages
                call mesh_cano%build_shaped_cart_arrays()
                call mesh_stag%build_shaped_cart_arrays()

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

        if (write_to_netcdf) then
            if (is_rank_info_writer) then
                write(get_stdout(),*)'Writing multigrid data to file (series): ' // nf90_fpath
            endif
            if (equi%is_axisymmetric()) then
                nf90_stat = nf90_create(nf90_fpath, NF90_NETCDF4, nf90_id)
                call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)
                call multigrid_cano%write_netcdf(nf90_id)
                nf90_stat = nf90_close(nf90_id)
                call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)
            else
                !For 3D equilibria all planes are written out, canonical and staggered to separate groups
                nf90_stat = nf90_create(nf90_fpath, NF90_NETCDF4, nf90_id)
                call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)
                nf90_stat = nf90_def_grp(nf90_id, 'canonical', nf90_grpid)
                call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)
                call multigrid_cano%write_netcdf(nf90_grpid)
                nf90_stat = nf90_def_grp(nf90_id, 'staggered', nf90_grpid)
                call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)
                call multigrid_stag%write_netcdf(nf90_grpid)
                nf90_stat = nf90_close(nf90_id)
                call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)
            endif
        endif

        if (is_rank_info_writer) then
            if ((enable_debug_output) .and. (multigrid_is_initialized)) then
                call mesh_cano%display()
            endif
            write(get_stdout(),*)repeat('-',80)
            write(get_stdout(),*)''
        endif

    end subroutine

end submodule