static_data_masks_s.f90 Source File


Contents


Source Code

submodule (static_data_m) static_data_masks_s
    !! Handles static data related with masks
    implicit none

contains

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

        logical :: write_to_netcdf = .false.
        character(len=3) :: plane_c
        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: masks'
        endif

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

        ! Associate staggered type
        if (equi%is_axisymmetric()) then
            masks_stag => masks_cano
        else
            allocate(masks_stag_mem)
            masks_stag => masks_stag_mem
        endif

        select case(cntrl_build_masks)
            case(BUILD_NOWRITE, BUILD_AND_WRITE)

                if (parbnd_immersed_is_initialized) then
                    call masks_cano%build(equi, mesh_cano, parbnd_immersed_cano)
                else
                    call masks_cano%build(equi, mesh_cano)
                endif
                if (.not.equi%is_axisymmetric()) then
                    if (parbnd_immersed_is_initialized) then
                        call masks_stag%build(equi, mesh_stag, parbnd_immersed_stag)
                    else
                        call masks_cano%build(equi, mesh_stag)
                    endif
                endif

                masks_is_initialized = .true.

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

                if (equi%is_axisymmetric()) then
                    call masks_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 masks_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 masks_stag%read_netcdf(nf90_grpid)
                endif
                nf90_stat = nf90_close(nf90_id)
                call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)

                masks_is_initialized = .true.

            case(NONE)
                if (is_rank_info_writer) then
                    write(get_stdout(),*) 'masks not built because it was not requested by the user.'
                endif
            case default
                call handle_error('cntrl_build_masks 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 masks 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 masks_cano%write_netcdf(nf90_id)
            else
                nf90_stat = nf90_def_grp(nf90_id, 'canonical', nf90_grpid)
                call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)
                call masks_cano%write_netcdf(nf90_grpid)
                nf90_stat = nf90_def_grp(nf90_id, 'staggered', nf90_grpid)
                call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)
                call masks_stag%write_netcdf(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