mask_data_netcdf_s.f90 Source File


Contents


Source Code

submodule(mask_data_m) mask_data_netcdf_s
    !! Routines related with NetCDF I/O of mask_data
    use netcdf
    use error_handling_grillix_m, only: handle_error_netcdf, handle_error
    implicit none

contains
    
    module subroutine write_netcdf_mask_data(self, fgid)
        class(mask_data_t), intent(in) :: self
        integer, intent(in) :: fgid

        integer :: npoints, l
        integer :: nf90_stat, iddim_npoints, id_is_inner, id_is_boundary, id_is_ghost, &
            id_is_in_vessel, id_is_parbnd_immersed
        integer, allocatable, dimension(:) :: wrk_int
        
        ! Define dimensions
        npoints = size(self%is_inner)
        nf90_stat = nf90_def_dim(fgid, 'npoints', npoints , iddim_npoints)
        call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)

        ! Define variables
        nf90_stat = nf90_def_var(fgid, 'is_inner', NF90_INT, iddim_npoints, id_is_inner)
        call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)

        nf90_stat = nf90_def_var(fgid, 'is_boundary', NF90_INT, iddim_npoints, id_is_boundary)
        call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)

        nf90_stat = nf90_def_var(fgid, 'is_ghost', NF90_INT, iddim_npoints, id_is_ghost)
        call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)

        nf90_stat = nf90_def_var(fgid, 'is_in_vessel', NF90_INT, iddim_npoints, id_is_in_vessel)
        call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)

        nf90_stat = nf90_def_var(fgid, 'is_parbnd_immersed', NF90_INT, iddim_npoints, id_is_parbnd_immersed)
        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
        allocate(wrk_int(npoints))
        !$omp parallel do private(l) shared(npoints, self, wrk_int)
        do l = 1, npoints
            wrk_int(l) = merge(1, 0, self%is_inner(l))
        enddo
        !$omp end parallel do
        nf90_stat = nf90_put_var(fgid, id_is_inner, wrk_int)
        call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)

        !$omp parallel do private(l) shared(npoints, self, wrk_int)
        do l = 1, npoints
            wrk_int(l) = merge(1, 0, self%is_boundary(l))
        enddo
        !$omp end parallel do
        nf90_stat = nf90_put_var(fgid, id_is_boundary, wrk_int)
        call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)

        !$omp parallel do private(l) shared(npoints, self, wrk_int)
        do l = 1, npoints
            wrk_int(l) = merge(1, 0, self%is_ghost(l))
        enddo
        !$omp end parallel do
        nf90_stat = nf90_put_var(fgid, id_is_ghost, wrk_int)
        call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)

        !$omp parallel do private(l) shared(npoints, self, wrk_int)
        do l = 1, npoints
            wrk_int(l) = merge(1, 0, self%is_in_vessel(l))
        enddo
        !$omp end parallel do
        nf90_stat = nf90_put_var(fgid, id_is_in_vessel, wrk_int)
        call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)

        !$omp parallel do private(l) shared(npoints, self, wrk_int)
        do l = 1, npoints
            wrk_int(l) = merge(1, 0, self%is_parbnd_immersed(l))
        enddo
        !$omp end parallel do
        nf90_stat = nf90_put_var(fgid, id_is_parbnd_immersed, wrk_int)
        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__)

        deallocate(wrk_int)

    end subroutine
    
    module subroutine read_netcdf_mask_data(self, fgid)
        class(mask_data_t), intent(inout) :: self
        integer, intent(in) :: fgid

        integer :: npoints, l
        integer, allocatable, dimension(:) :: wrk_int
        integer :: nf90_stat, iddim_npoints, id_is_inner, id_is_boundary, id_is_ghost, &
            id_is_in_vessel, id_is_parbnd_immersed

        ! Get dimensions
        nf90_stat = nf90_inq_dimid(fgid, 'npoints', iddim_npoints)
        call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)
        nf90_stat =  nf90_inquire_dimension(fgid, iddim_npoints, len=npoints)
        call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)

        ! Allocate fields
        allocate(self%is_inner(npoints))
        allocate(self%is_boundary(npoints))
        allocate(self%is_ghost(npoints))
        allocate(self%is_in_vessel(npoints))
        allocate(self%is_parbnd_immersed(npoints))

        allocate(wrk_int(npoints))

        ! Read variables
        nf90_stat = nf90_inq_varid(fgid, 'is_inner', id_is_inner)
        call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)
        nf90_stat = nf90_get_var(fgid, id_is_inner, wrk_int)
        call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)
        !$omp parallel do private(l) shared(npoints, self, wrk_int)
        do l = 1, npoints
            self%is_inner(l) = (wrk_int(l) == 1)
        enddo
        !$omp end parallel do

        nf90_stat = nf90_inq_varid(fgid, 'is_boundary', id_is_boundary)
        call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)
        nf90_stat = nf90_get_var(fgid, id_is_boundary, wrk_int)
        call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)
        !$omp parallel do private(l) shared(npoints, self, wrk_int)
        do l = 1, npoints
            self%is_boundary(l) = (wrk_int(l) == 1)
        enddo
        !$omp end parallel do

        nf90_stat = nf90_inq_varid(fgid, 'is_ghost', id_is_ghost)
        call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)
        nf90_stat = nf90_get_var(fgid, id_is_ghost, wrk_int)
        call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)
        !$omp parallel do private(l) shared(npoints, self, wrk_int)
        do l = 1, npoints
            self%is_ghost(l) = (wrk_int(l) == 1)
        enddo
        !$omp end parallel do

        nf90_stat = nf90_inq_varid(fgid, 'is_in_vessel', id_is_in_vessel)
        call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)
        nf90_stat = nf90_get_var(fgid, id_is_in_vessel, wrk_int)
        call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)
        !$omp parallel do private(l) shared(npoints, self, wrk_int)
        do l = 1, npoints
            self%is_in_vessel(l) = (wrk_int(l) == 1)
        enddo
        !$omp end parallel do

        nf90_stat = nf90_inq_varid(fgid, 'is_parbnd_immersed', id_is_parbnd_immersed)
        call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)
        nf90_stat = nf90_get_var(fgid, id_is_parbnd_immersed, wrk_int)
        call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)
        !$omp parallel do private(l) shared(npoints, self, wrk_int)
        do l = 1, npoints
            self%is_parbnd_immersed(l) = (wrk_int(l) == 1)
        enddo
        !$omp end parallel do

    end subroutine

end submodule