Adjusted test programs to ALLOCATABLEs.

psblas3-type-indexed
Salvatore Filippone 18 years ago
parent 656eba0524
commit f6655e698d

@ -339,9 +339,9 @@ Contains
function psb_isize1d(vin)
integer :: psb_isize1d
integer, pointer :: vin(:)
integer, allocatable, intent(in) :: vin(:)
if (.not.associated(vin)) then
if (.not.allocated(vin)) then
psb_isize1d = 0
else
psb_isize1d = size(vin)
@ -349,10 +349,10 @@ Contains
end function psb_isize1d
function psb_isize2d(vin,dim)
integer :: psb_isize2d
integer, pointer :: vin(:,:)
integer, allocatable, intent(in) :: vin(:,:)
integer, optional :: dim
if (.not.associated(vin)) then
if (.not.allocated(vin)) then
psb_isize2d = 0
else
if (present(dim)) then
@ -365,9 +365,9 @@ Contains
function psb_dsize1d(vin)
integer :: psb_dsize1d
real(kind(1.d0)), pointer :: vin(:)
real(kind(1.d0)), allocatable, intent(in) :: vin(:)
if (.not.associated(vin)) then
if (.not.allocated(vin)) then
psb_dsize1d = 0
else
psb_dsize1d = size(vin)
@ -375,10 +375,10 @@ Contains
end function psb_dsize1d
function psb_dsize2d(vin,dim)
integer :: psb_dsize2d
real(kind(1.d0)), pointer :: vin(:,:)
real(kind(1.d0)), allocatable, intent(in) :: vin(:,:)
integer, optional :: dim
if (.not.associated(vin)) then
if (.not.allocated(vin)) then
psb_dsize2d = 0
else
if (present(dim)) then
@ -392,9 +392,9 @@ Contains
function psb_zsize1d(vin)
integer :: psb_zsize1d
complex(kind(1.d0)), pointer :: vin(:)
complex(kind(1.d0)), allocatable, intent(in) :: vin(:)
if (.not.associated(vin)) then
if (.not.allocated(vin)) then
psb_zsize1d = 0
else
psb_zsize1d = size(vin)
@ -402,10 +402,10 @@ Contains
end function psb_zsize1d
function psb_zsize2d(vin,dim)
integer :: psb_zsize2d
complex(kind(1.d0)), pointer :: vin(:,:)
complex(kind(1.d0)), allocatable, intent(in) :: vin(:,:)
integer, optional :: dim
if (.not.associated(vin)) then
if (.not.allocated(vin)) then
psb_zsize2d = 0
else
if (present(dim)) then

@ -68,9 +68,10 @@ program df_sample
integer :: igsmth, matop, novr
! dense matrices
real(kind(1.d0)), pointer :: aux_b(:,:), d(:)
real(kind(1.d0)), pointer, save :: b_col(:), x_col(:), r_col(:), &
& b_col_glob(:), x_col_glob(:), r_col_glob(:)
real(kind(1.d0)), allocatable, target :: aux_b(:,:), d(:)
real(kind(1.d0)), allocatable , save :: b_col(:), x_col(:), r_col(:), &
& x_col_glob(:), r_col_glob(:)
real(kind(1.d0)), pointer :: b_col_glob(:)
! communications data structure
type(psb_desc_type):: desc_a
@ -125,7 +126,6 @@ program df_sample
nrhs = 1
if (amroot) then
nullify(aux_b)
call readmat(mtrx_file, aux_a, ictxt)
m_problem = aux_a%m
@ -136,14 +136,14 @@ program df_sample
call read_rhs(rhs_file,aux_b,ictxt)
end if
if (associated(aux_b).and.size(aux_b,1)==m_problem) then
if (psb_size(aux_b,dim=1)==m_problem) then
! if any rhs were present, broadcast the first one
write(0,'("Ok, got an rhs ")')
b_col_glob =>aux_b(:,1)
else
write(*,'("Generating an rhs...")')
write(*,'(" ")')
allocate(aux_b(m_problem,1), stat=ircode)
call psb_realloc(m_problem,1,aux_b,ircode)
if (ircode /= 0) then
call psb_errpush(4000,name)
goto 9999
@ -157,7 +157,7 @@ program df_sample
call psb_bcast(ictxt,b_col_glob(1:m_problem))
else
call psb_bcast(ictxt,m_problem)
allocate(aux_b(m_problem,1), stat=ircode)
call psb_realloc(m_problem,1,aux_b,ircode)
if (ircode /= 0) then
call psb_errpush(4000,name)
goto 9999
@ -274,8 +274,8 @@ program df_sample
call psb_barrier(ictxt)
t2 = psb_wtime() - t1
call psb_amx(ictxt,t2)
call psb_geaxpby(1.d0,b_col,0.d0,r_col,desc_a,info)
call psb_spmm(-1.d0,a,x_col,1.d0,r_col,desc_a,info)
call psb_geaxpby(done,b_col,dzero,r_col,desc_a,info)
call psb_spmm(-done,a,x_col,done,r_col,desc_a,info)
call psb_genrm2s(resmx,r_col,desc_a,info)
call psb_geamaxs(resmxp,r_col,desc_a,info)

@ -82,11 +82,11 @@ contains
! on entry: fresh variable.
! on exit : the updated array descriptor
!
! real(kind(1.d0)), pointer, optional :: b_glob(:)
! real(kind(1.d0)), optional :: b_glob(:)
! on entry: this contains right hand side.
! on exit :
!
! real(kind(1.d0)), pointer, optional :: b(:)
! real(kind(1.d0)), allocatable, optional :: b(:)
! on entry: fresh variable.
! on exit : this will contain the local right hand side.
!
@ -99,10 +99,10 @@ contains
! parameters
type(psb_dspmat_type) :: a_glob
real(kind(1.d0)), pointer :: b_glob(:)
real(kind(1.d0)) :: b_glob(:)
integer :: ictxt
type(psb_dspmat_type) :: a
real(kind(1.d0)), pointer :: b(:)
real(kind(1.d0)), allocatable :: b(:)
type (psb_desc_type) :: desc_a
integer, intent(out) :: info
integer, optional :: inroot
@ -123,7 +123,7 @@ contains
integer :: ircode, length_row, i_count, j_count,&
& k_count, blockdim, root, liwork, nrow, ncol, nnzero, nrhs,&
& i,j,k, ll, nz, isize, iproc, nnr, err, err_act, int_err(5)
integer, pointer :: iwork(:)
integer, allocatable :: iwork(:)
character :: afmt*5, atyp*5
integer, allocatable :: irow(:),icol(:)
real(kind(1.d0)), allocatable :: val(:)
@ -530,11 +530,11 @@ contains
! on entry: fresh variable.
! on exit : the updated array descriptor
!
! real(kind(1.d0)), pointer, optional :: b_glob(:)
! real(kind(1.d0)), optional :: b_glob(:)
! on entry: this contains right hand side.
! on exit :
!
! real(kind(1.d0)), pointer, optional :: b(:)
! real(kind(1.d0)), allocatable, optional :: b(:)
! on entry: fresh variable.
! on exit : this will contain the local right hand side.
!
@ -545,10 +545,10 @@ contains
use psb_sparse_mod
implicit none ! parameters
type(psb_dspmat_type) :: a_glob
real(kind(1.d0)), pointer :: b_glob(:)
real(kind(1.d0)) :: b_glob(:)
integer :: ictxt, v(:)
type(psb_dspmat_type) :: a
real(kind(1.d0)), pointer :: b(:)
real(kind(1.d0)), allocatable :: b(:)
type (psb_desc_type) :: desc_a
integer, intent(out) :: info
integer, optional :: inroot
@ -558,7 +558,7 @@ contains
integer :: ircode, length_row, i_count, j_count,&
& k_count, blockdim, root, liwork, nrow, ncol, nnzero, nrhs,&
& i,j,k, ll, nz, isize, iproc, nnr, err, err_act, int_err(5)
integer, pointer :: iwork(:)
integer, allocatable :: iwork(:)
character :: afmt*5, atyp*5
integer, allocatable :: irow(:),icol(:)
real(kind(1.d0)), allocatable :: val(:)
@ -857,11 +857,11 @@ contains
! on entry: fresh variable.
! on exit : the updated array descriptor
!
! real(kind(1.d0)), pointer, optional :: b_glob(:)
! real(kind(1.d0)), optional :: b_glob(:)
! on entry: this contains right hand side.
! on exit :
!
! real(kind(1.d0)), pointer, optional :: b(:)
! real(kind(1.d0)), allocatable, optional :: b(:)
! on entry: fresh variable.
! on exit : this will contain the local right hand side.
!
@ -874,10 +874,10 @@ contains
! parameters
type(psb_zspmat_type) :: a_glob
complex(kind(1.d0)), pointer :: b_glob(:)
complex(kind(1.d0)) :: b_glob(:)
integer :: ictxt
type(psb_zspmat_type) :: a
complex(kind(1.d0)), pointer :: b(:)
complex(kind(1.d0)), allocatable :: b(:)
type (psb_desc_type) :: desc_a
integer, intent(out) :: info
integer, optional :: inroot
@ -898,7 +898,7 @@ contains
integer :: ircode, length_row, i_count, j_count,&
& k_count, blockdim, root, liwork, nrow, ncol, nnzero, nrhs,&
& i,j,k, ll, nz, isize, iproc, nnr, err, err_act, int_err(5)
integer, pointer :: iwork(:)
integer, allocatable :: iwork(:)
character :: afmt*5, atyp*5
integer, allocatable :: irow(:),icol(:)
complex(kind(1.d0)), allocatable :: val(:)
@ -1302,11 +1302,11 @@ contains
! on entry: fresh variable.
! on exit : the updated array descriptor
!
! real(kind(1.d0)), pointer, optional :: b_glob(:)
! real(kind(1.d0)), optional :: b_glob(:)
! on entry: this contains right hand side.
! on exit :
!
! real(kind(1.d0)), pointer, optional :: b(:)
! real(kind(1.d0)), allocatable, optional :: b(:)
! on entry: fresh variable.
! on exit : this will contain the local right hand side.
!
@ -1317,10 +1317,10 @@ contains
use psb_sparse_mod
implicit none ! parameters
type(psb_zspmat_type) :: a_glob
complex(kind(1.d0)), pointer :: b_glob(:)
complex(kind(1.d0)) :: b_glob(:)
integer :: ictxt, v(:)
type(psb_zspmat_type) :: a
complex(kind(1.d0)), pointer :: b(:)
complex(kind(1.d0)), allocatable :: b(:)
type(psb_desc_type) :: desc_a
integer, intent(out) :: info
integer, optional :: inroot
@ -1330,7 +1330,7 @@ contains
integer :: ircode, length_row, i_count, j_count,&
& k_count, blockdim, root, liwork, nrow, ncol, nnzero, nrhs,&
& i,j,k, ll, nz, isize, iproc, nnr, err, err_act, int_err(5)
integer, pointer :: iwork(:)
integer, allocatable :: iwork(:)
character :: afmt*5, atyp*5
integer, allocatable :: irow(:),icol(:)
complex(kind(1.d0)), allocatable :: val(:)

@ -57,12 +57,6 @@
! On Entry: Index of root processor (default: 0)
! On Exit : unchanged.
!
! Real(Kind(1.D0)), Pointer, Optional :: indwork(:)
! On Entry/Exit: Double Precision Work Area.
!
! Integer, Pointer, Optional :: iniwork()
! On Entry/Exit: Integer Work Area.
!
module read_mat
interface readmat
module procedure dreadmat, zreadmat
@ -104,20 +98,17 @@ contains
end subroutine dreadmat
subroutine dread_rhs (filename, b, ictxt, inroot,&
& indwork, iniwork)
subroutine dread_rhs (filename, b, ictxt, inroot)
use psb_sparse_mod
implicit none
integer :: ictxt
character :: filename*(*)
integer, optional :: inroot
real(kind(1.0d0)), pointer, optional :: indwork(:)
integer, pointer, optional :: iniwork(:) ! local variables
integer, parameter :: infile = 2
integer :: nrow, ncol, i,root, nprow, npcol, myprow, mypcol, ircode, j
character :: mmheader*15, fmt*15, object*10, type*10, sym*15,&
& line*1024
real(kind(1.0d0)), pointer :: b(:,:)
real(kind(1.0d0)), allocatable :: b(:,:)
if (present(inroot)) then
root = inroot
else
@ -195,20 +186,17 @@ contains
end subroutine zreadmat
subroutine zread_rhs (filename, b, ictxt, inroot,&
& indwork, iniwork)
subroutine zread_rhs (filename, b, ictxt, inroot)
use psb_sparse_mod
implicit none
integer :: ictxt
character :: filename*(*)
integer, optional :: inroot
complex(kind(1.0d0)), pointer, optional :: indwork(:)
integer, pointer, optional :: iniwork(:) ! local variables
integer, parameter :: infile = 2
integer :: nrow, ncol, i,root, nprow, npcol, myprow, mypcol, ircode, j
character :: mmheader*15, fmt*15, object*10, type*10, sym*15,&
& line*1024
complex(kind(1.0d0)), pointer :: b(:,:)
complex(kind(1.0d0)), allocatable :: b(:,:)
if (present(inroot)) then
root = inroot
else

@ -68,9 +68,10 @@ program zf_sample
integer :: igsmth, matop, novr
! dense matrices
complex(kind(1.d0)), pointer :: aux_b(:,:), d(:)
complex(kind(1.d0)), pointer, save :: b_col(:), x_col(:), r_col(:), &
& b_col_glob(:), x_col_glob(:), r_col_glob(:)
complex(kind(1.d0)), allocatable, target :: aux_b(:,:), d(:)
complex(kind(1.d0)), allocatable , save :: b_col(:), x_col(:), r_col(:), &
& x_col_glob(:), r_col_glob(:)
complex(kind(1.d0)), pointer :: b_col_glob(:)
! communications data structure
type(psb_desc_type):: desc_a
@ -126,7 +127,6 @@ program zf_sample
nrhs = 1
if (amroot) then
nullify(aux_b)
call readmat(mtrx_file, aux_a, ictxt)
m_problem = aux_a%m
@ -137,14 +137,14 @@ program zf_sample
call read_rhs(rhs_file,aux_b,ictxt)
end if
if (associated(aux_b).and.size(aux_b,1)==m_problem) then
if (psb_size(aux_b,dim=1)==m_problem) then
! if any rhs were present, broadcast the first one
write(0,'("Ok, got an rhs ")')
b_col_glob =>aux_b(:,1)
else
write(*,'("Generating an rhs...")')
write(*,'(" ")')
allocate(aux_b(m_problem,1), stat=ircode)
call psb_realloc(m_problem,1,aux_b,ircode)
if (ircode /= 0) then
call psb_errpush(4000,name)
goto 9999
@ -158,7 +158,7 @@ program zf_sample
call psb_bcast(ictxt,b_col_glob(1:m_problem))
else
call psb_bcast(ictxt,m_problem)
allocate(aux_b(m_problem,1), stat=ircode)
call psb_realloc(m_problem,1,aux_b,ircode)
if (ircode /= 0) then
call psb_errpush(4000,name)
goto 9999
@ -192,10 +192,9 @@ program zf_sample
if (amroot) then
write(*,'("Partition type: graph")')
write(*,'(" ")')
write(0,'("Build type: graph")')
! write(0,'("Build type: graph")')
call build_grppart(aux_a%m,aux_a%fida,aux_a%ia1,aux_a%ia2,np)
endif
write(0,'("Done graph build")')
call psb_barrier(ictxt)
call distr_grppart(root,ictxt)
call getv_grppart(ivg)
@ -284,7 +283,7 @@ program zf_sample
if (amroot) then
call psb_prec_descr(6,pre)
write(*,'("Matrix: ",a)')mtrx_file
write(*,'("Computed solution on ",i4," processors")')np
write(*,'("Computed solution on ",i8," processors")')np
write(*,'("Iterations to convergence: ",i6)')iter
write(*,'("Error indicator on exit: ",f7.2)')err
write(*,'("Time to buil prec. : ",es10.4)')tprec

@ -97,8 +97,8 @@ program pde90
! descriptor
type(psb_desc_type) :: desc_a, desc_a_out
! dense matrices
real(kind(1.d0)), pointer :: b(:), x(:), d(:),ld(:)
integer, pointer :: work(:)
real(kind(1.d0)), allocatable :: b(:), x(:), d(:),ld(:)
integer, allocatable :: work(:)
! blacs parameters
integer :: ictxt, iam, np
@ -394,12 +394,12 @@ contains
use psb_sparse_mod
implicit none
integer :: idim
integer, parameter :: nbmax=10
real(kind(1.d0)),pointer :: b(:),t(:)
type(psb_desc_type) :: desc_a
integer :: ictxt, info
character :: afmt*5
integer :: idim
integer, parameter :: nbmax=10
real(kind(1.d0)), allocatable :: b(:),t(:)
type(psb_desc_type) :: desc_a
integer :: ictxt, info
character :: afmt*5
interface
! .....user passed subroutine.....
subroutine parts(global_indx,n,np,pv,nv)
@ -419,9 +419,6 @@ contains
integer, allocatable :: irow(:),icol(:)
real(kind(1.d0)), allocatable :: val(:)
integer, allocatable :: prv(:)
integer, pointer :: ierrv(:)
real(kind(1.d0)), pointer :: dwork(:)
integer,pointer :: iwork(:)
! deltah dimension of each grid cell
! deltat discretization time
real(kind(1.d0)) :: deltah

Loading…
Cancel
Save