*** empty log message ***

psblas3-type-indexed
Alfredo Buttari 20 years ago
parent 17ae68eacb
commit 1b407fa3ec

@ -128,7 +128,6 @@ subroutine psi_dswapdatam(flag,n,beta,y,desc_a,work,info)
sndbuf => work(1:idxs)
rcvbuf => work(idxs+1:idxs+idxr)
else
write(0,'(i2," allocating",3(i6,2x))')myrow,idxs,idxr,size(work)
allocate(sndbuf(idxs),rcvbuf(idxr), stat=info)
if(info.ne.0) then
call psb_errpush(4000,name)
@ -300,9 +299,7 @@ subroutine psi_dswapdatam(flag,n,beta,y,desc_a,work,info)
end do
do i=1, totxch
write(0,'(i2," waiting")')myrow
call mpi_waitany(nprow,rvhd,ixrec,p2pstat,iret)
write(0,'(i2," done")')myrow
if(iret.ne.mpi_success) then
int_err(1) = iret
info=400
@ -498,6 +495,7 @@ subroutine psi_dswapdatav(flag,beta,y,desc_a,work,info)
call blacs_get(icontxt,10,icomm)
allocate(sdsz(0:nprow-1), rvsz(0:nprow-1), bsdidx(0:nprow-1),&
& brvidx(0:nprow-1), rvhd(0:nprow-1), prcid(0:nprow-1),&
& ptp(0:nprow-1), stat=info)
@ -728,9 +726,9 @@ subroutine psi_dswapdatav(flag,beta,y,desc_a,work,info)
nesd = h_idx(point_to_proc+nerv+psb_n_elem_send_)
idx_pt = point_to_proc+psb_elem_recv_
snd_pt = bsdidx(proc_to_comm)
rcv_pt = brvidx(proc_to_comm)
call psi_sct(nerv,h_idx(idx_pt:idx_pt+nerv-1),&
& sndbuf(snd_pt:snd_pt+nesd-1),beta,y)
& rcvbuf(rcv_pt:rcv_pt+nerv-1),beta,y)
else
int_err(1) = ixrec
info=400

@ -77,7 +77,7 @@ Subroutine psb_dcgstab(a,prec,b,x,eps,desc_a,info,&
Logical, Parameter :: exchange=.True., noexchange=.False., debug1 = .False.
Integer, Parameter :: ione=1
Integer, Parameter :: irmax = 8
Integer :: itx, i, isvch, ich, icontxt, err_act, int_err(5)
Integer :: itx, i, isvch, ich, icontxt, err_act, int_err(5),ii
Integer :: listop
Logical :: do_renum_left
Real(Kind(1.d0)), Parameter :: one=1.d0, zero=0.d0, epstol=1.d-35
@ -175,9 +175,6 @@ Subroutine psb_dcgstab(a,prec,b,x,eps,desc_a,info,&
Else If (listop == 2) Then
bn2 = psb_nrm2(b,desc_a,info)
Endif
call blacs_barrier(icontxt,'All') ! to be removed
write(0,'(i2," ani bni bn2 ",3(f10.2,2x))')myrow,ani,bni,bn2
call blacs_barrier(icontxt,'All') ! to be removed
if (info /= 0) Then
info=4011
call psb_errpush(info,name)
@ -190,17 +187,9 @@ Subroutine psb_dcgstab(a,prec,b,x,eps,desc_a,info,&
!!$
If (itx >= itmax) Exit restart
it = 0
write(0,'(i2," b ",10(f10.2,2x))')myrow,b(1:10)
write(0,'(i2," r ",10(f10.2,2x))')myrow,r(1:10)
Call psb_axpby(one,b,zero,r,desc_a,info)
write(0,'(i2," b ",10(f10.2,2x))')myrow,b(1:10)
write(0,'(i2," r ",10(f10.2,2x))')myrow,r(1:10)
Call psb_spmm(-one,a,x,one,r,desc_a,info,work=aux)
write(0,'(i2," x ",10(f10.2,2x))')myrow,x(1:10)
write(0,'(i2," r ",10(f10.2,2x))')myrow,r(1:10)
Call psb_axpby(one,r,zero,q,desc_a,info)
write(0,'(i2," q ",10(f10.2,2x))')myrow,q(1:10)
write(0,'(i2," r ",10(f10.2,2x))')myrow,r(1:10)
if (info /= 0) Then
info=4011
call psb_errpush(info,name)
@ -220,9 +209,6 @@ Subroutine psb_dcgstab(a,prec,b,x,eps,desc_a,info,&
Else If (listop == 2) Then
rni = psb_nrm2(r,desc_a,info)
Endif
call blacs_barrier(icontxt,'All') ! to be removed
write(0,'(i2," rni xni ",2(f10.2,2x))')myrow,rni,xni
call blacs_barrier(icontxt,'All') ! to be removed
if (info /= 0) Then
info=4011
call psb_errpush(info,name)
@ -269,7 +255,9 @@ Subroutine psb_dcgstab(a,prec,b,x,eps,desc_a,info,&
If (debug) Write(*,*) 'Iteration: ',itx
rho_old = rho
rho = psb_dot(q,r,desc_a,info)
write(0,'(i2," rho old ",2(f10.2,2x))')myrow,rho,rho_old
!!$ call blacs_barrier(icontxt,'All') ! to be removed
!!$ write(0,'(i2," rho old ",2(f,2x))')myrow,rho,rho_old
!!$ call blacs_barrier(icontxt,'All') ! to be removed
If (rho==zero) Then
If (debug) Write(0,*) 'Bi-CGSTAB Itxation breakdown R',rho
Exit iteration
@ -296,13 +284,30 @@ Subroutine psb_dcgstab(a,prec,b,x,eps,desc_a,info,&
alpha = rho/sigma
Call psb_axpby(one,r,zero,s,desc_a,info)
if(info.ne.0) then
call psb_errpush(4010,name,a_err='psb_axpby')
goto 9999
end if
Call psb_axpby(-alpha,v,one,s,desc_a,info)
if(info.ne.0) then
call psb_errpush(4010,name,a_err='psb_axpby')
goto 9999
end if
Call psb_prcaply(prec,s,z,desc_a,info,work=aux)
if(info.ne.0) then
call psb_errpush(4010,name,a_err='psb_prcaply')
goto 9999
end if
Call psb_spmm(one,a,z,zero,t,desc_a,info,&
& work=aux)
if(info.ne.0) then
call psb_errpush(4010,name,a_err='psb_spmm')
goto 9999
end if
sigma = psb_dot(t,t,desc_a,info)
If (sigma==zero) Then
If (debug) Write(0,*) 'BI-CGSTAB ITERATION BREAKDOWN S2', sigma

@ -22,6 +22,7 @@
subroutine psb_daxpby(alpha, x, beta,y,desc_a,info, n, jx, jy)
use psb_descriptor_type
use psb_check_mod
use psb_const_mod
use psb_error_mod
implicit none
@ -147,6 +148,7 @@ end subroutine psb_daxpby
!
subroutine psb_daxpbyv(alpha, x, beta,y,desc_a,info)
use psb_descriptor_type
use psb_const_mod
use psb_check_mod
use psb_error_mod
implicit none
@ -200,8 +202,6 @@ subroutine psb_daxpbyv(alpha, x, beta,y,desc_a,info)
call psb_errpush(info,name)
end if
write(0,'(i2," before daxpby",2(i6,2x),2(f10.2,2x))')myrow,desc_a%matrix_data(psb_n_row_),&
& desc_a%matrix_data(psb_n_col_),alpha,beta
if(desc_a%matrix_data(psb_n_row_).gt.0) then
call daxpby(desc_a%matrix_data(psb_n_col_),ione,&
& alpha,x,size(x),beta,&

@ -368,6 +368,7 @@ subroutine psb_dspmv(alpha,a,x,beta,y,desc_a,info,&
use psb_serial_mod
use psb_descriptor_type
use psb_comm_mod
use psb_const_mod
use psi_mod
use psb_check_mod
use psb_error_mod
@ -526,7 +527,6 @@ subroutine psb_dspmv(alpha,a,x,beta,y,desc_a,info,&
& a%descra,a%aspk,a%ia1,a%ia2,a%infoa,a%pl,&
& x(iix),lldx,beta,y(iiy),lldy,&
& iwork,liwork,info)
if(info.ne.0) then
info = 4011
call psb_errpush(info,name)

@ -175,7 +175,7 @@ C .. Local Scalars ..
INTEGER LWORKM, LWORKB, LWORKC, LWORKS, P, ERR_ACT
LOGICAL LP, RP
C .. Local Array..
INTEGER INT_VAL(5)
INTEGER INT_VAL(5),iunit
CHARACTER*20 NAME
DOUBLE PRECISION REAL_VAL(5)
CHARACTER*30 STRINGS(2)

@ -4,7 +4,7 @@ NONE Preconditioner ILU DIAGSC NONE
CSR A Storage format CSR COO JAD
20 Domain size (acutal sistem is this**3)
1 Stopping criterion
2 MAXIT
80 MAXIT
00 ITRACE
02 ML

@ -217,7 +217,7 @@ program pde90
write(6,*) 'time per iteration : ',t2/iter
write(6,*) 'number of iterations : ',iter
write(6,*) 'error on exit : ',err
write(6,*) 'info on exit : ',ierr
write(6,*) 'info on exit : ',info
end if
!

Loading…
Cancel
Save