From 39d6bf01e0ac483d21c61c89a021c492db482be6 Mon Sep 17 00:00:00 2001 From: Salvatore Filippone Date: Thu, 23 Jul 2020 12:00:44 +0200 Subject: [PATCH] Reworked interface METISPART: use LPK matrix input. --- test/fileread/psb_cf_sample.f90 | 11 +-- test/fileread/psb_df_sample.f90 | 11 +-- test/fileread/psb_sf_sample.f90 | 11 +-- test/fileread/psb_zf_sample.f90 | 11 +-- util/Makefile | 2 +- util/psb_metis_int.c | 35 ++-------- util/psb_metispart_mod.F90 | 118 ++++++++++++++++---------------- util/psi_build_mtpart.F90 | 100 +++++++++++++++++++++------ 8 files changed, 174 insertions(+), 125 deletions(-) diff --git a/test/fileread/psb_cf_sample.f90 b/test/fileread/psb_cf_sample.f90 index 8170a446..699c4347 100644 --- a/test/fileread/psb_cf_sample.f90 +++ b/test/fileread/psb_cf_sample.f90 @@ -38,10 +38,11 @@ program psb_cf_sample implicit none ! input parameters - character(len=40) :: kmethd, ptype, mtrx_file, rhs_file,renum + character(len=40) :: kmethd, ptype, mtrx_file, rhs_file ! sparse matrices - type(psb_cspmat_type) :: a, aux_a + type(psb_cspmat_type) :: a + type(psb_lcspmat_type) :: aux_a ! preconditioner data type(psb_cprec_type) :: prec @@ -56,6 +57,7 @@ program psb_cf_sample type(psb_desc_type):: desc_a integer(psb_ipk_) :: ictxt, iam, np + integer(psb_lpk_) :: lnp ! solver paramters integer(psb_ipk_) :: iter, itmax, ierr, itrace, ircode,& @@ -140,7 +142,6 @@ program psb_cf_sample m_problem = aux_a%get_nrows() call psb_bcast(ictxt,m_problem) - call psb_mat_renum(psb_mat_renum_identity_,aux_a,info,perm) ! At this point aux_b may still be unallocated if (size(aux_b,dim=1) == m_problem) then @@ -180,7 +181,9 @@ program psb_cf_sample write(psb_out_unit,'("Partition type: graph vector")') write(psb_out_unit,'(" ")') ! write(psb_err_unit,'("Build type: graph")') - call build_mtpart(aux_a,np) + call aux_a%cscnv(info,type='csr') + lnp = np + call build_mtpart(aux_a,lnp) endif call psb_barrier(ictxt) diff --git a/test/fileread/psb_df_sample.f90 b/test/fileread/psb_df_sample.f90 index 17f821b3..8fbf8bc4 100644 --- a/test/fileread/psb_df_sample.f90 +++ b/test/fileread/psb_df_sample.f90 @@ -38,10 +38,11 @@ program psb_df_sample implicit none ! input parameters - character(len=40) :: kmethd, ptype, mtrx_file, rhs_file,renum + character(len=40) :: kmethd, ptype, mtrx_file, rhs_file ! sparse matrices - type(psb_dspmat_type) :: a, aux_a + type(psb_dspmat_type) :: a + type(psb_ldspmat_type) :: aux_a ! preconditioner data type(psb_dprec_type) :: prec @@ -56,6 +57,7 @@ program psb_df_sample type(psb_desc_type):: desc_a integer(psb_ipk_) :: ictxt, iam, np + integer(psb_lpk_) :: lnp ! solver paramters integer(psb_ipk_) :: iter, itmax, ierr, itrace, ircode,& @@ -140,7 +142,6 @@ program psb_df_sample m_problem = aux_a%get_nrows() call psb_bcast(ictxt,m_problem) - call psb_mat_renum(psb_mat_renum_identity_,aux_a,info,perm) ! At this point aux_b may still be unallocated if (size(aux_b,dim=1) == m_problem) then @@ -180,7 +181,9 @@ program psb_df_sample write(psb_out_unit,'("Partition type: graph vector")') write(psb_out_unit,'(" ")') ! write(psb_err_unit,'("Build type: graph")') - call build_mtpart(aux_a,np) + call aux_a%cscnv(info,type='csr') + lnp = np + call build_mtpart(aux_a,lnp) endif call psb_barrier(ictxt) diff --git a/test/fileread/psb_sf_sample.f90 b/test/fileread/psb_sf_sample.f90 index 05127088..71dd646b 100644 --- a/test/fileread/psb_sf_sample.f90 +++ b/test/fileread/psb_sf_sample.f90 @@ -38,10 +38,11 @@ program psb_sf_sample implicit none ! input parameters - character(len=40) :: kmethd, ptype, mtrx_file, rhs_file,renum + character(len=40) :: kmethd, ptype, mtrx_file, rhs_file ! sparse matrices - type(psb_sspmat_type) :: a, aux_a + type(psb_sspmat_type) :: a + type(psb_lsspmat_type) :: aux_a ! preconditioner data type(psb_sprec_type) :: prec @@ -56,6 +57,7 @@ program psb_sf_sample type(psb_desc_type):: desc_a integer(psb_ipk_) :: ictxt, iam, np + integer(psb_lpk_) :: lnp ! solver paramters integer(psb_ipk_) :: iter, itmax, ierr, itrace, ircode,& @@ -140,7 +142,6 @@ program psb_sf_sample m_problem = aux_a%get_nrows() call psb_bcast(ictxt,m_problem) - call psb_mat_renum(psb_mat_renum_identity_,aux_a,info,perm) ! At this point aux_b may still be unallocated if (size(aux_b,dim=1) == m_problem) then @@ -180,7 +181,9 @@ program psb_sf_sample write(psb_out_unit,'("Partition type: graph vector")') write(psb_out_unit,'(" ")') ! write(psb_err_unit,'("Build type: graph")') - call build_mtpart(aux_a,np) + call aux_a%cscnv(info,type='csr') + lnp = np + call build_mtpart(aux_a,lnp) endif call psb_barrier(ictxt) diff --git a/test/fileread/psb_zf_sample.f90 b/test/fileread/psb_zf_sample.f90 index a2d54e6a..860405ba 100644 --- a/test/fileread/psb_zf_sample.f90 +++ b/test/fileread/psb_zf_sample.f90 @@ -38,10 +38,11 @@ program psb_zf_sample implicit none ! input parameters - character(len=40) :: kmethd, ptype, mtrx_file, rhs_file,renum + character(len=40) :: kmethd, ptype, mtrx_file, rhs_file ! sparse matrices - type(psb_zspmat_type) :: a, aux_a + type(psb_zspmat_type) :: a + type(psb_lzspmat_type) :: aux_a ! preconditioner data type(psb_zprec_type) :: prec @@ -56,6 +57,7 @@ program psb_zf_sample type(psb_desc_type):: desc_a integer(psb_ipk_) :: ictxt, iam, np + integer(psb_lpk_) :: lnp ! solver paramters integer(psb_ipk_) :: iter, itmax, ierr, itrace, ircode,& @@ -140,7 +142,6 @@ program psb_zf_sample m_problem = aux_a%get_nrows() call psb_bcast(ictxt,m_problem) - call psb_mat_renum(psb_mat_renum_identity_,aux_a,info,perm) ! At this point aux_b may still be unallocated if (size(aux_b,dim=1) == m_problem) then @@ -180,7 +181,9 @@ program psb_zf_sample write(psb_out_unit,'("Partition type: graph vector")') write(psb_out_unit,'(" ")') ! write(psb_err_unit,'("Build type: graph")') - call build_mtpart(aux_a,np) + call aux_a%cscnv(info,type='csr') + lnp = np + call build_mtpart(aux_a,lnp) endif call psb_barrier(ictxt) diff --git a/util/Makefile b/util/Makefile index 3f2f3553..9809a34f 100644 --- a/util/Makefile +++ b/util/Makefile @@ -38,7 +38,7 @@ $(HERE)/$(LIBNAME): $(OBJS) $(OBJS): $(MODDIR)/$(BASEMODNAME)$(.mod) psb_util_mod.o: $(BASEOBJS) -psb_metispart_mod.o: metis_int.o +psb_metispart_mod.o: psb_metis_int.o psb_mat_dist_mod.o: psb_s_mat_dist_mod.o psb_d_mat_dist_mod.o psb_c_mat_dist_mod.o psb_z_mat_dist_mod.o $(IMPLOBJS): $(BASEOBJS) diff --git a/util/psb_metis_int.c b/util/psb_metis_int.c index ad3abec9..0b76c2a6 100644 --- a/util/psb_metis_int.c +++ b/util/psb_metis_int.c @@ -1,31 +1,10 @@ #include #if defined(HAVE_METIS_) -#include "metis.h" +#include "psb_metis_int.h" -typedef int32_t psb_m_t; - -#if defined(IPK4) && defined(LPK4) -typedef int32_t psb_i_t; -typedef int32_t psb_l_t; -#elif defined(IPK4) && defined(LPK8) -typedef int32_t psb_i_t; -typedef int64_t psb_l_t; -#elif defined(IPK8) && defined(LPK8) -typedef int64_t psb_i_t; -typedef int64_t psb_l_t; -#else -#endif -typedef int64_t psb_e_t; - -typedef float psb_s_t; -typedef double psb_d_t; -typedef float complex psb_c_t; -typedef double complex psb_z_t; - - -int metis_PartGraphKway_C(int *n, int *ixadj, int *iadj, int *ivwg, - int *iajw, int *nparts, float *weights, - int *graphpart) +int metis_PartGraphKway_C(idx_t *n, idx_t *ixadj, idx_t *iadj, idx_t *ivwg, + idx_t *iajw, idx_t *nparts, float *weights, + idx_t *graphpart) { int res = -1; idx_t objval = 0; @@ -58,9 +37,9 @@ int metis_PartGraphKway_C(int *n, int *ixadj, int *iadj, int *ivwg, #else -int metis_PartGraphKway_C(int *n, int *ixadj, int *iadj, int *ivwg, - int *iajw, int *nparts, float *weights, - int *graphpart) +int metis_PartGraphKway_C(idx_t *n, idx_t *ixadj, idx_t *iadj, idx_t *ivwg, + idx_t *iajw, idx_t *nparts, float *weights, + idx_t *graphpart) { return(-1); } diff --git a/util/psb_metispart_mod.F90 b/util/psb_metispart_mod.F90 index 425d8106..fa6581e9 100644 --- a/util/psb_metispart_mod.F90 +++ b/util/psb_metispart_mod.F90 @@ -54,31 +54,32 @@ ! uses information prepared by the previous two subroutines. ! module psb_metispart_mod - use psb_base_mod, only : psb_sspmat_type, psb_cspmat_type,& - & psb_dspmat_type, psb_zspmat_type, psb_err_unit, & - & psb_ipk_, psb_lpk_, psb_mpk_, psb_epk_, psb_spk_,& - & psb_s_csr_sparse_mat, psb_d_csr_sparse_mat, & - & psb_c_csr_sparse_mat, psb_z_csr_sparse_mat + use psb_base_mod, only : psb_ipk_, psb_lpk_, psb_mpk_, psb_epk_, & + & psb_err_unit, psb_spk_,& + & psb_lsspmat_type, psb_lcspmat_type,& + & psb_ldspmat_type, psb_lzspmat_type, & + & psb_ls_csr_sparse_mat, psb_ld_csr_sparse_mat, & + & psb_lc_csr_sparse_mat, psb_lz_csr_sparse_mat public part_graph, build_mtpart, distr_mtpart,& & getv_mtpart, free_part private - integer(psb_ipk_), allocatable, save :: graph_vect(:) + integer(psb_lpk_), allocatable, save :: graph_vect(:) interface build_mtpart - module procedure d_mat_build_mtpart, s_mat_build_mtpart,& - & z_mat_build_mtpart, c_mat_build_mtpart + module procedure ld_mat_build_mtpart, ls_mat_build_mtpart,& + & lz_mat_build_mtpart, lc_mat_build_mtpart end interface - interface - subroutine psi_build_mtpart(n,ja,irp,nparts,vect, weights) - import :: psb_ipk_, psb_spk_ + interface psi_build_mtpart + subroutine psi_l_build_mtpart(n,ja,irp,nparts,vect, weights) + import :: psb_lpk_, psb_spk_ implicit none - integer(psb_ipk_), intent(in) :: n, nparts - integer(psb_ipk_), intent(in) :: ja(:), irp(:) - integer(psb_ipk_), allocatable, intent(inout) :: vect(:) + integer(psb_lpk_), intent(in) :: n, nparts + integer(psb_lpk_), intent(in) :: ja(:), irp(:) + integer(psb_lpk_), allocatable, intent(inout) :: vect(:) real(psb_spk_),optional, intent(in) :: weights(:) - end subroutine psi_build_mtpart + end subroutine psi_l_build_mtpart end interface contains @@ -155,29 +156,30 @@ contains ivg(:) = graph_vect(:) end if end subroutine getv_mtpart - - subroutine d_mat_build_mtpart(a,nparts,weights) + + + subroutine ld_mat_build_mtpart(a,nparts,weights) use psb_base_mod implicit none - type(psb_dspmat_type), intent(in) :: a - integer(psb_ipk_) :: nparts + type(psb_ldspmat_type), intent(in) :: a + integer(psb_lpk_) :: nparts real(psb_dpk_), optional :: weights(:) select type (aa=>a%a) - type is (psb_d_csr_sparse_mat) - call d_csr_build_mtpart(aa,nparts,weights) + type is (psb_ld_csr_sparse_mat) + call ld_csr_build_mtpart(aa,nparts,weights) class default write(psb_err_unit,*) 'Sorry, right now we only take CSR input!' end select - end subroutine d_mat_build_mtpart + end subroutine ld_mat_build_mtpart - subroutine d_csr_build_mtpart(a,nparts,weights) + subroutine ld_csr_build_mtpart(a,nparts,weights) use psb_base_mod implicit none - type(psb_d_csr_sparse_mat), intent(in) :: a - integer(psb_ipk_) :: nparts + type(psb_ld_csr_sparse_mat), intent(in) :: a + integer(psb_lpk_) :: nparts real(psb_dpk_), optional :: weights(:) real(psb_spk_), allocatable :: wgh_(:) @@ -193,30 +195,30 @@ contains call psi_build_mtpart(a%get_nrows(),a%ja,a%irp,nparts,graph_vect) end if - end subroutine d_csr_build_mtpart + end subroutine ld_csr_build_mtpart - subroutine z_mat_build_mtpart(a,nparts,weights) + subroutine lz_mat_build_mtpart(a,nparts,weights) use psb_base_mod implicit none - type(psb_zspmat_type), intent(in) :: a - integer(psb_ipk_) :: nparts + type(psb_lzspmat_type), intent(in) :: a + integer(psb_lpk_) :: nparts real(psb_dpk_), optional :: weights(:) select type (aa=>a%a) - type is (psb_z_csr_sparse_mat) - call z_csr_build_mtpart(aa,nparts,weights) + type is (psb_lz_csr_sparse_mat) + call lz_csr_build_mtpart(aa,nparts,weights) class default write(psb_err_unit,*) 'Sorry, right now we only take CSR input!' end select - end subroutine z_mat_build_mtpart + end subroutine lz_mat_build_mtpart - subroutine z_csr_build_mtpart(a,nparts,weights) + subroutine lz_csr_build_mtpart(a,nparts,weights) use psb_base_mod implicit none - type(psb_z_csr_sparse_mat), intent(in) :: a - integer(psb_ipk_) :: nparts + type(psb_lz_csr_sparse_mat), intent(in) :: a + integer(psb_lpk_) :: nparts real(psb_dpk_), optional :: weights(:) real(psb_spk_), allocatable :: wgh_(:) @@ -232,66 +234,66 @@ contains call psi_build_mtpart(a%get_nrows(),a%ja,a%irp,nparts,graph_vect) end if - end subroutine z_csr_build_mtpart + end subroutine lz_csr_build_mtpart - subroutine s_mat_build_mtpart(a,nparts,weights) + subroutine ls_mat_build_mtpart(a,nparts,weights) use psb_base_mod implicit none - type(psb_sspmat_type), intent(in) :: a - integer(psb_ipk_) :: nparts + type(psb_lsspmat_type), intent(in) :: a + integer(psb_lpk_) :: nparts real(psb_spk_), optional :: weights(:) select type (aa=>a%a) - type is (psb_s_csr_sparse_mat) - call s_csr_build_mtpart(aa,nparts,weights) + type is (psb_ls_csr_sparse_mat) + call ls_csr_build_mtpart(aa,nparts,weights) class default write(psb_err_unit,*) 'Sorry, right now we only take CSR input!' end select - end subroutine s_mat_build_mtpart + end subroutine ls_mat_build_mtpart - subroutine c_mat_build_mtpart(a,nparts,weights) + subroutine lc_mat_build_mtpart(a,nparts,weights) use psb_base_mod implicit none - type(psb_cspmat_type), intent(in) :: a - integer(psb_ipk_) :: nparts + type(psb_lcspmat_type), intent(in) :: a + integer(psb_lpk_) :: nparts real(psb_spk_), optional :: weights(:) select type (aa=>a%a) - type is (psb_c_csr_sparse_mat) - call c_csr_build_mtpart(aa,nparts,weights) + type is (psb_lc_csr_sparse_mat) + call lc_csr_build_mtpart(aa,nparts,weights) class default write(psb_err_unit,*) 'Sorry, right now we only take CSR input!' end select - end subroutine c_mat_build_mtpart + end subroutine lc_mat_build_mtpart - subroutine c_csr_build_mtpart(a,nparts,weights) + subroutine lc_csr_build_mtpart(a,nparts,weights) use psb_base_mod implicit none - type(psb_c_csr_sparse_mat), intent(in) :: a - integer(psb_ipk_) :: nparts + type(psb_lc_csr_sparse_mat), intent(in) :: a + integer(psb_lpk_) :: nparts real(psb_spk_), optional :: weights(:) call psi_build_mtpart(a%get_nrows(),a%ja,a%irp,nparts,graph_vect,weights) - end subroutine c_csr_build_mtpart + end subroutine lc_csr_build_mtpart - subroutine s_csr_build_mtpart(a,nparts,weights) + subroutine ls_csr_build_mtpart(a,nparts,weights) use psb_base_mod implicit none - type(psb_s_csr_sparse_mat), intent(in) :: a - integer(psb_ipk_) :: nparts + type(psb_ls_csr_sparse_mat), intent(in) :: a + integer(psb_lpk_) :: nparts real(psb_spk_), optional :: weights(:) - + call psi_build_mtpart(a%get_nrows(),a%ja,a%irp,nparts,graph_vect,weights) - end subroutine s_csr_build_mtpart - + end subroutine ls_csr_build_mtpart + ! ! WARNING: called IRET otherwise Intel compiler complains, ! methinks it's a compiler bug, will need to report. diff --git a/util/psi_build_mtpart.F90 b/util/psi_build_mtpart.F90 index b8974185..d01c3d85 100644 --- a/util/psi_build_mtpart.F90 +++ b/util/psi_build_mtpart.F90 @@ -1,31 +1,21 @@ - -subroutine psi_build_mtpart(n,ja,irp,nparts,graph_vect,weights) +subroutine psi_l_build_mtpart(n,ja,irp,nparts,graph_vect,weights) use psb_base_mod use iso_c_binding implicit none - integer(psb_ipk_), intent(in) :: n, nparts - integer(psb_ipk_), intent(in) :: ja(:), irp(:) - integer(psb_ipk_), allocatable, intent(inout) :: graph_vect(:) + integer(psb_lpk_), intent(in) :: n, nparts + integer(psb_lpk_), intent(in) :: ja(:), irp(:) + integer(psb_lpk_), allocatable, intent(inout) :: graph_vect(:) real(psb_spk_),optional, intent(in) :: weights(:) ! local variables - integer(psb_ipk_) :: i,numflag, nedc,wgflag - integer(psb_ipk_) :: iopt(10),idummy(2),jdummy(2), info - integer(psb_ipk_) :: nl,nptl - integer(psb_ipk_), allocatable :: irpl(:),jal(:),gvl(:) + integer(psb_lpk_) :: i,numflag, nedc,wgflag + integer(psb_lpk_) :: iopt(10),idummy(2),jdummy(2) + integer(psb_ipk_) :: info + integer(psb_lpk_) :: nl,nptl + integer(psb_lpk_), allocatable :: irpl(:),jal(:),gvl(:) real(psb_spk_),allocatable :: wgh_(:) -#if defined(HAVE_METIS) && defined(IPK4) +#if defined(HAVE_METIS) && defined(LPK4) && defined(METIS_32) interface - ! subroutine METIS_PartGraphKway(n,ixadj,iadj,ivwg,iajw,& - ! & wgflag,numflag,nparts,weights,iopt,nedc,part) bind(c) - ! use iso_c_binding - ! integer(c_int) :: n,wgflag,numflag,nparts,nedc - ! integer(c_int) :: ixadj(*),iadj(*),ivwg(*),iajw(*),iopt(*),part(*) - ! real(c_float) :: weights(*) - ! !integer(psb_ipk_) :: n,wgflag,numflag,nparts,nedc - ! !integer(psb_ipk_) :: ixadj(*),iadj(*),ivwg(*),iajw(*),iopt(*),part(*) - ! end subroutine METIS_PartGraphKway - function METIS_PartGraphKway(n,ixadj,iadj,ivwg,iajw,& & nparts,weights,part) bind(c,name="metis_PartGraphKway_C") result(res) use iso_c_binding @@ -86,10 +76,76 @@ subroutine psi_build_mtpart(n,ja,irp,nparts,graph_vect,weights) graph_vect(i) = 0 enddo endif + +#elif defined(HAVE_METIS) && defined(LPK8) && defined(METIS_64) + + interface + function METIS_PartGraphKway(n,ixadj,iadj,ivwg,iajw,& + & nparts,weights,part) bind(c,name="metis_PartGraphKway_C") result(res) + use iso_c_binding + integer(c_long_long) :: res + integer(c_long_long) :: n,nparts + integer(c_long_long) :: ixadj(*),iadj(*),ivwg(*),iajw(*),part(*) + real(c_float) :: weights(*) + !integer(psb_ipk_) :: n,wgflag,numflag,nparts,nedc + !integer(psb_ipk_) :: ixadj(*),iadj(*),ivwg(*),iajw(*),iopt(*),part(*) + end function METIS_PartGraphKway + end interface + + call psb_realloc(n,graph_vect,info) + if (info == psb_success_) allocate(gvl(n),wgh_(nparts),stat=info) + + if (info /= psb_success_) then + write(psb_err_unit,*) 'Fatal error in BUILD_MTPART: memory allocation ',& + & ' failure.' + return + endif + if (nparts > 1) then + iopt(1) = 0 + numflag = 1 + wgflag = 0 + +!!$ write(*,*) 'Before allocation',nparts + + irpl=irp + jal = ja + nl = n + nptl = nparts + wgh_ = -1.0 + if(present(weights)) then + if (size(weights) == nptl) then +!!$ write(*,*) 'weights present',weights + ! call METIS_PartGraphKway(n,irp,ja,idummy,jdummy,& + ! & wgflag,numflag,nparts,weights,iopt,nedc,graph_vect) + info = METIS_PartGraphKway(nl,irpl,jal,idummy,jdummy,& + & nptl,weights,gvl) + + else +!!$ write(*,*) 'weights absent',wgh_ + info = METIS_PartGraphKway(nl,irpl,jal,idummy,jdummy,& + & nptl,wgh_,gvl) + end if + else +!!$ write(*,*) 'weights absent',wgh_ + info = METIS_PartGraphKway(nl,irpl,jal,idummy,jdummy,& + & nptl,wgh_,gvl) + endif +!!$ write(*,*) 'after allocation',info + + do i=1, n + graph_vect(i) = gvl(i) - 1 + enddo + else + do i=1, n + graph_vect(i) = 0 + enddo + endif + #else - write(psb_err_unit,*) 'Warning: METIS was not configured at PSBLAS compile time !' + + write(psb_err_unit,*) 'Warning: no suitable METIS interface for LPK indices' #endif return -end subroutine psi_build_mtpart +end subroutine psi_l_build_mtpart