Renaming: moved PSB_DSCFREE into PSB_CDFREE.

psblas3-type-indexed
Salvatore Filippone 19 years ago
parent 49f80f671d
commit 1a6aa93f75

@ -439,7 +439,7 @@ contains
end if
if (associated(p%desc_data)) then
if (associated(p%desc_data%matrix_data)) then
call psb_dscfree(p%desc_data,info)
call psb_cdfree(p%desc_data,info)
end if
deallocate(p%desc_data)
endif

@ -297,21 +297,21 @@ Module psb_tools_mod
end subroutine psb_dptins
end interface
interface psb_dscall
subroutine psb_dscall(m, n, parts, icontxt, desc_a, info)
interface psb_cdall
subroutine psb_cdall(m, n, parts, icontxt, desc_a, info)
use psb_descriptor_type
include 'parts.fh'
Integer, intent(in) :: m,n,icontxt
Type(psb_desc_type), intent(out) :: desc_a
integer, intent(out) :: info
end subroutine psb_dscall
subroutine psb_dscalv(m, v, icontxt, desc_a, info, flag)
end subroutine psb_cdall
subroutine psb_cdalv(m, v, icontxt, desc_a, info, flag)
use psb_descriptor_type
Integer, intent(in) :: m,icontxt, v(:)
integer, intent(in), optional :: flag
integer, intent(out) :: info
Type(psb_desc_type), intent(out) :: desc_a
end subroutine psb_dscalv
end subroutine psb_cdalv
end interface
@ -335,12 +335,12 @@ Module psb_tools_mod
end interface
interface psb_dscfree
subroutine psb_dscfree(desc_a,info)
interface psb_cdfree
subroutine psb_cdfree(desc_a,info)
use psb_descriptor_type
type(psb_desc_type), intent(inout) :: desc_a
integer, intent(out) :: info
end subroutine psb_dscfree
end subroutine psb_cdfree
end interface
interface psb_dscins
@ -500,22 +500,22 @@ Module psb_tools_mod
end interface
interface psb_dscrep
subroutine psb_dscrep(m, icontxt, desc_a,info)
interface psb_cdrep
subroutine psb_cdrep(m, icontxt, desc_a,info)
use psb_descriptor_type
Integer, intent(in) :: m,icontxt
Type(psb_desc_type), intent(out) :: desc_a
integer, intent(out) :: info
end subroutine psb_dscrep
end subroutine psb_cdrep
end interface
interface psb_dscdec
subroutine psb_dscdec(nloc, icontxt, desc_a,info)
interface psb_cddec
subroutine psb_cddec(nloc, icontxt, desc_a,info)
use psb_descriptor_type
Integer, intent(in) :: nloc,icontxt
Type(psb_desc_type), intent(out) :: desc_a
integer, intent(out) :: info
end subroutine psb_dscdec
end subroutine psb_cddec
end interface

@ -239,7 +239,7 @@ contains
if (p%iprcparm(coarse_mat_) == mat_repl_) then
call psb_dscrep(ntaggr,icontxt,p%desc_data,info)
call psb_cdrep(ntaggr,icontxt,p%desc_data,info)
nzbr(:) = 0
nzbr(myprow+1) = irs
@ -288,7 +288,7 @@ contains
else if (p%iprcparm(coarse_mat_) == mat_distr_) then
call psb_dscdec(naggr,icontxt,p%desc_data,info)
call psb_cddec(naggr,icontxt,p%desc_data,info)
call psb_spclone(b,bg,info)
if(info /= 0) then
call psb_errpush(4010,name,a_err='spclone')
@ -715,9 +715,9 @@ contains
i = i + 1
end do
end do
call psb_dscall(ntaggr,ivall,icontxt,p%desc_data,info,flag=1)
call psb_cdall(ntaggr,ivall,icontxt,p%desc_data,info,flag=1)
if(info /= 0) then
call psb_errpush(4010,name,a_err='psb_dscall')
call psb_errpush(4010,name,a_err='psb_cdall')
goto 9999
end if
@ -833,7 +833,7 @@ contains
nzbr(:) = 0
nzbr(myprow+1) = b%infoa(psb_nnz_)
call psb_dscrep(ntaggr,icontxt,p%desc_data,info)
call psb_cdrep(ntaggr,icontxt,p%desc_data,info)
call igsum2d(icontxt,'All',' ',np,1,nzbr,np,-1,-1)
nzbg = sum(nzbr)
@ -887,7 +887,7 @@ contains
call psb_errpush(4010,name,a_err='spclone')
goto 9999
end if
call psb_dscdec(naggr,icontxt,p%desc_data,info)
call psb_cddec(naggr,icontxt,p%desc_data,info)
call psb_spfree(b,info)
if(info /= 0) then
@ -902,7 +902,7 @@ contains
nzbr(:) = 0
nzbr(myprow+1) = b%infoa(psb_nnz_)
call psb_dscrep(ntaggr,icontxt,p%desc_data,info)
call psb_cdrep(ntaggr,icontxt,p%desc_data,info)
call igsum2d(icontxt,'All',' ',np,1,nzbr,np,-1,-1)

@ -2,9 +2,9 @@ include ../../Make.inc
FOBJS = psb_dallc.o psb_dasb.o psb_dcsrp.o psb_cdprt.o \
psb_dfree.o psb_dgelp.o psb_dins.o \
psb_dscall.o psb_dscalv.o psb_dscasb.o psb_dsccpy.o \
psb_dscdec.o psb_dscfree.o psb_dscins.o psb_cdovr.o \
psb_dscren.o psb_dscrep.o psb_dspalloc.o psb_dspasb.o \
psb_cdall.o psb_cdalv.o psb_dscasb.o psb_dsccpy.o \
psb_cddec.o psb_cdfree.o psb_dscins.o psb_cdovr.o \
psb_dscren.o psb_cdrep.o psb_dspalloc.o psb_dspasb.o \
psb_dspcnv.o psb_dspfree.o psb_dspins.o psb_dsprn.o \
psb_glob_to_loc.o psb_ialloc.o psb_iasb.o \
psb_ifree.o psb_iins.o psb_loc_to_glob.o

@ -28,9 +28,9 @@
!!$ POSSIBILITY OF SUCH DAMAGE.
!!$
!!$
! File: psb_dscall.f90
! File: psb_cdall.f90
!
! Subroutine: psb_dscall
! Subroutine: psb_cdall
! Allocate descriptor
! and checks correctness of PARTS subroutine
!
@ -41,7 +41,7 @@
! icontxt - integer. The communication context.
! desc_a - type(<psb_desc_type>). The communication descriptor.
! info - integer. Eventually returns an error code
subroutine psb_dscall(m, n, parts, icontxt, desc_a, info)
subroutine psb_cdall(m, n, parts, icontxt, desc_a, info)
use psb_error_mod
use psb_descriptor_type
use psb_realloc_mod
@ -67,11 +67,11 @@ subroutine psb_dscall(m, n, parts, icontxt, desc_a, info)
if(psb_get_errstatus().ne.0) return
info=0
err=0
name = 'psb_dscall'
name = 'psb_cdall'
call psb_erractionsave(err_act)
call blacs_gridinfo(icontxt, nprow, npcol, myrow, mycol)
if (debug) write(*,*) 'psb_dscall: ',nprow,npcol,myrow,mycol
if (debug) write(*,*) 'psb_cdall: ',nprow,npcol,myrow,mycol
! ....verify blacs grid correctness..
if (npcol /= 1) then
info = 2030
@ -100,7 +100,7 @@ subroutine psb_dscall(m, n, parts, icontxt, desc_a, info)
endif
if (debug) write(*,*) 'psb_dscall: doing global checks'
if (debug) write(*,*) 'psb_cdall: doing global checks'
!global check on m and n parameters
if (myrow.eq.psb_root_) then
exch(1)=m
@ -137,7 +137,7 @@ subroutine psb_dscall(m, n, parts, icontxt, desc_a, info)
endif
if (debug) write(*,*) 'PSB_DSCALL: starting main loop' ,info
if (debug) write(*,*) 'PSB_CDALL: starting main loop' ,info
counter = 0
itmpov = 0
temp_ovrlap(:) = -1
@ -210,10 +210,10 @@ subroutine psb_dscall(m, n, parts, icontxt, desc_a, info)
loc_row=counter
! check on parts function
if (debug) write(*,*) 'PSB_DSCALL: End main loop:' ,loc_row,itmpov,info
if (debug) write(*,*) 'PSB_CDALL: End main loop:' ,loc_row,itmpov,info
if (debug) write(*,*) 'PSB_DSCALL: error check:' ,err
if (debug) write(*,*) 'PSB_CDALL: error check:' ,err
l_ov_ix=0
l_ov_el=0
@ -231,7 +231,7 @@ subroutine psb_dscall(m, n, parts, icontxt, desc_a, info)
l_ov_ix = l_ov_ix+3
l_ov_el = l_ov_el+3
if (debug) write(*,*) 'PSB_DSCALL: Ov len',l_ov_ix,l_ov_el
if (debug) write(*,*) 'PSB_CDALL: Ov len',l_ov_ix,l_ov_el
allocate(ov_idx(l_ov_ix),ov_el(l_ov_el), stat=info)
if (info /= no_err) then
info=4010
@ -295,7 +295,7 @@ subroutine psb_dscall(m, n, parts, icontxt, desc_a, info)
enddo
nullify(desc_a%bnd_elem,desc_a%halo_index)
!!$ if (debug) write(*,*) 'PSB_DSCALL: Last bits in desc_a', loc_row,k
!!$ if (debug) write(*,*) 'PSB_CDALL: Last bits in desc_a', loc_row,k
! set fields in desc_a%MATRIX_DATA....
desc_a%matrix_data(psb_n_row_) = loc_row
desc_a%matrix_data(psb_n_col_) = loc_row
@ -328,4 +328,4 @@ subroutine psb_dscall(m, n, parts, icontxt, desc_a, info)
end if
return
end subroutine psb_dscall
end subroutine psb_cdall

@ -28,9 +28,9 @@
!!$ POSSIBILITY OF SUCH DAMAGE.
!!$
!!$
! File: psb_dscalv.f90
! File: psb_cdalv.f90
!
! Subroutine: psb_dscalv
! Subroutine: psb_cdalv
! Allocate descriptor
! and checks correctness of PARTS subroutine
!
@ -41,7 +41,7 @@
! desc_a - type(<psb_desc_type>). The communication descriptor.
! info - integer. Eventually returns an error code
! flag - integer. ???
subroutine psb_dscalv(m, v, icontxt, desc_a, info, flag)
subroutine psb_cdalv(m, v, icontxt, desc_a, info, flag)
use psb_descriptor_type
use psb_serial_mod
use psb_const_mod
@ -66,10 +66,10 @@ subroutine psb_dscalv(m, v, icontxt, desc_a, info, flag)
if(psb_get_errstatus().ne.0) return
info=0
err=0
name = 'psb_dscalv'
name = 'psb_cdalv'
call blacs_gridinfo(icontxt, nprow, npcol, myrow, mycol)
if (debug) write(*,*) 'psb_dscall: ',nprow,npcol,myrow,mycol
if (debug) write(*,*) 'psb_cdall: ',nprow,npcol,myrow,mycol
! ....verify blacs grid correctness..
if (npcol /= 1) then
info = 2030
@ -99,7 +99,7 @@ subroutine psb_dscalv(m, v, icontxt, desc_a, info, flag)
goto 9999
end if
if (debug) write(*,*) 'psb_dscall: doing global checks'
if (debug) write(*,*) 'psb_cdall: doing global checks'
!global check on m and n parameters
if (myrow.eq.psb_root_) then
exch(1)=m
@ -150,7 +150,7 @@ subroutine psb_dscalv(m, v, icontxt, desc_a, info, flag)
endif
if (debug) write(*,*) 'PSB_DSCALL: starting main loop' ,info
if (debug) write(*,*) 'PSB_CDALL: starting main loop' ,info
counter = 0
itmpov = 0
temp_ovrlap(:) = -1
@ -175,13 +175,13 @@ subroutine psb_dscalv(m, v, icontxt, desc_a, info, flag)
loc_row=counter
! check on parts function
if (debug) write(*,*) 'PSB_DSCALL: End main loop:' ,loc_row,itmpov,info
if (debug) write(*,*) 'PSB_CDALL: End main loop:' ,loc_row,itmpov,info
if (info /= 0) then
call psb_errpush(info,name,i_err=int_err)
goto 9999
end if
if (debug) write(*,*) 'PSB_DSCALL: error check:' ,err
if (debug) write(*,*) 'PSB_CDALL: error check:' ,err
l_ov_ix=0
l_ov_el=0
@ -199,7 +199,7 @@ subroutine psb_dscalv(m, v, icontxt, desc_a, info, flag)
l_ov_ix = l_ov_ix+3
l_ov_el = l_ov_el+3
if (debug) write(*,*) 'PSB_DSCALL: Ov len',l_ov_ix,l_ov_el
if (debug) write(*,*) 'PSB_CDALL: Ov len',l_ov_ix,l_ov_el
allocate(ov_idx(l_ov_ix),ov_el(l_ov_el), stat=info)
if (info /= 0) then
info=2025
@ -264,7 +264,7 @@ subroutine psb_dscalv(m, v, icontxt, desc_a, info, flag)
enddo
nullify(desc_a%bnd_elem,desc_a%halo_index)
!!$ if (debug) write(*,*) 'PSB_DSCALL: Last bits in desc_a', loc_row,k
!!$ if (debug) write(*,*) 'PSB_CDALL: Last bits in desc_a', loc_row,k
! set fields in desc_a%MATRIX_DATA....
desc_a%matrix_data(psb_n_row_) = loc_row
desc_a%matrix_data(psb_n_col_) = loc_row
@ -296,4 +296,4 @@ subroutine psb_dscalv(m, v, icontxt, desc_a, info, flag)
end if
return
end subroutine psb_dscalv
end subroutine psb_cdalv

@ -28,7 +28,7 @@
!!$ POSSIBILITY OF SUCH DAMAGE.
!!$
!!$
subroutine psb_dscdec(nloc, icontxt, desc_a, info)
subroutine psb_cddec(nloc, icontxt, desc_a, info)
! Purpose
! =======
@ -127,10 +127,10 @@ subroutine psb_dscdec(nloc, icontxt, desc_a, info)
if(psb_get_errstatus().ne.0) return
info=0
err=0
name = 'psb_dscdec'
name = 'psb_cddec'
call blacs_gridinfo(icontxt, nprow, npcol, me, mypcol)
if (debug) write(*,*) 'psb_dscall: ',nprow,npcol,me,mypcol
if (debug) write(*,*) 'psb_cdalll: ',nprow,npcol,me,mypcol
! ....verify blacs grid correctness..
if (npcol /= 1) then
info = 2030
@ -156,7 +156,7 @@ subroutine psb_dscdec(nloc, icontxt, desc_a, info)
goto 9999
end if
if (debug) write(*,*) 'psb_dscall: doing global checks'
if (debug) write(*,*) 'psb_cdall: doing global checks'
!global check on m and n parameters
allocate(nlv(0:nprow-1))
@ -230,4 +230,4 @@ subroutine psb_dscdec(nloc, icontxt, desc_a, info)
end if
return
end subroutine psb_dscdec
end subroutine psb_cddec

@ -28,15 +28,15 @@
!!$ POSSIBILITY OF SUCH DAMAGE.
!!$
!!$
! File: psb_dscfree.f90
! File: psb_cdfree.f90
!
! Subroutine: psb_dscfree
! Subroutine: psb_cdfree
! Frees a descriptor data structure.
!
! Parameters:
! desc_a - type(<psb_desc_type>). The communication descriptor to be freed.
! info - integer. Eventually returns an error code.
subroutine psb_dscfree(desc_a,info)
subroutine psb_cdfree(desc_a,info)
!...free descriptor structure...
use psb_descriptor_type
use psb_const_mod
@ -55,7 +55,7 @@ subroutine psb_dscfree(desc_a,info)
if(psb_get_errstatus().ne.0) return
info=0
call psb_erractionsave(err_act)
name = 'psb_dscfree'
name = 'psb_cdfree'
if (.not.associated(desc_a%matrix_data)) then
@ -181,4 +181,4 @@ subroutine psb_dscfree(desc_a,info)
end if
return
end subroutine psb_dscfree
end subroutine psb_cdfree

@ -28,7 +28,7 @@
!!$ POSSIBILITY OF SUCH DAMAGE.
!!$
!!$
subroutine psb_dscrep(m, icontxt, desc_a, info)
subroutine psb_cdrep(m, icontxt, desc_a, info)
! Purpose
! =======
@ -126,10 +126,10 @@ subroutine psb_dscrep(m, icontxt, desc_a, info)
if(psb_get_errstatus().ne.0) return
info=0
err=0
name = 'psb_dscrep'
name = 'psb_cdrep'
call blacs_gridinfo(icontxt, nprow, npcol, myrow, mycol)
if (debug) write(*,*) 'psb_dscall: ',nprow,npcol,myrow,mycol
if (debug) write(*,*) 'psb_cdrep: ',nprow,npcol,myrow,mycol
! ....verify blacs grid correctness..
if (npcol /= 1) then
info = 2030
@ -226,4 +226,4 @@ subroutine psb_dscrep(m, icontxt, desc_a, info)
end if
return
end subroutine psb_dscrep
end subroutine psb_cdrep

@ -24,7 +24,7 @@ read_mat.o: mmio.o
df_sample: $(DFOBJS)
$(F90LINK) $(LINKOPT) $(DFOBJS) -o df_sample\
$(PSBLAS_LIB) $(METIS_LIB) $(BLACS)
$(PSBLAS_LIB) $(METIS_LIB) $(BLACS) $(SLU) $(UMF) $(BLAS)
/bin/mv df_sample $(EXEDIR)

@ -350,7 +350,7 @@ program df_sample
call psb_free(x_col, desc_a,info)
call psb_spfree(a, desc_a,info)
call psb_precfree(pre,info)
call psb_dscfree(desc_a,info)
call psb_cdfree(desc_a,info)
9999 continue
if(info /= 0) then

@ -187,15 +187,15 @@ contains
&nrow, ncol, nnzero,nrhs
endif
if (newt) then
call psb_dscall(nrow,nrow,parts,icontxt,desc_a,info)
call psb_cdall(nrow,nrow,parts,icontxt,desc_a,info)
if(info/=0) then
info=4010
ch_err='psb_dscall'
ch_err='psb_cdall'
call psb_errpush(info,name,a_err=ch_err)
goto 9999
end if
else
call psb_dscall(nrow,nrow,parts,icontxt,desc_a,info)
call psb_cdall(nrow,nrow,parts,icontxt,desc_a,info)
if(info/=0) then
info=4010
ch_err='psb_psdscall'
@ -633,10 +633,10 @@ contains
goto 9999
endif
call psb_dscall(nrow,v,icontxt,desc_a,info)
call psb_cdall(nrow,v,icontxt,desc_a,info)
if(info/=0) then
info=4010
ch_err='psb_dscall'
ch_err='psb_cdall'
call psb_errpush(info,name,a_err=ch_err)
goto 9999
end if

@ -255,7 +255,7 @@ program pde90
call psb_free(b,desc_a,info)
call psb_free(x,desc_a,info)
call psb_spfree(a,desc_a,info)
call psb_dscfree(desc_a,info)
call psb_cdfree(desc_a,info)
if(info.ne.0) then
info=4010
ch_err='free routine'
@ -485,7 +485,7 @@ contains
nnz = ((n*9)/(nprow*npcol))
if(myprow.eq.psb_root_) write(0,'("Generating Matrix (size=",i0x,")...")')n
call psb_dscall(n,n,parts,icontxt,desc_a,info)
call psb_cdall(n,n,parts,icontxt,desc_a,info)
call psb_spalloc(a,desc_a,info,nnz=nnz)
! define rhs from boundary conditions; also build initial guess
call psb_alloc(n,b,desc_a,info)

Loading…
Cancel
Save