From d00d9cc517d5708b5155356914ec50f88e670393 Mon Sep 17 00:00:00 2001 From: Salvatore Filippone Date: Thu, 24 Jan 2008 16:09:53 +0000 Subject: [PATCH] In psblas/prec fixed some constant names. Some more to be fixed. --- prec/psb_dgprec_aply.f90 | 8 ++++---- prec/psb_dprecbld.f90 | 8 ++++---- prec/psb_dprecinit.f90 | 6 +++--- prec/psb_dprecset.f90 | 4 ++-- prec/psb_prec_type.f90 | 22 +++++++++++----------- prec/psb_zgprec_aply.f90 | 14 +++++++++----- prec/psb_zprecbld.f90 | 8 ++++---- prec/psb_zprecinit.f90 | 6 +++--- prec/psb_zprecset.f90 | 4 ++-- 9 files changed, 42 insertions(+), 38 deletions(-) diff --git a/prec/psb_dgprec_aply.f90 b/prec/psb_dgprec_aply.f90 index 30205867..f7d756a0 100644 --- a/prec/psb_dgprec_aply.f90 +++ b/prec/psb_dgprec_aply.f90 @@ -75,11 +75,11 @@ subroutine psb_dgprec_aply(alpha,prec,x,beta,y,desc_data,trans,work,info) select case(prec%iprcparm(p_type_)) - case(noprec_) + case(psb_noprec_) call psb_geaxpby(alpha,x,beta,y,desc_data,info) - case(diag_) + case(psb_diag_) if (size(work) >= size(x)) then ww => work @@ -103,9 +103,9 @@ subroutine psb_dgprec_aply(alpha,prec,x,beta,y,desc_data,trans,work,info) end if end if - case(bjac_) + case(psb_bjac_) - call psb_bjac_aply(alpha,prec,x,beta,y,desc_data,trans,work,info) + call psb_bjac_aply(alpha,prec,x,beta,y,desc_data,trans_,work,info) if(info /= 0) then info=4010 ch_err='psb_bjac_aply' diff --git a/prec/psb_dprecbld.f90 b/prec/psb_dprecbld.f90 index 9ba26257..461a4b0e 100644 --- a/prec/psb_dprecbld.f90 +++ b/prec/psb_dprecbld.f90 @@ -81,12 +81,12 @@ subroutine psb_dprecbld(a,desc_a,p,info,upd) ! call psb_check_def(p%iprcparm(p_type_),'base_prec',& - & diag_,is_legal_prec) + & psb_diag_,is_legal_prec) call psb_nullify_desc(p%desc_data) select case(p%iprcparm(p_type_)) - case (noprec_) + case (psb_noprec_) ! Do nothing. call psb_cdcpy(desc_a,p%desc_data,info) if(info /= 0) then @@ -96,7 +96,7 @@ subroutine psb_dprecbld(a,desc_a,p,info,upd) goto 9999 end if - case (diag_) + case (psb_diag_) call psb_diagsc_bld(a,desc_a,p,upd_,info) if(info /= 0) then @@ -106,7 +106,7 @@ subroutine psb_dprecbld(a,desc_a,p,info,upd) goto 9999 end if - case (bjac_) + case (psb_bjac_) call psb_check_def(p%iprcparm(f_type_),'fact',& & f_ilu_n_,is_legal_ml_fact) diff --git a/prec/psb_dprecinit.f90 b/prec/psb_dprecinit.f90 index 44934cd4..9b69fe8d 100644 --- a/prec/psb_dprecinit.f90 +++ b/prec/psb_dprecinit.f90 @@ -47,17 +47,17 @@ subroutine psb_dprecinit(p,ptype,info) select case(toupper(ptype(1:len_trim(ptype)))) case ('NONE','NOPREC') p%iprcparm(:) = 0 - p%iprcparm(p_type_) = noprec_ + p%iprcparm(p_type_) = psb_noprec_ p%iprcparm(f_type_) = f_none_ case ('DIAG') p%iprcparm(:) = 0 - p%iprcparm(p_type_) = diag_ + p%iprcparm(p_type_) = psb_diag_ p%iprcparm(f_type_) = f_none_ case ('BJAC') p%iprcparm(:) = 0 - p%iprcparm(p_type_) = bjac_ + p%iprcparm(p_type_) = psb_bjac_ p%iprcparm(f_type_) = f_ilu_n_ p%iprcparm(ilu_fill_in_) = 0 diff --git a/prec/psb_dprecset.f90 b/prec/psb_dprecset.f90 index e82f8fb8..0f501fd8 100644 --- a/prec/psb_dprecset.f90 +++ b/prec/psb_dprecset.f90 @@ -41,7 +41,7 @@ subroutine psb_dprecseti(p,what,val,info) select case(what) case (f_type_) - if (p%iprcparm(p_type_) /= bjac_) then + if (p%iprcparm(p_type_) /= psb_bjac_) then write(0,*) 'WHAT is invalid for current preconditioner ',p%iprcparm(p_type_),& & 'ignoring user specification' return @@ -49,7 +49,7 @@ subroutine psb_dprecseti(p,what,val,info) p%iprcparm(f_type_) = val case (ilu_fill_in_) - if ((p%iprcparm(p_type_) /= bjac_).or.(p%iprcparm(f_type_) /= f_ilu_n_)) then + if ((p%iprcparm(p_type_) /= psb_bjac_).or.(p%iprcparm(f_type_) /= f_ilu_n_)) then write(0,*) 'WHAT is invalid for current preconditioner ',p%iprcparm(p_type_),& & 'ignoring user specification' return diff --git a/prec/psb_prec_type.f90 b/prec/psb_prec_type.f90 index 5c95601a..9e27b221 100644 --- a/prec/psb_prec_type.f90 +++ b/prec/psb_prec_type.f90 @@ -39,8 +39,8 @@ module psb_prec_type use psb_base_mod, only : psb_dspmat_type, psb_zspmat_type, psb_desc_type,& & psb_sizeof - integer, parameter :: min_prec_=0, noprec_=0, diag_=1, bjac_=2,& - & max_prec_=2 + integer, parameter :: psb_min_prec_=0, psb_noprec_=0, psb_diag_=1, & + & psb_bjac_=2, psb_max_prec_=2 ! Entries in iprcparm: preconditioner type, factorization type, ! prolongation type, restriction type, renumbering algorithm, @@ -129,11 +129,11 @@ contains write(iout,*) 'Preconditioner description' select case(p%iprcparm(p_type_)) - case(noprec_) + case(psb_noprec_) write(iout,*) 'No preconditioning' - case(diag_) + case(psb_diag_) write(iout,*) 'Diagonal scaling' - case(bjac_) + case(psb_bjac_) write(iout,*) 'Block Jacobi with: ',& & fact_names(p%iprcparm(f_type_)) end select @@ -147,11 +147,11 @@ contains write(iout,*) 'Preconditioner description' select case(p%iprcparm(p_type_)) - case(noprec_) + case(psb_noprec_) write(iout,*) 'No preconditioning' - case(diag_) + case(psb_diag_) write(iout,*) 'Diagonal scaling' - case(bjac_) + case(psb_bjac_) write(iout,*) 'Block Jacobi with: ',& & fact_names(p%iprcparm(f_type_)) end select @@ -367,11 +367,11 @@ contains character(len=10) :: pr_to_str select case(iprec) - case(noprec_) + case(psb_noprec_) pr_to_str='NOPREC' - case(diag_) + case(psb_diag_) pr_to_str='DIAG' - case(bjac_) + case(psb_bjac_) pr_to_str='BJAC' case default pr_to_str='???' diff --git a/prec/psb_zgprec_aply.f90 b/prec/psb_zgprec_aply.f90 index dc2c4a01..08ce7d6f 100644 --- a/prec/psb_zgprec_aply.f90 +++ b/prec/psb_zgprec_aply.f90 @@ -76,11 +76,11 @@ subroutine psb_zgprec_aply(alpha,prec,x,beta,y,desc_data,trans,work,info) select case(prec%iprcparm(p_type_)) - case(noprec_) + case(psb_noprec_) call psb_geaxpby(alpha,x,beta,y,desc_data,info) - case(diag_) + case(psb_diag_) if (size(work) >= size(x)) then ww => work @@ -93,7 +93,11 @@ subroutine psb_zgprec_aply(alpha,prec,x,beta,y,desc_data,trans,work,info) end if n_row=desc_data%matrix_data(psb_n_row_) - ww(1:n_row) = x(1:n_row)*prec%d(1:n_row) + if (trans_=='C') then + ww(1:n_row) = x(1:n_row)*conjg(prec%d(1:n_row)) + else + ww(1:n_row) = x(1:n_row)*prec%d(1:n_row) + endif call psb_geaxpby(alpha,ww,beta,y,desc_data,info) if (size(work) < size(x)) then @@ -104,9 +108,9 @@ subroutine psb_zgprec_aply(alpha,prec,x,beta,y,desc_data,trans,work,info) end if end if - case(bjac_) + case(psb_bjac_) - call psb_bjac_aply(alpha,prec,x,beta,y,desc_data,trans,work,info) + call psb_bjac_aply(alpha,prec,x,beta,y,desc_data,trans_,work,info) if(info /= 0) then info=4010 ch_err='psb_bjac_aply' diff --git a/prec/psb_zprecbld.f90 b/prec/psb_zprecbld.f90 index 06313e39..e6e47a34 100644 --- a/prec/psb_zprecbld.f90 +++ b/prec/psb_zprecbld.f90 @@ -82,12 +82,12 @@ subroutine psb_zprecbld(a,desc_a,p,info,upd) ! call psb_check_def(p%iprcparm(p_type_),'base_prec',& - & diag_,is_legal_prec) + & psb_diag_,is_legal_prec) call psb_nullify_desc(p%desc_data) select case(p%iprcparm(p_type_)) - case (noprec_) + case (psb_noprec_) ! Do nothing. call psb_cdcpy(desc_a,p%desc_data,info) if(info /= 0) then @@ -97,7 +97,7 @@ subroutine psb_zprecbld(a,desc_a,p,info,upd) goto 9999 end if - case (diag_) + case (psb_diag_) call psb_diagsc_bld(a,desc_a,p,upd_,info) if(info /= 0) then @@ -107,7 +107,7 @@ subroutine psb_zprecbld(a,desc_a,p,info,upd) goto 9999 end if - case (bjac_) + case (psb_bjac_) call psb_check_def(p%iprcparm(f_type_),'fact',& & f_ilu_n_,is_legal_ml_fact) diff --git a/prec/psb_zprecinit.f90 b/prec/psb_zprecinit.f90 index d6ff3431..53dfd646 100644 --- a/prec/psb_zprecinit.f90 +++ b/prec/psb_zprecinit.f90 @@ -48,17 +48,17 @@ subroutine psb_zprecinit(p,ptype,info) select case(toupper(ptype(1:len_trim(ptype)))) case ('NONE','NOPREC') p%iprcparm(:) = 0 - p%iprcparm(p_type_) = noprec_ + p%iprcparm(p_type_) = psb_noprec_ p%iprcparm(f_type_) = f_none_ case ('DIAG') p%iprcparm(:) = 0 - p%iprcparm(p_type_) = diag_ + p%iprcparm(p_type_) = psb_diag_ p%iprcparm(f_type_) = f_none_ case ('BJAC') p%iprcparm(:) = 0 - p%iprcparm(p_type_) = bjac_ + p%iprcparm(p_type_) = psb_bjac_ p%iprcparm(f_type_) = f_ilu_n_ p%iprcparm(ilu_fill_in_) = 0 diff --git a/prec/psb_zprecset.f90 b/prec/psb_zprecset.f90 index 551ec04c..b13f8b8e 100644 --- a/prec/psb_zprecset.f90 +++ b/prec/psb_zprecset.f90 @@ -41,7 +41,7 @@ subroutine psb_zprecseti(p,what,val,info) select case(what) case (f_type_) - if (p%iprcparm(p_type_) /= bjac_) then + if (p%iprcparm(p_type_) /= psb_bjac_) then write(0,*) 'WHAT is invalid for current preconditioner ',p%iprcparm(p_type_),& & 'ignoring user specification' return @@ -49,7 +49,7 @@ subroutine psb_zprecseti(p,what,val,info) p%iprcparm(f_type_) = val case (ilu_fill_in_) - if ((p%iprcparm(p_type_) /= bjac_).or.(p%iprcparm(f_type_) /= f_ilu_n_)) then + if ((p%iprcparm(p_type_) /= psb_bjac_).or.(p%iprcparm(f_type_) /= f_ilu_n_)) then write(0,*) 'WHAT is invalid for current preconditioner ',p%iprcparm(p_type_),& & 'ignoring user specification' return