!!$ !!$ !!$ MD2P4 !!$ Multilevel Domain Decomposition Parallel Preconditioner Package for PSBLAS !!$ for !!$ Parallel Sparse BLAS v2.0 !!$ (C) Copyright 2006 Salvatore Filippone University of Rome Tor Vergata !!$ Alfredo Buttari University of Rome Tor Vergata !!$ Daniela di Serafino Second University of Naples !!$ Pasqua D'Ambra ICAR-CNR !!$ !!$ Redistribution and use in source and binary forms, with or without !!$ modification, are permitted provided that the following conditions !!$ are met: !!$ 1. Redistributions of source code must retain the above copyright !!$ notice, this list of conditions and the following disclaimer. !!$ 2. Redistributions in binary form must reproduce the above copyright !!$ notice, this list of conditions, and the following disclaimer in the !!$ documentation and/or other materials provided with the distribution. !!$ 3. The name of the MD2P4 group or the names of its contributors may !!$ not be used to endorse or promote products derived from this !!$ software without specific written permission. !!$ !!$ THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS !!$ ``AS IS'' AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED !!$ TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR !!$ PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE MD2P4 GROUP OR ITS CONTRIBUTORS !!$ BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR !!$ CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF !!$ SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS !!$ INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN !!$ CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) !!$ ARISING IN ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE !!$ POSSIBILITY OF SUCH DAMAGE. !!$ !!$ subroutine mld_zaggrmat_smth_asb(a,desc_a,ac,desc_ac,p,info) use psb_base_mod use mld_prec_mod, mld_protect_name => mld_zaggrmat_smth_asb #ifdef MPI_MOD use mpi #endif implicit none #ifdef MPI_H include 'mpif.h' #endif type(psb_zspmat_type), intent(in), target :: a type(mld_zbaseprc_type), intent(inout),target :: p type(psb_zspmat_type), intent(inout), target :: ac type(psb_desc_type), intent(in) :: desc_a type(psb_desc_type), intent(inout) :: desc_ac integer, intent(out) :: info type(psb_zspmat_type) :: b integer, pointer :: nzbr(:), idisp(:) integer :: nrow, nglob, ncol, ntaggr, nzac, ip, ndx,& & naggr, nzt,jl,nzl,nlr,naggrm1,naggrp1, i, j, k logical, parameter :: aggr_dump=.false. integer ::ictxt,np,me, err_act, icomm character(len=20) :: name, ch_err type(psb_zspmat_type), pointer :: am1,am2 type(psb_zspmat_type) :: am3,am4 logical :: ml_global_nmb logical, parameter :: test_dump=.false., debug=.false. integer, parameter :: ncmax=16 real(kind(1.d0)) :: omega, anorm, tmp, dg name='mld_aggrmat_smth_asb' if(psb_get_errstatus().ne.0) return info=0 call psb_erractionsave(err_act) ictxt = psb_cd_get_context(desc_a) icomm = psb_cd_get_mpic(desc_a) ictxt=psb_cd_get_context(desc_a) call psb_info(ictxt, me, np) call psb_nullify_sp(b) call psb_nullify_sp(am3) call psb_nullify_sp(am4) am2 => p%av(sm_pr_t_) am1 => p%av(sm_pr_) call psb_nullify_sp(am1) call psb_nullify_sp(am2) nglob = psb_cd_get_global_rows(desc_a) nrow = psb_cd_get_local_rows(desc_a) ncol = psb_cd_get_local_cols(desc_a) naggr = p%nlaggr(me+1) ntaggr = sum(p%nlaggr) allocate(nzbr(np), idisp(np),stat=info) if (info /= 0) then info=4025 call psb_errpush(info,name,i_err=(/2*np,0,0,0,0/),& & a_err='integer') goto 9999 end if naggrm1 = sum(p%nlaggr(1:me)) naggrp1 = sum(p%nlaggr(1:me+1)) ml_global_nmb = ( (p%iprcparm(aggr_kind_) == smooth_prol_).or.& & ( (p%iprcparm(aggr_kind_) == biz_prol_).and.& & (p%iprcparm(coarse_mat_) == repl_mat_)) ) if (ml_global_nmb) then p%mlia(1:nrow) = p%mlia(1:nrow) + naggrm1 call psb_halo(p%mlia,desc_a,info) if(info /= 0) then call psb_errpush(4010,name,a_err='psb_halo') goto 9999 end if end if if (aggr_dump) then open(30+me) write(30+me,*) '% Aggregation map' do i=1,ncol write(30+me,*) i,p%mlia(i) end do close(30+me) end if ! naggr: number of local aggregates ! nrow: local rows. ! allocate(p%dorig(nrow),stat=info) if (info /= 0) then info=4025 call psb_errpush(info,name,i_err=(/nrow,0,0,0,0/),& & a_err='real(kind(1.d0))') goto 9999 end if ! Get diagonal D call psb_sp_getdiag(a,p%dorig,info) if(info /= 0) then call psb_errpush(4010,name,a_err='sp_getdiag') goto 9999 end if do i=1,size(p%dorig) if (p%dorig(i) /= zzero) then p%dorig(i) = zone / p%dorig(i) else p%dorig(i) = zone end if end do ! where (p%dorig /= zzero) ! p%dorig = zone / p%dorig ! elsewhere ! p%dorig = zone ! end where ! 1. Allocate Ptilde in sparse matrix form am4%fida='COO' am4%m=ncol if (ml_global_nmb) then am4%k=ntaggr call psb_sp_all(ncol,ntaggr,am4,ncol,info) else am4%k=naggr call psb_sp_all(ncol,naggr,am4,ncol,info) endif if(info /= 0) then call psb_errpush(4010,name,a_err='spall') goto 9999 end if if (ml_global_nmb) then do i=1,ncol am4%aspk(i) = zone am4%ia1(i) = i am4%ia2(i) = p%mlia(i) end do am4%infoa(psb_nnz_) = ncol else do i=1,nrow am4%aspk(i) = zone am4%ia1(i) = i am4%ia2(i) = p%mlia(i) end do am4%infoa(psb_nnz_) = nrow endif call psb_spcnv(am4,info,afmt='csr',dupl=psb_dupl_add_) if(info /= 0) then call psb_errpush(4010,name,a_err='spcnv') goto 9999 end if call psb_sp_clone(a,am3,info) if(info /= 0) then call psb_errpush(4010,name,a_err='spclone') goto 9999 end if ! ! WARNING: the cycles below assume that AM3 does have ! its diagonal elements stored explicitly!!! ! Should we switch to something safer? ! call psb_sp_scal(am3,p%dorig,info) if(info /= 0) goto 9999 if (p%iprcparm(aggr_eig_) == max_norm_) then if (p%iprcparm(aggr_kind_) == biz_prol_) then ! ! This only works with CSR. ! anorm = dzero dg = done do i=1,am3%m tmp = dzero do j=am3%ia2(i),am3%ia2(i+1)-1 if (am3%ia1(j) <= am3%m) then tmp = tmp + abs(am3%aspk(j)) endif if (am3%ia1(j) == i ) then dg = abs(am3%aspk(j)) end if end do anorm = max(anorm,tmp/dg) enddo call psb_amx(ictxt,anorm) else anorm = psb_spnrmi(am3,desc_a,info) endif omega = 4.d0/(3.d0*anorm) p%dprcparm(aggr_damp_) = omega else if (p%iprcparm(aggr_eig_) == user_choice_) then omega = p%dprcparm(aggr_damp_) else if (p%iprcparm(aggr_eig_) /= user_choice_) then write(0,*) me,'Error: invalid choice for OMEGA in blaggrmat?? ',& & p%iprcparm(aggr_eig_) end if if (toupper(am3%fida)=='CSR') then do i=1,am3%m do j=am3%ia2(i),am3%ia2(i+1)-1 if (am3%ia1(j) == i) then am3%aspk(j) = zone - omega*am3%aspk(j) else am3%aspk(j) = - omega*am3%aspk(j) end if end do end do else if (toupper(am3%fida)=='COO') then do j=1,am3%infoa(psb_nnz_) if (am3%ia1(j) /= am3%ia2(j)) then am3%aspk(j) = - omega*am3%aspk(j) else am3%aspk(j) = zone - omega*am3%aspk(j) endif end do call psb_spcnv(am3,info,afmt='csr',dupl=psb_dupl_add_) if (info /=0) then call psb_errpush(4010,name,a_err='spcnv am3') goto 9999 end if else write(0,*) 'Missing implementation of I sum' call psb_errpush(4010,name) goto 9999 end if if (test_dump) then open(30+me) write(30+me,*) 'OMEGA: ',omega do i=1,size(p%dorig) write(30+me,*) p%dorig(i) end do close(30+me) end if if (test_dump) call & & psb_csprt(20+me,am4,head='% Operator Ptilde.',ivr=desc_a%loc_to_glob) if (test_dump) call psb_csprt(40+me,am3,head='% (I-wDA)',ivr=desc_a%loc_to_glob,& & ivc=desc_a%loc_to_glob) if (debug) write(0,*) me,'Done gather, going for SYMBMM 1' ! ! Symbmm90 does the allocation for its result. ! ! am1 = (i-wDA)Ptilde ! Doing it this way means to consider diag(Ai) ! ! call psb_symbmm(am3,am4,am1,info) if(info /= 0) then call psb_errpush(4010,name,a_err='symbmm 1') goto 9999 end if call psb_numbmm(am3,am4,am1) if (debug) write(0,*) me,'Done NUMBMM 1' call psb_sp_free(am4,info) if(info /= 0) then call psb_errpush(4010,name,a_err='sp_free') goto 9999 end if if (ml_global_nmb) then ! ! Now we have to gather the halo of am1, and add it to itself ! to multiply it by A, ! call psb_sphalo(am1,desc_a,am4,info,& & colcnv=.false.,rowscale=.true.) if(info /= 0) then call psb_errpush(4010,name,a_err='psb_sphalo') goto 9999 end if call psb_rwextd(ncol,am1,info,b=am4) if(info /= 0) then call psb_errpush(4010,name,a_err='psb_rwextd') goto 9999 end if call psb_sp_free(am4,info) if(info /= 0) then call psb_errpush(4010,name,a_err='psb_sp_free') goto 9999 end if else call psb_rwextd(ncol,am1,info) if(info /= 0) then call psb_errpush(4010,name,a_err='rwextd') goto 9999 end if endif if (test_dump) & & call psb_csprt(60+me,am1,head='% (I-wDA)Pt',ivr=desc_a%loc_to_glob) call psb_symbmm(a,am1,am3,info) if(info /= 0) then call psb_errpush(4010,name,a_err='symbmm 2') goto 9999 end if call psb_numbmm(a,am1,am3) if (debug) write(0,*) me,'Done NUMBMM 2' if (p%iprcparm(aggr_kind_) == smooth_prol_) then call psb_transp(am1,am2,fmt='COO') nzl = am2%infoa(psb_nnz_) i=0 ! ! Now we have to fix this. The only rows of B that are correct ! are those corresponding to "local" aggregates, i.e. indices in p%mlia(:) ! do k=1, nzl if ((naggrm1 < am2%ia1(k)) .and.(am2%ia1(k) <= naggrp1)) then i = i+1 am2%aspk(i) = am2%aspk(k) am2%ia1(i) = am2%ia1(k) am2%ia2(i) = am2%ia2(k) end if end do am2%infoa(psb_nnz_) = i call psb_spcnv(am2,info,afmt='csr',dupl=psb_dupl_add_) if (info /=0) then call psb_errpush(4010,name,a_err='spcnv am2') goto 9999 end if else call psb_transp(am1,am2) endif if (debug) write(0,*) me,'starting sphalo/ rwxtd' if (p%iprcparm(aggr_kind_) == smooth_prol_) then ! am2 = ((i-wDA)Ptilde)^T call psb_sphalo(am3,desc_a,am4,info,& & colcnv=.false.,rowscale=.true.) if(info /= 0) then call psb_errpush(4010,name,a_err='psb_sphalo') goto 9999 end if call psb_rwextd(ncol,am3,info,b=am4) if(info /= 0) then call psb_errpush(4010,name,a_err='psb_rwextd') goto 9999 end if call psb_sp_free(am4,info) if(info /= 0) then call psb_errpush(4010,name,a_err='psb_sp_free') goto 9999 end if else if (p%iprcparm(aggr_kind_) == biz_prol_) then call psb_rwextd(ncol,am3,info) if(info /= 0) then call psb_errpush(4010,name,a_err='psb_rwextd') goto 9999 end if endif if (debug) write(0,*) me,'starting symbmm 3' call psb_symbmm(am2,am3,b,info) if(info /= 0) then call psb_errpush(4010,name,a_err='symbmm 3') goto 9999 end if if (debug) write(0,*) me,'starting numbmm 3' call psb_numbmm(am2,am3,b) if (debug) write(0,*) me,'Done NUMBMM 3' !!$ if (aggr_dump) call csprt(50+me,am1,head='% Operator PTrans.') call psb_sp_free(am3,info) if(info /= 0) then call psb_errpush(4010,name,a_err='psb_sp_free') goto 9999 end if call psb_spcnv(b,info,afmt='coo',dupl=psb_dupl_add_) if (info /=0) then call psb_errpush(4010,name,a_err='spcnv b') goto 9999 end if if (test_dump) call psb_csprt(80+me,b,head='% Smoothed aggregate AC.') select case(p%iprcparm(aggr_kind_)) case(smooth_prol_) select case(p%iprcparm(coarse_mat_)) case(distr_mat_) call psb_sp_clone(b,ac,info) if(info /= 0) goto 9999 nzac = ac%infoa(psb_nnz_) nzl = ac%infoa(psb_nnz_) call psb_cdall(ictxt,desc_ac,info,nl=p%nlaggr(me+1)) if(info /= 0) then call psb_errpush(4010,name,a_err='psb_cdall') goto 9999 end if call psb_cdins(nzl,ac%ia1,ac%ia2,desc_ac,info) if(info /= 0) then call psb_errpush(4010,name,a_err='psb_cdins') goto 9999 end if if (debug) write(0,*) me,'Created aux descr. distr.' call psb_cdasb(desc_ac,info) if(info /= 0) then call psb_errpush(4010,name,a_err='psb_cdasb') goto 9999 end if if (debug) write(0,*) me,'Asmbld aux descr. distr.' call psb_glob_to_loc(ac%ia1(1:nzl),desc_ac,info,iact='I') if(info /= 0) then call psb_errpush(4010,name,a_err='psglob_to_loc') goto 9999 end if call psb_glob_to_loc(ac%ia2(1:nzl),desc_ac,info,iact='I') if(info /= 0) then call psb_errpush(4010,name,a_err='psglob_to_loc') goto 9999 end if ac%m=desc_ac%matrix_data(psb_n_row_) ac%k=desc_ac%matrix_data(psb_n_col_) ac%fida='COO' ac%descra='G' call psb_sp_free(b,info) if (info == 0) deallocate(nzbr,idisp,stat=info) if(info /= 0) then call psb_errpush(4010,name,a_err='psb_sp_free') goto 9999 end if if (np>1) then nzl = psb_sp_get_nnzeros(am1) call psb_glob_to_loc(am1%ia1(1:nzl),desc_ac,info,'I') if(info /= 0) then call psb_errpush(4010,name,a_err='psb_glob_to_loc') goto 9999 end if endif am1%k=desc_ac%matrix_data(psb_n_col_) if (np>1) then call psb_spcnv(am2,info,afmt='coo',dupl=psb_dupl_add_) if(info /= 0) then call psb_errpush(4010,name,a_err='spcnv') goto 9999 end if nzl = am2%infoa(psb_nnz_) call psb_glob_to_loc(am2%ia1(1:nzl),desc_ac,info,'I') if(info /= 0) then call psb_errpush(4010,name,a_err='psb_glob_to_loc') goto 9999 end if call psb_spcnv(am2,info,afmt='csr',dupl=psb_dupl_add_) if(info /= 0) then call psb_errpush(4010,name,a_err='spcnv') goto 9999 end if end if am2%m=desc_ac%matrix_data(psb_n_col_) if (debug) write(0,*) me,'Done ac ' case(repl_mat_) ! ! call psb_cdrep(ntaggr,ictxt,desc_ac,info) if(info /= 0) then call psb_errpush(4010,name,a_err='psb_cdrep') goto 9999 end if nzbr(:) = 0 nzbr(me+1) = b%infoa(psb_nnz_) call psb_sum(ictxt,nzbr(1:np)) nzac = sum(nzbr) call psb_sp_all(ntaggr,ntaggr,ac,nzac,info) if(info /= 0) goto 9999 do ip=1,np idisp(ip) = sum(nzbr(1:ip-1)) enddo ndx = nzbr(me+1) call mpi_allgatherv(b%aspk,ndx,mpi_double_complex,ac%aspk,nzbr,idisp,& & mpi_double_complex,icomm,info) call mpi_allgatherv(b%ia1,ndx,mpi_integer,ac%ia1,nzbr,idisp,& & mpi_integer,icomm,info) call mpi_allgatherv(b%ia2,ndx,mpi_integer,ac%ia2,nzbr,idisp,& & mpi_integer,icomm,info) if(info /= 0) goto 9999 ac%m = ntaggr ac%k = ntaggr ac%infoa(psb_nnz_) = nzac ac%fida='COO' ac%descra='G' call psb_spcnv(ac,info,afmt='coo',dupl=psb_dupl_add_) if(info /= 0) goto 9999 call psb_sp_free(b,info) if(info /= 0) goto 9999 if (me==0) then if (test_dump) call psb_csprt(80+me,ac,head='% Smoothed aggregate AC.') endif deallocate(nzbr,idisp) case default write(0,*) 'Inconsistent input in smooth_new_aggregate' end select case(biz_prol_) select case(p%iprcparm(coarse_mat_)) case(distr_mat_) call psb_sp_clone(b,ac,info) if(info /= 0) then call psb_errpush(4010,name,a_err='spclone') goto 9999 end if call psb_cdall(ictxt,desc_ac,info,nl=naggr) if(info /= 0) then call psb_errpush(4010,name,a_err='psb_cdall') goto 9999 end if call psb_cdasb(desc_ac,info) if(info /= 0) then call psb_errpush(4010,name,a_err='psb_cdasb') goto 9999 end if call psb_sp_free(b,info) if(info /= 0) then call psb_errpush(4010,name,a_err='sp_free') goto 9999 end if case(repl_mat_) ! ! call psb_cdrep(ntaggr,ictxt,desc_ac,info) if(info /= 0) then call psb_errpush(4010,name,a_err='psb_cdrep') goto 9999 end if nzbr(:) = 0 nzbr(me+1) = b%infoa(psb_nnz_) call psb_sum(ictxt,nzbr(1:np)) nzac = sum(nzbr) call psb_sp_all(ntaggr,ntaggr,ac,nzac,info) if(info /= 0) then call psb_errpush(4010,name,a_err='psb_sp_all') goto 9999 end if do ip=1,np idisp(ip) = sum(nzbr(1:ip-1)) enddo ndx = nzbr(me+1) call mpi_allgatherv(b%aspk,ndx,mpi_double_complex,ac%aspk,nzbr,idisp,& & mpi_double_complex,icomm,info) call mpi_allgatherv(b%ia1,ndx,mpi_integer,ac%ia1,nzbr,idisp,& & mpi_integer,icomm,info) call mpi_allgatherv(b%ia2,ndx,mpi_integer,ac%ia2,nzbr,idisp,& & mpi_integer,icomm,info) if(info /= 0) then info=-1 call psb_errpush(info,name) goto 9999 end if ac%m = ntaggr ac%k = ntaggr ac%infoa(psb_nnz_) = nzac ac%fida='COO' ac%descra='G' call psb_spcnv(ac,info,afmt='coo',dupl=psb_dupl_add_) if(info /= 0) then call psb_errpush(4010,name,a_err='spcnv') goto 9999 end if call psb_sp_free(b,info) if(info /= 0) then call psb_errpush(4010,name,a_err='psb_sp_free') goto 9999 end if end select deallocate(nzbr,idisp) end select call psb_spcnv(ac,info,afmt='csr',dupl=psb_dupl_add_) if(info /= 0) then call psb_errpush(4010,name,a_err='spcnv') goto 9999 end if if (debug) write(0,*) me,'Done smooth_aggregate ' call psb_erractionrestore(err_act) return 9999 continue call psb_errpush(info,name) call psb_erractionrestore(err_act) if (err_act.eq.psb_act_abort_) then call psb_error() return end if return end subroutine mld_zaggrmat_smth_asb