test_parallel_map3d_m.pf Source File


Contents


Source Code

module test_parallel_map3d_m
    !! Tests for the map routines in 3d stellarator geometries
    use pfunit
    use helpers_m, only : almost_equal, digit_sum_csr
    use grillix_build_info_m, only : git_repo_path
    use MPI
    use comm_handler_m, only : comm_handler_t
    use precision_grillix_m, only : GP
    use constants_m, only : TWO_PI
    use csrmat_m, only : csrmat_t
    use equilibrium_m, only : equilibrium_t
    use equilibrium_factory_m, only : create_equilibrium, DOMMASCHK
    use mesh_cart_m, only : mesh_cart_t
    use mpi_mapping_auxiliaries_m, only : getdata_fwdbwdplane
    use parallel_map_m, only : parallel_map_t
    use netcdf
    implicit none
    
    real(GP), parameter :: trace_prec = 10E-7_GP
    
    ! Reference data
    real(GP), dimension(4), parameter :: ref_dpar_fwd_half = &
        [8.109940915299E-01_GP, 8.044066191614E-01_GP, 8.298382641725E-01_GP, 8.371604916492E-01_GP]
    real(GP), dimension(4), parameter :: ref_dpar_bwd_half = &
        [8.097832649010E-01_GP, 8.202903105972E-01_GP, 8.471039311781E-01_GP, 8.355766086507E-01_GP]
    real(GP), dimension(4), parameter :: ref_dpar_fwd_full = &
        [1.596287027802E+00_GP, 1.611192453990E+00_GP, 1.647061165634E+00_GP, 1.629495083576E+00_GP]
    real(GP), dimension(4), parameter :: ref_dpar_bwd_full = &
        [1.646852278182E+00_GP, 1.633834783669E+00_GP, 1.675962788927E+00_GP, 1.693809193333E+00_GP]
    real(GP), dimension(4), parameter :: ref_fluxbox_vol = &
        [1.604679369248E-04_GP, 1.602771975793E-04_GP, 1.725651938550E-04_GP, 1.726492995507E-04_GP]
    
    real(GP), dimension(4), parameter :: ref_dpar_fwd_half_stag = &
        [8.391447748292E-01_GP, 8.252207747489E-01_GP, 8.026059378764E-01_GP, 8.154286539574E-01_GP]
    real(GP), dimension(4), parameter :: ref_dpar_bwd_half_stag = &
        [8.471918318741E-01_GP, 8.200521764416E-01_GP, 8.099458741252E-01_GP, 8.355591111122E-01_GP]
    real(GP), dimension(4), parameter :: ref_dpar_fwd_full_stag = &
        [1.645662688617E+00_GP, 1.609375378967E+00_GP, 1.596802539322E+00_GP, 1.632191638072E+00_GP]
    real(GP), dimension(4), parameter :: ref_dpar_bwd_full_stag = &
        [1.694765877082E+00_GP, 1.672732839240E+00_GP, 1.631689585973E+00_GP, 1.651271938455E+00_GP]
    real(GP), dimension(4), parameter :: ref_fluxbox_vol_stag = &
        [1.753256596058E-04_GP, 1.664289926310E-04_GP, 1.580065973086E-04_GP, 1.662010377293E-04_GP]

    integer, dimension(4), parameter :: ref_map_stag_cano_fwd_nnz = &
        [3032, 3008, 3029, 3011]    
    integer, dimension(4), parameter :: ref_map_stag_cano_fwd_j = &
        [257,  223,  161,  194]
    real(GP), dimension(4), parameter :: ref_map_stag_cano_fwd_val = &
        [3.374995633528E-04_GP, 1.222234982036E-03_GP, 1.598124000418E-03_GP, 1.436909739511E-03_GP]
    real(GP), dimension(4), parameter :: ref_map_stag_cano_fwd_dsum = &
        [2.690000000000E+02_GP, 2.690000000000E+02_GP, 2.750000000000E+02_GP, 2.630000000000E+02_GP]
    integer, dimension(4), parameter :: ref_map_stag_cano_bwd_nnz = &
        [3011, 3029, 3008, 3032]    
    integer, dimension(4), parameter :: ref_map_stag_cano_bwd_j = &
        [184,  253,  144,  172]
    real(GP), dimension(4), parameter :: ref_map_stag_cano_bwd_val = &
        [9.730549633343E-02_GP, 8.423970861761E-01_GP, 4.096881745275E-03_GP, 3.426701368566E-03_GP]
    real(GP), dimension(4), parameter :: ref_map_stag_cano_bwd_dsum = &
        [2.630000000000E+02_GP, 2.750000000000E+02_GP, 2.690000000000E+02_GP, 2.690000000000E+02_GP]

    integer, dimension(4), parameter :: ref_map_cano_stag_fwd_nnz = &
        [3061, 3108, 3061, 2981]    
    integer, dimension(4), parameter :: ref_map_cano_stag_fwd_j = &
        [210,  143,  245,  264]
    real(GP), dimension(4), parameter :: ref_map_cano_stag_fwd_val = &
        [3.460249705519E-03_GP, 1.998407315066E-03_GP, 4.888856946138E-04_GP, 1.613797275680E-01_GP]
    real(GP), dimension(4), parameter :: ref_map_cano_stag_fwd_dsum = &
        [2.710000000000E+02_GP, 2.700000000000E+02_GP, 2.710000000000E+02_GP, 2.630000000000E+02_GP]
    integer, dimension(4), parameter :: ref_map_cano_stag_bwd_nnz = &
        [3061, 2981, 3061, 3108]    
    integer, dimension(4), parameter :: ref_map_cano_stag_bwd_j = &
        [239,  125,  173,  249]
    real(GP), dimension(4), parameter :: ref_map_cano_stag_bwd_val = &
        [1.366779434239E-03_GP, 1.492377500038E-03_GP, 4.108573911257E-03_GP, 2.134135175149E-01_GP]
    real(GP), dimension(4), parameter :: ref_map_cano_stag_bwd_dsum = &
        [2.710000000000E+02_GP, 2.630000000000E+02_GP, 2.710000000000E+02_GP, 2.700000000000E+02_GP]    
        
    integer, dimension(4), parameter :: ref_grad_stag_cano_fwd_nnz = &
        [4521, 4442, 4567, 4460]    
    integer, dimension(4), parameter :: ref_grad_stag_cano_fwd_j = &
        [ 257,  215,  160,  193]
    real(GP), dimension(4), parameter :: ref_grad_stag_cano_fwd_val = &
        [3.982309700978E-04_GP, 5.728191225148E-04_GP, 6.346928619128E-05_GP, 5.296650445312E-04_GP]
    real(GP), dimension(4), parameter :: ref_grad_stag_cano_fwd_dsum = &
        [1.716142759309E+02_GP, 1.729354151250E+02_GP, 1.753824705662E+02_GP, 1.674894722847E+02_GP]
    integer, dimension(4), parameter :: ref_grad_stag_cano_bwd_nnz = &
        [4460, 4567, 4442, 4521]    
    integer, dimension(4), parameter :: ref_grad_stag_cano_bwd_j = &
        [174,  239,  136,  166]
    real(GP), dimension(4), parameter :: ref_grad_stag_cano_bwd_val = &
        [8.796162159344E-05_GP, 4.580986364616E-04_GP, 3.281858157403E-04_GP, 2.749827551295E-04_GP]
    real(GP), dimension(4), parameter :: ref_grad_stag_cano_bwd_dsum = &
        [1.674894727016E+02_GP, 1.753824710060E+02_GP, 1.729354147046E+02_GP, 1.716142754954E+02_GP]
       
    integer, dimension(4), parameter :: ref_pdiv_cano_stag_fwd_nnz = &
        [4460, 4567, 4442, 4521]    
    integer, dimension(4), parameter :: ref_pdiv_cano_stag_fwd_j = &
        [156,  151,  199,  196]
    real(GP), dimension(4), parameter :: ref_pdiv_cano_stag_fwd_val = &
        [-5.845057707572E-07_GP, 1.316895418855E-04_GP, -1.827118536358E-04_GP, -9.238020016211E-05_GP]
    real(GP), dimension(4), parameter :: ref_pdiv_cano_stag_fwd_dsum = &
        [-1.674437972331E+02_GP, -1.766337369756E+02_GP,  -1.744762490743E+02_GP, -1.727737730022E+02_GP]   
    integer, dimension(4), parameter :: ref_pdiv_cano_stag_bwd_nnz = &
        [4460, 4521, 4442, 4567]    
    integer, dimension(4), parameter :: ref_pdiv_cano_stag_bwd_j = &
        [ 142, 54,  176,  200]
    real(GP), dimension(4), parameter :: ref_pdiv_cano_stag_bwd_val = &
        [3.034447258607E-05_GP, -3.201678639284E-04_GP, -4.445158439019E-02_GP, -6.503608420823E-04_GP]
    real(GP), dimension(4), parameter :: ref_pdiv_cano_stag_bwd_dsum = &
        [-1.674437971647E+02_GP, -1.727737720537E+02_GP, -1.744762482897E+02_GP, -1.766337370860E+02_GP]
    
    integer, dimension(4), parameter :: ref_map_cano_cano_fwd_nnz = &
        [3045, 3017, 3085, 3040]        
    integer, dimension(4), parameter :: ref_map_cano_cano_fwd_j = &
        [202,  182,  223,  230]
    real(GP), dimension(4), parameter :: ref_map_cano_cano_fwd_val = &
        [1.407753994887E-03_GP, 2.714439376066E-06_GP, 1.418522821409E-03_GP, 3.744379122314E-01_GP]
    real(GP), dimension(4), parameter :: ref_map_cano_cano_fwd_dsum = &
        [2.730000000000E+02_GP, 2.750000000000E+02_GP, 2.710000000000E+02_GP, 2.770000000000E+02_GP] 
    integer, dimension(4), parameter :: ref_map_cano_cano_bwd_nnz = &
        [3045, 3040, 3085, 3017]      
    integer, dimension(4), parameter :: ref_map_cano_cano_bwd_j = &
        [178,  169,  102,  158]
    real(GP), dimension(4), parameter :: ref_map_cano_cano_bwd_val = &
        [1.124360392877E-03_GP, 9.364405067968E-04_GP, 1.358226668233E-03_GP, 4.304396154699E-01_GP]
    real(GP), dimension(4), parameter :: ref_map_cano_cano_bwd_dsum = &
        [2.730000000000E+02_GP, 2.770000000000E+02_GP, 2.710000000000E+02_GP, 2.750000000000E+02_GP]
    
    integer, dimension(4), parameter :: ref_map_stag_stag_fwd_nnz = &
        [3126, 3088, 2990, 3017]      
    integer, dimension(4), parameter :: ref_map_stag_stag_fwd_j = &
        [235,  218,  193,  203]
    real(GP), dimension(4), parameter :: ref_map_stag_stag_fwd_val = &
        [6.474938589139E-04_GP, 3.485006669940E-03_GP, 9.738190800335E-04_GP, 2.998704589256E-03_GP]
    real(GP), dimension(4), parameter :: ref_map_stag_stag_fwd_dsum = &
        [2.700000000000E+02_GP, 2.800000000000E+02_GP, 2.750000000000E+02_GP, 2.690000000000E+02_GP]
    integer, dimension(4), parameter :: ref_map_stag_stag_bwd_nnz = &
        [3017, 2990, 3088, 3126]   
    integer, dimension(4), parameter :: ref_map_stag_stag_bwd_j = &
        [94,  169,  185,  161]
    real(GP), dimension(4), parameter :: ref_map_stag_stag_bwd_val = &
        [2.717561639101E-03_GP, 4.551374492498E-01_GP, 2.022067797877E-04_GP, 1.668248076998E-03_GP]
    real(GP), dimension(4), parameter :: ref_map_stag_stag_bwd_dsum = &
        [2.690000000000E+02_GP, 2.750000000000E+02_GP, 2.800000000000E+02_GP, 2.700000000000E+02_GP]

    integer, dimension(4), parameter :: ref_map_cano_cano_fwd2_nnz = &
        [2990, 3038, 3112, 3116]
    integer, dimension(4), parameter :: ref_map_cano_cano_fwd2_j = &
        [187,  212,  209,  207]
    real(GP), dimension(4), parameter :: ref_map_cano_cano_fwd2_val = &
        [1.388401141315758E-003_GP, 3.237923035015976E-003_GP, 2.109140503404123E-003_GP, 8.092166928255919E-002_GP]
    real(GP), dimension(4), parameter :: ref_map_cano_cano_fwd2_dsum = &
        [269.0_GP, 272.0_GP, 271.0_GP, 278.0_GP]
    integer, dimension(4), parameter :: ref_map_cano_cano_bwd2_nnz = &
        [2990, 3116, 3112, 3038]
    integer, dimension(4), parameter :: ref_map_cano_cano_bwd2_j = &
        [41,  22,  8,  27]
    real(GP), dimension(4), parameter :: ref_map_cano_cano_bwd2_val = &
        [2.022131891398691E-003_GP, 3.394727559428792E-004_GP, 2.387225561809455E-003_GP, 3.119472216713043E-002_GP]
    real(GP), dimension(4), parameter :: ref_map_cano_cano_bwd2_dsum = &
        [269.0_GP, 278.0_GP, 271.0_GP, 272.0_GP]

contains

    @test(npes = [4])
    subroutine parallel_map_dommaschk(this)
        class (MpiTestMethod), intent(inout) :: this
        integer :: comm_world
        type(comm_handler_t) :: comm_handler
        integer :: rank, ierr
        class(equilibrium_t), allocatable :: equi
        character(len=:), allocatable :: equi_params_file
        type(mesh_cart_t) :: mesh_cano, mesh_stag
        type(parallel_map_t), target  :: map
        real(GP) :: phi_cano, phi_stag, dphi, spacing_f, rret
        integer :: ic, jc, iz, iret, np, np_stag, nrecv, isum, jsum
        integer :: l, ls, lpnt
        integer, allocatable, dimension(:) :: np_recv
        integer, dimension(-4:4) :: nps
        logical :: fine
        integer :: imat
        type(csrmat_t), pointer :: test_mat => null()
        integer :: ref_ndim, ref_ncol, ref_nnz, ref_j
        real(GP) :: ref_val, ref_dsum
                     
        comm_world = this%getMpiCommunicator()
        call MPI_comm_rank(comm_world, rank, ierr)

        if (rank == 0) then     
            write(*,*)''
            write(*,*)''
            write(*,'(A80)')'test_parallel_map_dommaschk -----------------------------------------------'
        endif
        
        call comm_handler%init(comm_world, 4, 1)
        dphi = TWO_PI / comm_handler%get_nplanes()
        phi_cano = dphi * comm_handler%get_plane()
        phi_stag = phi_cano + dphi / 2.0_GP

        ! Create equilibrium
        equi_params_file = git_repo_path // &
                           '/cicd/res/params_equi_dommaschk.in'
        call create_equilibrium(equi, DOMMASCHK, equi_params_file)

        ! Create mesh
        spacing_f = 1E-2
        call mesh_cano%init(equi, phi_cano, lvl=1, spacing_f=spacing_f, &
                            size_neighbor=2, size_ghost_layer=2, &
                            extend_beyond_wall=.false., dbgout=0)
        call mesh_cano%reorder(grid_reorder_block=8, dbgout=0)
        call mesh_stag%init(equi, phi_stag, lvl=1, spacing_f=spacing_f, &
                            size_neighbor=2, size_ghost_layer=2,  &
                            extend_beyond_wall=.false., dbgout=0)
        call mesh_stag%reorder(grid_reorder_block=8, dbgout=0)
        
        ! Number of points of adjacent meshes and staggered meshes
        np = mesh_cano%get_n_points()
        np_stag = mesh_stag%get_n_points()
        call getdata_fwdbwdplane(comm_world, -2, 1, [np], nrecv, np_recv)
        nps(-4) = np_recv(1)
        deallocate(np_recv)
        call getdata_fwdbwdplane(comm_world, -2, 1, [np_stag], nrecv, np_recv)
        nps(-3) = np_recv(1)
        deallocate(np_recv)
        call getdata_fwdbwdplane(comm_world, -1, 1, [np], nrecv, np_recv)
        nps(-2) = np_recv(1)
        deallocate(np_recv)
        call getdata_fwdbwdplane(comm_world, -1, 1, [np_stag], nrecv, np_recv)
        nps(-1) = np_recv(1)
        deallocate(np_recv)
        call getdata_fwdbwdplane(comm_world, 0, 1, [np], nrecv, np_recv)
        nps(0) = np_recv(1)
        deallocate(np_recv)
        call getdata_fwdbwdplane(comm_world, 0, 1, [np_stag], nrecv, np_recv)
        nps(1) = np_recv(1)
        deallocate(np_recv)
        call getdata_fwdbwdplane(comm_world, 1, 1, [np], nrecv, np_recv)
        nps(2) = np_recv(1)
        deallocate(np_recv)
        call getdata_fwdbwdplane(comm_world, 1, 1, [np_stag], nrecv, np_recv)
        nps(3) = np_recv(1)
        deallocate(np_recv)
        call getdata_fwdbwdplane(comm_world, 2, 1, [np], nrecv, np_recv)
        nps(4) = np_recv(1)
        deallocate(np_recv)

        ! Build map
        call map%set_parameters(comm_handler, intorder_in=3, xorder_in=1, &
                                use_fixed_stencil_in=.false., &
                                use_gauss_quadrature_in=.false.)
        call map%build(comm_handler, equi, mesh_cano, mesh_stag, dbgout=1)
                
        ! Test metadata
        iret = map%get_nplanes()
        @assertEqual(iret, 4)
        iret = map%get_iplane()
        @assertEqual(iret, rank)
        rret = map%get_phi_cano()
        fine = almost_equal(rret, phi_cano, rtol=0.0_GP, atol=1.0E-12_GP) 
        @assertTrue(fine)
        rret = map%get_phi_stag()
        fine = almost_equal(rret, phi_stag, rtol=0.0_GP, atol=1.0E-12_GP) 
        @assertTrue(fine)
        rret = map%get_dphi()
        fine = almost_equal(rret, dphi, rtol=0.0_GP, atol=1.0E-12_GP) 
        @assertTrue(fine)
        iret = map%get_intorder()
        @assertEqual(iret, 3)
        iret = map%get_xorder()
        @assertEqual(iret, 1)
        
        ! Select test points
        ic = 25
        jc = 24
        l = mesh_cano%cart_to_mesh_index(ic, jc)
        ls = mesh_stag%cart_to_mesh_index(ic, jc)
        
        ! Check trace data (mesh_cano)
        rret = map%dpar_cano_stag_fwd(l)
        fine = almost_equal(rret, ref_dpar_fwd_half(rank+1), &
                            rtol=trace_prec, atol=0.0_GP) 
        @assertTrue(fine)
                
        rret = map%dpar_cano_stag_bwd(l)
        fine = almost_equal(rret, ref_dpar_bwd_half(rank+1), &
                            rtol=trace_prec, atol=0.0_GP) 
        @assertTrue(fine)
        
        rret = map%dpar_cano_cano_fwd(l)
        fine = almost_equal(rret, ref_dpar_fwd_full(rank+1), &
                            rtol=trace_prec, atol=0.0_GP) 
        @assertTrue(fine)
        
        rret = map%dpar_cano_cano_bwd(l)
        fine = almost_equal(rret, ref_dpar_bwd_full(rank+1), &
                            rtol=trace_prec, atol=0.0_GP) 
        @assertTrue(fine)
        
        rret = map%fluxbox_vol_cano(l)
        fine = almost_equal(rret, ref_fluxbox_vol(rank+1), &
                            rtol=trace_prec, atol=0.0_GP) 
        @assertTrue(fine)
                
        ! Check trace data (mesh_staggered)
        rret = map%dpar_stag_cano_fwd(ls)
        fine = almost_equal(rret, ref_dpar_fwd_half_stag(rank+1), &
                            rtol=trace_prec, atol=0.0_GP) 
        @assertTrue(fine)       
        
        rret = map%dpar_stag_cano_bwd(ls)
        fine = almost_equal(rret, ref_dpar_bwd_half_stag(rank+1), &
                            rtol=trace_prec, atol=0.0_GP) 
        @assertTrue(fine)
        
        rret = map%dpar_stag_stag_fwd(ls)
        fine = almost_equal(rret, ref_dpar_fwd_full_stag(rank+1), &
                            rtol=trace_prec, atol=0.0_GP) 
        @assertTrue(fine)
        
        rret = map%dpar_stag_stag_bwd(ls)
        fine = almost_equal(rret, ref_dpar_bwd_full_stag(rank+1), &
                            rtol=trace_prec, atol=0.0_GP) 
        @assertTrue(fine)
              
        rret = map%fluxbox_vol_stag(ls)
        fine = almost_equal(rret, ref_fluxbox_vol_stag(rank+1), &
                            rtol=trace_prec, atol=0.0_GP) 
        @assertTrue(fine)
                     
        ! Check matrices: Dimensions, specific entry and digit sum of values
        do imat = 1, 14
            select case(imat)
                case(1)
                    lpnt = ls
                    test_mat => map%map_stag_cano_fwd
                    ref_ndim = nps(1)
                    ref_ncol = nps(2)
                    ref_nnz = ref_map_stag_cano_fwd_nnz(rank+1)
                    ref_j = ref_map_stag_cano_fwd_j(rank+1)
                    ref_val =  ref_map_stag_cano_fwd_val(rank+1)
                    ref_dsum = 1.0_GP * ref_map_stag_cano_fwd_dsum(rank+1)
                case(2)
                    lpnt = ls
                    test_mat => map%map_stag_cano_bwd
                    ref_ndim = nps(1)
                    ref_ncol = nps(0)
                    ref_nnz = ref_map_stag_cano_bwd_nnz(rank+1)
                    ref_j = ref_map_stag_cano_bwd_j(rank+1)
                    ref_val =  ref_map_stag_cano_bwd_val(rank+1)
                    ref_dsum = ref_map_stag_cano_bwd_dsum(rank+1)
                case(3)
                    lpnt = l
                    test_mat => map%map_cano_stag_fwd
                    ref_ndim = nps(0)
                    ref_ncol = nps(1)
                    ref_nnz = ref_map_cano_stag_fwd_nnz(rank+1)
                    ref_j = ref_map_cano_stag_fwd_j(rank+1)
                    ref_val =  ref_map_cano_stag_fwd_val(rank+1)
                    ref_dsum = ref_map_cano_stag_fwd_dsum(rank+1)
                case(4)
                    lpnt = l
                    test_mat => map%map_cano_stag_bwd
                    ref_ndim = nps(0)
                    ref_ncol = nps(-1)
                    ref_nnz = ref_map_cano_stag_bwd_nnz(rank+1)
                    ref_j = ref_map_cano_stag_bwd_j(rank+1)
                    ref_val =  ref_map_cano_stag_bwd_val(rank+1)
                    ref_dsum = ref_map_cano_stag_bwd_dsum(rank+1)
                case(5)
                    lpnt = ls
                    test_mat => map%grad_stag_cano_fwd
                    ref_ndim = nps(1)
                    ref_ncol = nps(2)
                    ref_nnz = ref_grad_stag_cano_fwd_nnz(rank+1)
                    ref_j = ref_grad_stag_cano_fwd_j(rank+1)
                    ref_val =  ref_grad_stag_cano_fwd_val(rank+1)
                    ref_dsum = ref_grad_stag_cano_fwd_dsum(rank+1)
                case(6)
                    lpnt = ls
                    test_mat => map%grad_stag_cano_bwd
                    ref_ndim = nps(1)
                    ref_ncol = nps(0)
                    ref_nnz = ref_grad_stag_cano_bwd_nnz(rank+1)
                    ref_j = ref_grad_stag_cano_bwd_j(rank+1)
                    ref_val =  ref_grad_stag_cano_bwd_val(rank+1)
                    ref_dsum = ref_grad_stag_cano_bwd_dsum(rank+1)
                case(7)
                    lpnt = l
                    test_mat => map%pdiv_cano_stag_fwd
                    ref_ndim = nps(0)
                    ref_ncol = nps(1)
                    ref_nnz = ref_pdiv_cano_stag_fwd_nnz(rank+1)
                    ref_j = ref_pdiv_cano_stag_fwd_j(rank+1)
                    ref_val =  ref_pdiv_cano_stag_fwd_val(rank+1)
                    ref_dsum = ref_pdiv_cano_stag_fwd_dsum(rank+1)
                case(8)
                    lpnt = l
                    test_mat => map%pdiv_cano_stag_bwd
                    ref_ndim = nps(0)
                    ref_ncol = nps(-1)
                    ref_nnz = ref_pdiv_cano_stag_bwd_nnz(rank+1)
                    ref_j = ref_pdiv_cano_stag_bwd_j(rank+1)
                    ref_val =  ref_pdiv_cano_stag_bwd_val(rank+1)
                    ref_dsum = ref_pdiv_cano_stag_bwd_dsum(rank+1)
                case(9)
                    lpnt = l
                    test_mat => map%map_cano_cano_fwd
                    ref_ndim = nps(0)
                    ref_ncol = nps(2)
                    ref_nnz = ref_map_cano_cano_fwd_nnz(rank+1)
                    ref_j = ref_map_cano_cano_fwd_j(rank+1)
                    ref_val =  ref_map_cano_cano_fwd_val(rank+1)
                    ref_dsum = ref_map_cano_cano_fwd_dsum(rank+1)
                case(10)
                    lpnt = l
                    test_mat => map%map_cano_cano_bwd
                    ref_ndim = nps(0)
                    ref_ncol = nps(-2)
                    ref_nnz = ref_map_cano_cano_bwd_nnz(rank+1)
                    ref_j = ref_map_cano_cano_bwd_j(rank+1)
                    ref_val =  ref_map_cano_cano_bwd_val(rank+1)
                    ref_dsum = ref_map_cano_cano_bwd_dsum(rank+1)
                case(11)
                    lpnt = ls
                    test_mat => map%map_stag_stag_fwd
                    ref_ndim = nps(1)
                    ref_ncol = nps(3)
                    ref_nnz = ref_map_stag_stag_fwd_nnz(rank+1)
                    ref_j = ref_map_stag_stag_fwd_j(rank+1)
                    ref_val =  ref_map_stag_stag_fwd_val(rank+1)
                    ref_dsum = ref_map_stag_stag_fwd_dsum(rank+1)
                case(12)
                    lpnt = ls
                    test_mat => map%map_stag_stag_bwd
                    ref_ndim = nps(1)
                    ref_ncol = nps(-1)
                    ref_nnz = ref_map_stag_stag_bwd_nnz(rank+1)
                    ref_j = ref_map_stag_stag_bwd_j(rank+1)
                    ref_val =  ref_map_stag_stag_bwd_val(rank+1)
                    ref_dsum = ref_map_stag_stag_bwd_dsum(rank+1)
                case(13)
                    lpnt = l
                    test_mat => map%map_cano_cano_fwd2
                    ref_ndim = nps(0)
                    ref_ncol = nps(4)
                    ref_nnz = ref_map_cano_cano_fwd2_nnz(rank+1)
                    ref_j = ref_map_cano_cano_fwd2_j(rank+1)
                    ref_val =  ref_map_cano_cano_fwd2_val(rank+1)
                    ref_dsum = ref_map_cano_cano_fwd2_dsum(rank+1)
                case(14)
                    lpnt = l
                    test_mat => map%map_cano_cano_bwd2
                    ref_ndim = nps(0)
                    ref_ncol = nps(-4)
                    ref_nnz = ref_map_cano_cano_bwd2_nnz(rank+1)
                    ref_j = ref_map_cano_cano_bwd2_j(rank+1)
                    ref_val =  ref_map_cano_cano_bwd2_val(rank+1)
                    ref_dsum = ref_map_cano_cano_bwd2_dsum(rank+1)
                case default
                    @assertTrue(.false.)                
            end select
            
            iz = test_mat%i(lpnt)
            @assertEqual(test_mat%ndim, ref_ndim)
            @assertEqual(test_mat%ncol, ref_ncol)
            @assertEqual(test_mat%nnz, ref_nnz)
            @assertEqual(test_mat%j(iz), ref_j)

            fine = almost_equal(test_mat%val(iz), ref_val, rtol=0.0_GP, atol=1.0E-10_GP)
            @assertTrue(fine)

            call digit_sum_csr(test_mat, isum, jsum, rret)
            fine = almost_equal(rret, ref_dsum, rtol=0.0_GP, atol=1.0E-9_GP)
            @assertTrue(fine)
                        
        enddo
               
        if (rank == 0) then     
            write(*,'(A80)')'test_parallel_map_dommaschk complete --------------------------------------'
            write(*,*)''
            write(*,*)''
        endif  

    end subroutine
    
end module