static_data_polars_s.f90 Source File


Contents


Source Code

submodule (static_data_m) static_data_polars_s
    !! Handles static data related with polars
    implicit none

contains

    module subroutine static_polars_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, nf90_grpid
        character(len=3) :: plane_c
        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: polars'
        endif

        if (.not.equi_is_initialized) then
            call handle_error('equi must be initialized at polars initialization', &
                              GRILLIX_ERR_STATIC_DATA, __LINE__, __FILE__)
        endif
        if (.not.multigrid_is_initialized) then
            call handle_error('multigrid must be initialized at polars 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)//'/polars.nc'
                if ((cntrl_build_polars == BUILD_AND_WRITE) .and. (is_rank_info_writer)) then
                    write_to_netcdf = .true.
                endif
            else
                nf90_fpath = trim(static_data_dirpath)//'/polars_plane'//plane_c//'.nc'
                if (cntrl_build_polars == BUILD_AND_WRITE) then
                    write_to_netcdf = .true.
                endif
            endif
        endif

        ! Associate staggered type
        if (equi%is_axisymmetric()) then
            polars_stag => polars_cano
        else
            allocate(polars_stag_mem)
            polars_stag => polars_stag_mem
        endif

        select case(cntrl_build_polars)
            case(BUILD_NOWRITE, BUILD_AND_WRITE)
                if (is_rank_info_writer) then
                    write(get_stdout(),*)'Building polars from scratch'
                endif

                call polars_cano%set_parameters(filename=fpath)
                call polars_cano%build(equi, mesh_cano, dbgout=idbg)
                if (.not.equi%is_axisymmetric()) then
                    call polars_stag%set_parameters(filename=fpath)
                    call polars_stag%build(equi, mesh_stag, dbgout=idbg)
                endif

                polars_is_initialized = .true.

            case(READ)
                if (is_rank_info_writer) then
                    write(get_stdout(),*)'Reading polars data from file: ' // nf90_fpath
                endif
                nf90_stat = nf90_open(nf90_fpath, NF90_NETCDF4, nf90_id)
                call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)
                
                ! Read polars from trunk
                if (equi%is_axisymmetric()) then
                    call polars_cano%read_netcdf(nf90_id)
                else
                    nf90_stat = nf90_inq_grp_ncid(nf90_id, 'canonical', nf90_grpid)
                    call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)
                    call polars_cano%read_netcdf(nf90_grpid)

                    nf90_stat = nf90_inq_grp_ncid(nf90_id, 'staggered', nf90_grpid)
                    call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)
                    call polars_stag_mem%read_netcdf(nf90_grpid)
                endif

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

                polars_is_initialized = .true.

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

        if ((enable_debug_output) .and. (polars_is_initialized)) then
            call polars_cano%display()
        endif

        if (write_to_netcdf) then
            if (is_rank_info_writer) then
                write(get_stdout(),*)'Writing polars data to file (series): ' // nf90_fpath
            endif

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

            if (equi%is_axisymmetric()) then
                call polars_cano%write_netcdf(equi, nf90_id)
            else
                nf90_stat = nf90_def_grp(nf90_id, 'canonical', nf90_grpid)
                call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)
                call polars_cano%write_netcdf(equi, nf90_grpid)

                nf90_stat = nf90_def_grp(nf90_id, 'staggered', nf90_grpid)
                call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)
                call polars_stag%write_netcdf(equi, nf90_grpid)
            endif

            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