polygon_in_mesh_netcdf_s.f90 Source File


Contents


Source Code

submodule (polygon_in_mesh_m) polygon_in_mesh_netcdf_s
    !! Routines related with I/O of polygon_in_mesh
    use netcdf
    use error_handling_grillix_m, only: handle_error_netcdf
    implicit none

contains

    module subroutine write_netcdf_polygon_in_mesh(self, fgid)
        class(polygon_in_mesh_t), intent(in) :: self
        integer, intent(in) :: fgid

         integer :: nf90_stat, iddim_np, iddim_nseg_pone, id_inds, id_ki_seg
         integer :: is_closed_int, is_cc_wise_int, warning_self_intersection_int

        ! define dimensions
        nf90_stat = nf90_def_dim(fgid, 'np', self%np , iddim_np)
        call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)

        nf90_stat = nf90_def_dim(fgid, 'nseg_plus_one', self%nsegs+1 , iddim_nseg_pone)
        call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)

        ! define variables
        nf90_stat = nf90_def_var(fgid, 'inds', NF90_INT, iddim_np, id_inds)
        call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)
        
        nf90_stat = nf90_def_var(fgid, 'ki_seg', NF90_INT, iddim_nseg_pone, id_ki_seg)
        call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)

        ! write global attributes (metadata)
        if (self%is_closed) then
            is_closed_int = 1
        else
            is_closed_int = 0
        endif
        nf90_stat = nf90_put_att(fgid, NF90_GLOBAL, 'is_closed', is_closed_int)
        call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)

        if (self%is_cc_wise) then
            is_cc_wise_int = 1
        else
            is_cc_wise_int = 0
        endif
        nf90_stat = nf90_put_att(fgid, NF90_GLOBAL, 'is_cc_wise', is_cc_wise_int)
        call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)

        if (self%warning_self_intersection) then
            warning_self_intersection_int = 1
        else
            warning_self_intersection_int = 0
        endif
        nf90_stat = nf90_put_att(fgid, NF90_GLOBAL, 'warning_self_intersection', &
                                 warning_self_intersection_int)
        call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)

        ! end of definition
        nf90_stat = nf90_enddef(fgid)
        call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)

        ! put variables
        nf90_stat = nf90_put_var(fgid, id_inds, self%inds)
        call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)

        nf90_stat = nf90_put_var(fgid, id_ki_seg, self%ki_seg)
        call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)

        ! put into redefinition state
        nf90_stat = nf90_redef(fgid)
        call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)

    end subroutine


    module subroutine read_netcdf_polygon_in_mesh(self, fgid)
        class(polygon_in_mesh_t), intent(inout) :: self
        integer, intent(in) :: fgid
    
        integer :: nf90_stat, iddim_np, iddim_nseg_pone, id_inds, id_ki_seg
        integer :: is_closed_int, is_cc_wise_int, warning_self_intersection_int, nsegs_pone

        ! Make sure that properties of instance are deallocated
        call destructor(self)

        ! read global attributes
        nf90_stat = nf90_get_att(fgid, NF90_GLOBAL, 'is_closed', is_closed_int)
        call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)
        if (is_closed_int == 1) then
            self%is_closed = .true.
        else
            self%is_closed = .false.
        endif

        nf90_stat = nf90_get_att(fgid, NF90_GLOBAL, 'is_cc_wise', is_cc_wise_int)
        call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)
        if (is_cc_wise_int == 1) then
            self%is_cc_wise = .true.
        else
            self%is_cc_wise = .false.
        endif

        nf90_stat = nf90_get_att(fgid, NF90_GLOBAL, 'warning_self_intersection', &
                                 warning_self_intersection_int)
        call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)
        if (warning_self_intersection_int== 1) then
            self%warning_self_intersection = .true.
        else
            self%warning_self_intersection = .false.
        endif

        ! get dimensions
        nf90_stat = nf90_inq_dimid(fgid, 'np', iddim_np)
        call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)
        nf90_stat = nf90_inquire_dimension(fgid, iddim_np, len=self%np)
        call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)

        nf90_stat = nf90_inq_dimid(fgid, 'nseg_plus_one', iddim_nseg_pone)
        call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)
        nf90_stat = nf90_inquire_dimension(fgid, iddim_nseg_pone, len=nsegs_pone)
        call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)
        self%nsegs = nsegs_pone - 1

        ! allocate fields
        allocate(self%inds(self%np))
        allocate(self%ki_seg(self%nsegs+1))

         ! read fields
        nf90_stat = nf90_inq_varid(fgid, 'inds', id_inds)
        call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)
        nf90_stat = nf90_get_var(fgid, id_inds, self%inds)
        call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)

        nf90_stat = nf90_inq_varid(fgid, 'ki_seg', id_ki_seg)
        call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)
        nf90_stat = nf90_get_var(fgid, id_ki_seg, self%ki_seg)
        call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)
    
    end subroutine

end submodule