added weighted 2-norm function

merge-paraggr-newops
Cirdans-Home 5 years ago
parent 63885c36c7
commit 1b214f3f00

@ -273,11 +273,21 @@ module psb_c_psblas_mod
integer(psb_ipk_), intent(out) :: info
logical, intent(in), optional :: global
end function psb_cnrm2_vect
function psb_cnrm2_weight_vect(x,w, desc_a, info,global) result(res)
import :: psb_desc_type, psb_spk_, psb_ipk_, &
& psb_c_vect_type, psb_cspmat_type
real(psb_spk_) :: res
type(psb_c_vect_type), intent (inout) :: x
type(psb_c_vect_type), intent (inout) :: w
type(psb_desc_type), intent (in) :: desc_a
integer(psb_ipk_), intent(out) :: info
logical, intent(in), optional :: global
end function psb_cnrm2_weight_vect
end interface
#if ! defined(HAVE_BUGGY_GENERICS)
interface psb_norm2
procedure psb_cnrm2, psb_cnrm2v, psb_cnrm2_vect
procedure psb_cnrm2, psb_cnrm2v, psb_cnrm2_vect, psb_cnrm2_weight_vect
end interface
#endif

@ -273,11 +273,21 @@ module psb_d_psblas_mod
integer(psb_ipk_), intent(out) :: info
logical, intent(in), optional :: global
end function psb_dnrm2_vect
function psb_dnrm2_weight_vect(x,w, desc_a, info,global) result(res)
import :: psb_desc_type, psb_dpk_, psb_ipk_, &
& psb_d_vect_type, psb_dspmat_type
real(psb_dpk_) :: res
type(psb_d_vect_type), intent (inout) :: x
type(psb_d_vect_type), intent (inout) :: w
type(psb_desc_type), intent (in) :: desc_a
integer(psb_ipk_), intent(out) :: info
logical, intent(in), optional :: global
end function psb_dnrm2_weight_vect
end interface
#if ! defined(HAVE_BUGGY_GENERICS)
interface psb_norm2
procedure psb_dnrm2, psb_dnrm2v, psb_dnrm2_vect
procedure psb_dnrm2, psb_dnrm2v, psb_dnrm2_vect, psb_dnrm2_weight_vect
end interface
#endif

@ -273,11 +273,21 @@ module psb_s_psblas_mod
integer(psb_ipk_), intent(out) :: info
logical, intent(in), optional :: global
end function psb_snrm2_vect
function psb_snrm2_weight_vect(x,w, desc_a, info,global) result(res)
import :: psb_desc_type, psb_spk_, psb_ipk_, &
& psb_s_vect_type, psb_sspmat_type
real(psb_spk_) :: res
type(psb_s_vect_type), intent (inout) :: x
type(psb_s_vect_type), intent (inout) :: w
type(psb_desc_type), intent (in) :: desc_a
integer(psb_ipk_), intent(out) :: info
logical, intent(in), optional :: global
end function psb_snrm2_weight_vect
end interface
#if ! defined(HAVE_BUGGY_GENERICS)
interface psb_norm2
procedure psb_snrm2, psb_snrm2v, psb_snrm2_vect
procedure psb_snrm2, psb_snrm2v, psb_snrm2_vect, psb_snrm2_weight_vect
end interface
#endif

@ -273,11 +273,21 @@ module psb_z_psblas_mod
integer(psb_ipk_), intent(out) :: info
logical, intent(in), optional :: global
end function psb_znrm2_vect
function psb_znrm2_weight_vect(x,w, desc_a, info,global) result(res)
import :: psb_desc_type, psb_dpk_, psb_ipk_, &
& psb_z_vect_type, psb_zspmat_type
real(psb_dpk_) :: res
type(psb_z_vect_type), intent (inout) :: x
type(psb_z_vect_type), intent (inout) :: w
type(psb_desc_type), intent (in) :: desc_a
integer(psb_ipk_), intent(out) :: info
logical, intent(in), optional :: global
end function psb_znrm2_weight_vect
end interface
#if ! defined(HAVE_BUGGY_GENERICS)
interface psb_norm2
procedure psb_znrm2, psb_znrm2v, psb_znrm2_vect
procedure psb_znrm2, psb_znrm2v, psb_znrm2_vect, psb_znrm2_weight_vect
end interface
#endif

@ -108,7 +108,9 @@ module psb_c_vect_mod
procedure, pass(x) :: absval1 => c_vect_absval1
procedure, pass(x) :: absval2 => c_vect_absval2
generic, public :: absval => absval1, absval2
procedure, pass(x) :: nrm2 => c_vect_nrm2
procedure, pass(x) :: nrm2std => c_vect_nrm2
procedure, pass(x) :: nrm2weight => c_vect_nrm2_weight
generic, public :: nrm2 => nrm2std, nrm2weight
procedure, pass(x) :: amax => c_vect_amax
procedure, pass(x) :: asum => c_vect_asum
end type psb_c_vect_type
@ -903,6 +905,23 @@ contains
end function c_vect_nrm2
function c_vect_nrm2_weight(n,x,w) result(res)
implicit none
class(psb_c_vect_type), intent(inout) :: x
class(psb_c_vect_type), intent(inout) :: w
integer(psb_ipk_), intent(in) :: n
real(psb_spk_) :: res
integer(psb_ipk_) :: info
if (allocated(x%v).and.allocated(w%v)) then
call w%v%mlt(x%v,info)
res = w%v%nrm2(n)
else
res = szero
end if
end function c_vect_nrm2_weight
function c_vect_amax(n,x) result(res)
implicit none
class(psb_c_vect_type), intent(inout) :: x

@ -108,7 +108,9 @@ module psb_d_vect_mod
procedure, pass(x) :: absval1 => d_vect_absval1
procedure, pass(x) :: absval2 => d_vect_absval2
generic, public :: absval => absval1, absval2
procedure, pass(x) :: nrm2 => d_vect_nrm2
procedure, pass(x) :: nrm2std => d_vect_nrm2
procedure, pass(x) :: nrm2weight => d_vect_nrm2_weight
generic, public :: nrm2 => nrm2std, nrm2weight
procedure, pass(x) :: amax => d_vect_amax
procedure, pass(x) :: asum => d_vect_asum
procedure, pass(z) :: cmp_a2 => d_vect_cmp_a2
@ -933,6 +935,23 @@ contains
end function d_vect_nrm2
function d_vect_nrm2_weight(n,x,w) result(res)
implicit none
class(psb_d_vect_type), intent(inout) :: x
class(psb_d_vect_type), intent(inout) :: w
integer(psb_ipk_), intent(in) :: n
real(psb_dpk_) :: res
integer(psb_ipk_) :: info
if (allocated(x%v).and.allocated(w%v)) then
call w%v%mlt(x%v,info)
res = w%v%nrm2(n)
else
res = dzero
end if
end function d_vect_nrm2_weight
function d_vect_amax(n,x) result(res)
implicit none
class(psb_d_vect_type), intent(inout) :: x

@ -108,7 +108,9 @@ module psb_s_vect_mod
procedure, pass(x) :: absval1 => s_vect_absval1
procedure, pass(x) :: absval2 => s_vect_absval2
generic, public :: absval => absval1, absval2
procedure, pass(x) :: nrm2 => s_vect_nrm2
procedure, pass(x) :: nrm2std => s_vect_nrm2
procedure, pass(x) :: nrm2weight => s_vect_nrm2_weight
generic, public :: nrm2 => nrm2std, nrm2weight
procedure, pass(x) :: amax => s_vect_amax
procedure, pass(x) :: asum => s_vect_asum
procedure, pass(z) :: cmp_a2 => s_vect_cmp_a2
@ -933,6 +935,23 @@ contains
end function s_vect_nrm2
function s_vect_nrm2_weight(n,x,w) result(res)
implicit none
class(psb_s_vect_type), intent(inout) :: x
class(psb_s_vect_type), intent(inout) :: w
integer(psb_ipk_), intent(in) :: n
real(psb_spk_) :: res
integer(psb_ipk_) :: info
if (allocated(x%v).and.allocated(w%v)) then
call w%v%mlt(x%v,info)
res = w%v%nrm2(n)
else
res = szero
end if
end function s_vect_nrm2_weight
function s_vect_amax(n,x) result(res)
implicit none
class(psb_s_vect_type), intent(inout) :: x

@ -108,7 +108,9 @@ module psb_z_vect_mod
procedure, pass(x) :: absval1 => z_vect_absval1
procedure, pass(x) :: absval2 => z_vect_absval2
generic, public :: absval => absval1, absval2
procedure, pass(x) :: nrm2 => z_vect_nrm2
procedure, pass(x) :: nrm2std => z_vect_nrm2
procedure, pass(x) :: nrm2weight => z_vect_nrm2_weight
generic, public :: nrm2 => nrm2std, nrm2weight
procedure, pass(x) :: amax => z_vect_amax
procedure, pass(x) :: asum => z_vect_asum
end type psb_z_vect_type
@ -903,6 +905,23 @@ contains
end function z_vect_nrm2
function z_vect_nrm2_weight(n,x,w) result(res)
implicit none
class(psb_z_vect_type), intent(inout) :: x
class(psb_z_vect_type), intent(inout) :: w
integer(psb_ipk_), intent(in) :: n
real(psb_dpk_) :: res
integer(psb_ipk_) :: info
if (allocated(x%v).and.allocated(w%v)) then
call w%v%mlt(x%v,info)
res = w%v%nrm2(n)
else
res = dzero
end if
end function z_vect_nrm2_weight
function z_vect_amax(n,x) result(res)
implicit none
class(psb_z_vect_type), intent(inout) :: x

@ -370,6 +370,113 @@ function psb_cnrm2_vect(x, desc_a, info,global) result(res)
return
end function psb_cnrm2_vect
! Function: psb_cnrm2_weight_vect
! Computes the weighted norm2 of a distributed vector,
!
! norm2 := sqrt ( (w.*X)**C * (w.*X))
!
! Arguments:
! x - type(psb_c_vect_type) The input vector containing the entries of X.
! w - type(psb_c_vect_type) The input vector containing the entries of W.
! desc_a - type(psb_desc_type). The communication descriptor.
! info - integer. Return code
! global - logical(optional) Whether to perform the global reduction, default: .true.
!
function psb_cnrm2_weight_vect(x,w, desc_a, info,global) result(res)
use psb_desc_mod
use psb_check_mod
use psb_error_mod
use psb_penv_mod
use psb_c_vect_mod
implicit none
real(psb_spk_) :: res
type(psb_c_vect_type), intent (inout) :: x
type(psb_c_vect_type), intent (inout) :: w
type(psb_desc_type), intent(in) :: desc_a
integer(psb_ipk_), intent(out) :: info
logical, intent(in), optional :: global
! locals
integer(psb_ipk_) :: ictxt, np, me,&
& err_act, iix, jjx, ndim, i, id, idx, ndm, ldx
integer(psb_lpk_) :: ix, jx, iy, ijy, m
logical :: global_
real(psb_spk_) :: snrm2, dd
character(len=20) :: name, ch_err
name='psb_cnrm2v_weight'
if (psb_errstatus_fatal()) then
info = psb_err_internal_error_ ; goto 9999
end if
info=psb_success_
call psb_erractionsave(err_act)
ictxt=desc_a%get_context()
call psb_info(ictxt, me, np)
if (np == -1) then
info=psb_err_context_error_
call psb_errpush(info,name)
goto 9999
endif
if (.not.allocated(x%v)) then
info = psb_err_invalid_vect_state_
call psb_errpush(info,name)
goto 9999
endif
if (present(global)) then
global_ = global
else
global_ = .true.
end if
ix = 1
jx = 1
m = desc_a%get_global_rows()
ldx = x%get_nrows()
call psb_chkvect(m,lone,ldx,ix,jx,desc_a,info,iix,jjx)
if(info /= psb_success_) then
info=psb_err_from_subroutine_
ch_err='psb_chkvect'
call psb_errpush(info,name,a_err=ch_err)
end if
if (iix /= 1) then
info=psb_err_ix_n1_iy_n1_unsupported_
call psb_errpush(info,name)
goto 9999
end if
if (desc_a%get_local_rows() > 0) then
ndim = desc_a%get_local_rows()
res = x%nrm2(ndim,w)
! adjust because overlapped elements are computed more than once
if (size(desc_a%ovrlap_elem,1)>0) then
if (x%is_dev()) call x%sync()
do i=1,size(desc_a%ovrlap_elem,1)
idx = desc_a%ovrlap_elem(i,1)
ndm = desc_a%ovrlap_elem(i,2)
dd = dble(ndm-1)/dble(ndm)
res = res - sqrt(cone - dd*(abs(x%v%v(idx))/res)**2)
end do
end if
else
res = szero
end if
if (global_) call psb_nrm2(ictxt,res)
call psb_erractionrestore(err_act)
return
9999 call psb_error_handler(ictxt,err_act)
return
end function psb_cnrm2_weight_vect
!!$
!!$ Parallel Sparse BLAS version 3.5

@ -370,6 +370,113 @@ function psb_dnrm2_vect(x, desc_a, info,global) result(res)
return
end function psb_dnrm2_vect
! Function: psb_dnrm2_weight_vect
! Computes the weighted norm2 of a distributed vector,
!
! norm2 := sqrt ( (w.*X)**C * (w.*X))
!
! Arguments:
! x - type(psb_d_vect_type) The input vector containing the entries of X.
! w - type(psb_d_vect_type) The input vector containing the entries of W.
! desc_a - type(psb_desc_type). The communication descriptor.
! info - integer. Return code
! global - logical(optional) Whether to perform the global reduction, default: .true.
!
function psb_dnrm2_weight_vect(x,w, desc_a, info,global) result(res)
use psb_desc_mod
use psb_check_mod
use psb_error_mod
use psb_penv_mod
use psb_d_vect_mod
implicit none
real(psb_dpk_) :: res
type(psb_d_vect_type), intent (inout) :: x
type(psb_d_vect_type), intent (inout) :: w
type(psb_desc_type), intent(in) :: desc_a
integer(psb_ipk_), intent(out) :: info
logical, intent(in), optional :: global
! locals
integer(psb_ipk_) :: ictxt, np, me,&
& err_act, iix, jjx, ndim, i, id, idx, ndm, ldx
integer(psb_lpk_) :: ix, jx, iy, ijy, m
logical :: global_
real(psb_dpk_) :: snrm2, dd
character(len=20) :: name, ch_err
name='psb_dnrm2v_weight'
if (psb_errstatus_fatal()) then
info = psb_err_internal_error_ ; goto 9999
end if
info=psb_success_
call psb_erractionsave(err_act)
ictxt=desc_a%get_context()
call psb_info(ictxt, me, np)
if (np == -1) then
info=psb_err_context_error_
call psb_errpush(info,name)
goto 9999
endif
if (.not.allocated(x%v)) then
info = psb_err_invalid_vect_state_
call psb_errpush(info,name)
goto 9999
endif
if (present(global)) then
global_ = global
else
global_ = .true.
end if
ix = 1
jx = 1
m = desc_a%get_global_rows()
ldx = x%get_nrows()
call psb_chkvect(m,lone,ldx,ix,jx,desc_a,info,iix,jjx)
if(info /= psb_success_) then
info=psb_err_from_subroutine_
ch_err='psb_chkvect'
call psb_errpush(info,name,a_err=ch_err)
end if
if (iix /= 1) then
info=psb_err_ix_n1_iy_n1_unsupported_
call psb_errpush(info,name)
goto 9999
end if
if (desc_a%get_local_rows() > 0) then
ndim = desc_a%get_local_rows()
res = x%nrm2(ndim,w)
! adjust because overlapped elements are computed more than once
if (size(desc_a%ovrlap_elem,1)>0) then
if (x%is_dev()) call x%sync()
do i=1,size(desc_a%ovrlap_elem,1)
idx = desc_a%ovrlap_elem(i,1)
ndm = desc_a%ovrlap_elem(i,2)
dd = dble(ndm-1)/dble(ndm)
res = res - sqrt(done - dd*(abs(x%v%v(idx))/res)**2)
end do
end if
else
res = dzero
end if
if (global_) call psb_nrm2(ictxt,res)
call psb_erractionrestore(err_act)
return
9999 call psb_error_handler(ictxt,err_act)
return
end function psb_dnrm2_weight_vect
!!$
!!$ Parallel Sparse BLAS version 3.5

@ -370,6 +370,113 @@ function psb_snrm2_vect(x, desc_a, info,global) result(res)
return
end function psb_snrm2_vect
! Function: psb_snrm2_weight_vect
! Computes the weighted norm2 of a distributed vector,
!
! norm2 := sqrt ( (w.*X)**C * (w.*X))
!
! Arguments:
! x - type(psb_s_vect_type) The input vector containing the entries of X.
! w - type(psb_s_vect_type) The input vector containing the entries of W.
! desc_a - type(psb_desc_type). The communication descriptor.
! info - integer. Return code
! global - logical(optional) Whether to perform the global reduction, default: .true.
!
function psb_snrm2_weight_vect(x,w, desc_a, info,global) result(res)
use psb_desc_mod
use psb_check_mod
use psb_error_mod
use psb_penv_mod
use psb_s_vect_mod
implicit none
real(psb_spk_) :: res
type(psb_s_vect_type), intent (inout) :: x
type(psb_s_vect_type), intent (inout) :: w
type(psb_desc_type), intent(in) :: desc_a
integer(psb_ipk_), intent(out) :: info
logical, intent(in), optional :: global
! locals
integer(psb_ipk_) :: ictxt, np, me,&
& err_act, iix, jjx, ndim, i, id, idx, ndm, ldx
integer(psb_lpk_) :: ix, jx, iy, ijy, m
logical :: global_
real(psb_spk_) :: snrm2, dd
character(len=20) :: name, ch_err
name='psb_snrm2v_weight'
if (psb_errstatus_fatal()) then
info = psb_err_internal_error_ ; goto 9999
end if
info=psb_success_
call psb_erractionsave(err_act)
ictxt=desc_a%get_context()
call psb_info(ictxt, me, np)
if (np == -1) then
info=psb_err_context_error_
call psb_errpush(info,name)
goto 9999
endif
if (.not.allocated(x%v)) then
info = psb_err_invalid_vect_state_
call psb_errpush(info,name)
goto 9999
endif
if (present(global)) then
global_ = global
else
global_ = .true.
end if
ix = 1
jx = 1
m = desc_a%get_global_rows()
ldx = x%get_nrows()
call psb_chkvect(m,lone,ldx,ix,jx,desc_a,info,iix,jjx)
if(info /= psb_success_) then
info=psb_err_from_subroutine_
ch_err='psb_chkvect'
call psb_errpush(info,name,a_err=ch_err)
end if
if (iix /= 1) then
info=psb_err_ix_n1_iy_n1_unsupported_
call psb_errpush(info,name)
goto 9999
end if
if (desc_a%get_local_rows() > 0) then
ndim = desc_a%get_local_rows()
res = x%nrm2(ndim,w)
! adjust because overlapped elements are computed more than once
if (size(desc_a%ovrlap_elem,1)>0) then
if (x%is_dev()) call x%sync()
do i=1,size(desc_a%ovrlap_elem,1)
idx = desc_a%ovrlap_elem(i,1)
ndm = desc_a%ovrlap_elem(i,2)
dd = dble(ndm-1)/dble(ndm)
res = res - sqrt(sone - dd*(abs(x%v%v(idx))/res)**2)
end do
end if
else
res = szero
end if
if (global_) call psb_nrm2(ictxt,res)
call psb_erractionrestore(err_act)
return
9999 call psb_error_handler(ictxt,err_act)
return
end function psb_snrm2_weight_vect
!!$
!!$ Parallel Sparse BLAS version 3.5

@ -370,6 +370,113 @@ function psb_znrm2_vect(x, desc_a, info,global) result(res)
return
end function psb_znrm2_vect
! Function: psb_znrm2_weight_vect
! Computes the weighted norm2 of a distributed vector,
!
! norm2 := sqrt ( (w.*X)**C * (w.*X))
!
! Arguments:
! x - type(psb_z_vect_type) The input vector containing the entries of X.
! w - type(psb_z_vect_type) The input vector containing the entries of W.
! desc_a - type(psb_desc_type). The communication descriptor.
! info - integer. Return code
! global - logical(optional) Whether to perform the global reduction, default: .true.
!
function psb_znrm2_weight_vect(x,w, desc_a, info,global) result(res)
use psb_desc_mod
use psb_check_mod
use psb_error_mod
use psb_penv_mod
use psb_z_vect_mod
implicit none
real(psb_dpk_) :: res
type(psb_z_vect_type), intent (inout) :: x
type(psb_z_vect_type), intent (inout) :: w
type(psb_desc_type), intent(in) :: desc_a
integer(psb_ipk_), intent(out) :: info
logical, intent(in), optional :: global
! locals
integer(psb_ipk_) :: ictxt, np, me,&
& err_act, iix, jjx, ndim, i, id, idx, ndm, ldx
integer(psb_lpk_) :: ix, jx, iy, ijy, m
logical :: global_
real(psb_dpk_) :: snrm2, dd
character(len=20) :: name, ch_err
name='psb_znrm2v_weight'
if (psb_errstatus_fatal()) then
info = psb_err_internal_error_ ; goto 9999
end if
info=psb_success_
call psb_erractionsave(err_act)
ictxt=desc_a%get_context()
call psb_info(ictxt, me, np)
if (np == -1) then
info=psb_err_context_error_
call psb_errpush(info,name)
goto 9999
endif
if (.not.allocated(x%v)) then
info = psb_err_invalid_vect_state_
call psb_errpush(info,name)
goto 9999
endif
if (present(global)) then
global_ = global
else
global_ = .true.
end if
ix = 1
jx = 1
m = desc_a%get_global_rows()
ldx = x%get_nrows()
call psb_chkvect(m,lone,ldx,ix,jx,desc_a,info,iix,jjx)
if(info /= psb_success_) then
info=psb_err_from_subroutine_
ch_err='psb_chkvect'
call psb_errpush(info,name,a_err=ch_err)
end if
if (iix /= 1) then
info=psb_err_ix_n1_iy_n1_unsupported_
call psb_errpush(info,name)
goto 9999
end if
if (desc_a%get_local_rows() > 0) then
ndim = desc_a%get_local_rows()
res = x%nrm2(ndim,w)
! adjust because overlapped elements are computed more than once
if (size(desc_a%ovrlap_elem,1)>0) then
if (x%is_dev()) call x%sync()
do i=1,size(desc_a%ovrlap_elem,1)
idx = desc_a%ovrlap_elem(i,1)
ndm = desc_a%ovrlap_elem(i,2)
dd = dble(ndm-1)/dble(ndm)
res = res - sqrt(zone - dd*(abs(x%v%v(idx))/res)**2)
end do
end if
else
res = dzero
end if
if (global_) call psb_nrm2(ictxt,res)
call psb_erractionrestore(err_act)
return
9999 call psb_error_handler(ictxt,err_act)
return
end function psb_znrm2_weight_vect
!!$
!!$ Parallel Sparse BLAS version 3.5

@ -53,7 +53,7 @@ program vecoperation
integer(psb_ipk_) :: nr, nlr, info, i, ii, ib=1
integer(psb_lpk_) :: nt
integer(psb_lpk_), allocatable :: myidx(:)
real(psb_dpk_) :: zt(1), dotresult, norm2, norm1, norminf
real(psb_dpk_) :: zt(1), dotresult, norm2, norm1, norminf, norm2w
character(len=20) :: name,ch_err,readinput
real(psb_dpk_), allocatable :: vx(:), vy(:), vz(:)
real(psb_dpk_) :: c
@ -160,10 +160,6 @@ program vecoperation
!
dotresult = psb_gedot(x,y,desc_a,info) ! Dot-product
if (iam == psb_root_) write(psb_out_unit,'("Dot product result : ",es12.5)')dotresult
norm1 = psb_norm1(x,desc_a,info)
norm2 = psb_norm2(x,desc_a,info)
norminf = psb_normi(x,desc_a,info)
if (iam == psb_root_) write(psb_out_unit,'("\|x\|_inf : ",es12.5," \|x\|_1 :",es12.5," \|x\|_2",es12.5)')norminf,norm1,norm2
call psb_geaxpby(1.0_psb_dpk_, x, 1.0_psb_dpk_, y, desc_a, info) ! \alpha x + \beta y
if (iam == psb_root_) then
@ -218,7 +214,7 @@ program vecoperation
write(psb_out_unit,'("|x| = ",es12.1)')vz(:)
end if
c = 1.0/2.0;
c = 1.0/5.0;
call psb_gecmp(x,c,z,desc_a,info);
if (iam == psb_root_) then
@ -229,6 +225,16 @@ program vecoperation
write(psb_out_unit,'("z = ",es12.1)')vz(:)
end if
write(psb_out_unit,'("Computation of vector norms:")')
norm1 = psb_norm1(x,desc_a,info)
norm2 = psb_norm2(x,desc_a,info)
norm2w = psb_norm2(x,absz,desc_a,info)
norminf = psb_normi(x,desc_a,info)
write(psb_out_unit,'("\|x\|_inf : ",es12.5," \|x\|_1 :",es12.5)')norminf,norm1
write(psb_out_unit,'(" \|x\|_2 : ",es12.5," \|x\|_2,w : ",es12.5)')norm2,norm2w
!
! cleanup storage and exit
!

Loading…
Cancel
Save