psblas3-matasb:

base/modules/psb_base_mat_mod.f90
 base/modules/psb_d_base_vect_mod.f90
 base/modules/psb_d_mat_mod.f90
 base/modules/psb_i_base_vect_mod.f90

Redefine get_vect for multivector to use n_rows/n_cols.
psblas-3.3.1-1
Salvatore Filippone 11 years ago
parent cf679ebc1a
commit 85774c22ec

@ -124,6 +124,8 @@ module psb_base_mat_mod
logical, private :: unitd
!> Are the coefficients sorted ?
logical, private :: sorted
logical, private :: repeatable_updates=.false.
contains
! == = =================================
@ -152,6 +154,7 @@ module psb_base_mat_mod
procedure, pass(a) :: is_unit => psb_base_is_unit
procedure, pass(a) :: is_by_rows => psb_base_is_by_rows
procedure, pass(a) :: is_by_cols => psb_base_is_by_cols
procedure, pass(a) :: is_repeatable_updates => psb_base_is_repeatable_updates
! == = =================================
!
@ -172,6 +175,8 @@ module psb_base_mat_mod
procedure, pass(a) :: set_triangle => psb_base_set_triangle
procedure, pass(a) :: set_unit => psb_base_set_unit
procedure, pass(a) :: set_repeatable_updates => psb_base_set_repeatable_updates
! == = =================================
!
@ -614,6 +619,18 @@ contains
end if
end subroutine psb_base_set_upper
subroutine psb_base_set_repeatable_updates(a,val)
implicit none
class(psb_base_sparse_mat), intent(inout) :: a
logical, intent(in), optional :: val
if (present(val)) then
a%repeatable_updates = val
else
a%repeatable_updates = .true.
end if
end subroutine psb_base_set_repeatable_updates
function psb_base_is_triangle(a) result(res)
implicit none
class(psb_base_sparse_mat), intent(in) :: a
@ -692,6 +709,13 @@ contains
res = .false.
end function psb_base_is_by_cols
function psb_base_is_repeatable_updates(a) result(res)
implicit none
class(psb_base_sparse_mat), intent(in) :: a
logical :: res
res = a%repeatable_updates
end function psb_base_is_repeatable_updates
!
! TRANSP: note sorted=.false.
@ -712,6 +736,7 @@ contains
b%unitd = a%unitd
b%upper = .not.a%upper
b%sorted = .false.
b%repeatable_updates = .false.
end subroutine psb_base_transp_2mat
@ -730,6 +755,7 @@ contains
b%unitd = a%unitd
b%upper = .not.a%upper
b%sorted = .false.
b%repeatable_updates = .false.
end subroutine psb_base_transc_2mat
@ -748,6 +774,7 @@ contains
a%unitd = a%unitd
a%upper = .not.a%upper
a%sorted = .false.
a%repeatable_updates = .false.
end subroutine psb_base_transp_1mat

@ -1704,16 +1704,17 @@ contains
function d_base_mv_get_vect(x) result(res)
class(psb_d_base_multivect_type), intent(inout) :: x
real(psb_dpk_), allocatable :: res(:,:)
integer(psb_ipk_) :: info
integer(psb_ipk_) :: info,m,n
m = x%get_nrows()
n = x%get_ncols()
if (.not.allocated(x%v)) return
call x%sync()
allocate(res(x%get_nrows(),x%get_ncols()),stat=info)
allocate(res(m,n),stat=info)
if (info /= 0) then
call psb_errpush(psb_err_alloc_dealloc_,'base_mv_get_vect')
return
end if
res(:,:) = x%v(:,:)
res(1:m,1:n) = x%v(1:m,1:n)
end function d_base_mv_get_vect
!

@ -98,6 +98,7 @@ module psb_d_mat_mod
procedure, pass(a) :: is_lower => psb_d_is_lower
procedure, pass(a) :: is_triangle => psb_d_is_triangle
procedure, pass(a) :: is_unit => psb_d_is_unit
procedure, pass(a) :: is_repeatable_updates => psb_d_is_repeatable_updates
procedure, pass(a) :: get_fmt => psb_d_get_fmt
procedure, pass(a) :: sizeof => psb_d_sizeof
@ -115,6 +116,9 @@ module psb_d_mat_mod
procedure, pass(a) :: set_triangle => psb_d_set_triangle
procedure, pass(a) :: set_unit => psb_d_set_unit
procedure, pass(a) :: set_repeatable_updates => psb_d_set_repeatable_updates
! Memory/data management
procedure, pass(a) :: csall => psb_d_csall
procedure, pass(a) :: free => psb_d_free
@ -1109,6 +1113,31 @@ contains
end function psb_d_is_by_cols
function psb_d_is_repeatable_updates(a) result(res)
implicit none
class(psb_dspmat_type), intent(in) :: a
logical :: res
if (allocated(a%a)) then
res = a%a%is_repeatable_updates()
else
res = .false.
end if
end function psb_d_is_repeatable_updates
subroutine psb_d_set_repeatable_updates(a,val)
implicit none
class(psb_dspmat_type), intent(inout) :: a
logical, intent(in), optional :: val
if (allocated(a%a)) then
call a%a%set_repeatable_updates(val)
end if
end subroutine psb_d_set_repeatable_updates
function psb_d_get_nzeros(a) result(res)
implicit none

@ -1670,16 +1670,19 @@ contains
function i_base_mv_get_vect(x) result(res)
class(psb_i_base_multivect_type), intent(inout) :: x
integer(psb_ipk_), allocatable :: res(:,:)
integer(psb_ipk_) :: info
integer(psb_ipk_) :: info, m, n
m = x%get_nrows()
n = x%get_ncols()
if (.not.allocated(x%v)) return
call x%sync()
allocate(res(x%get_nrows(),x%get_ncols()),stat=info)
allocate(res(m,n),stat=info)
if (info /= 0) then
call psb_errpush(psb_err_alloc_dealloc_,'base_mv_get_vect')
return
end if
res(:,:) = x%v(:,:)
res(1:m,1:n) = x%v(1:m,1:n)
end function i_base_mv_get_vect
!

Loading…
Cancel
Save