psblas3-matasb:

base/modules/psb_base_mat_mod.f90
 base/modules/psb_d_base_mat_mod.f90
 base/modules/psb_d_base_vect_mod.f90
 base/modules/psb_d_mat_mod.f90
 base/serial/impl/psb_d_mat_impl.F90
 base/tools/psb_dins.f90
 base/tools/psb_dspasb.f90

Mods to get imaging application up & running
psblas-3.3.1-1
Salvatore Filippone 10 years ago
parent 047d928ed5
commit cf679ebc1a

@ -786,7 +786,7 @@ contains
!
subroutine psb_base_mat_sync(a)
implicit none
class(psb_base_sparse_mat), intent(inout) :: a
class(psb_base_sparse_mat), target, intent(in) :: a
end subroutine psb_base_mat_sync

@ -117,7 +117,7 @@ module psb_d_base_mat_mod
private :: d_base_mat_sync, d_base_mat_is_host, d_base_mat_is_dev, &
& d_base_mat_is_sync, d_base_mat_set_host, d_base_mat_set_dev,&
& d_base_mat_set_sync, d_base_mat_asb
& d_base_mat_set_sync
!> \namespace psb_base_mod \class psb_d_coo_sparse_mat
!! \extends psb_d_base_mat_mod::psb_d_base_sparse_mat

@ -349,6 +349,7 @@ contains
! !$ goto 9999
end select
end if
call x%set_host()
if (info /= 0) then
call psb_errpush(info,'base_vect_ins')
return
@ -370,8 +371,9 @@ contains
info = 0
if (psb_errstatus_fatal()) return
call irl%sync()
call val%sync()
if (irl%is_dev()) call irl%sync()
if (val%is_dev()) call val%sync()
if (x%is_dev()) call x%sync()
call x%ins(n,irl%v,val%v,dupl,info)
if (info /= 0) then
@ -425,7 +427,7 @@ contains
& call psb_realloc(n,x%v,info)
if (info /= 0) &
& call psb_errpush(psb_err_alloc_dealloc_,'vect_asb')
call x%sync()
end subroutine d_base_asb
@ -609,7 +611,7 @@ contains
real(psb_dpk_), allocatable :: res(:)
integer(psb_ipk_) :: info
if (.not.allocated(x%v)) return
!if (.not.allocated(x%v)) return
call x%sync()
allocate(res(x%get_nrows()),stat=info)
if (info /= 0) then

@ -511,9 +511,10 @@ module psb_d_mat_mod
end interface
interface
subroutine psb_d_asb(a)
import :: psb_ipk_, psb_dspmat_type
subroutine psb_d_asb(a,mold)
import :: psb_ipk_, psb_dspmat_type, psb_d_base_sparse_mat
class(psb_dspmat_type), intent(inout) :: a
class(psb_d_base_sparse_mat), optional, intent(in) :: mold
end subroutine psb_d_asb
end interface

@ -1917,12 +1917,14 @@ subroutine psb_d_transc_2mat(a,b)
end subroutine psb_d_transc_2mat
subroutine psb_d_asb(a)
subroutine psb_d_asb(a,mold)
use psb_d_mat_mod, psb_protect_name => psb_d_asb
use psb_error_mod
implicit none
class(psb_dspmat_type), intent(inout) :: a
class(psb_d_base_sparse_mat), optional, intent(in) :: mold
class(psb_d_base_sparse_mat), allocatable :: tmp
integer(psb_ipk_) :: err_act, info
character(len=20) :: name='reinit'
@ -1934,6 +1936,15 @@ subroutine psb_d_asb(a)
endif
call a%a%asb()
if (present(mold)) then
if (.not.same_type_as(a%a,mold)) then
allocate(tmp,mold=mold)
call tmp%mv_from_fmt(a%a,info)
call a%a%free()
call move_alloc(tmp,a%a)
end if
end if
call psb_erractionrestore(err_act)
return

@ -336,7 +336,7 @@ subroutine psb_dins_vect_v(m, irw, val, x, desc_a, info, dupl,local)
if (psb_errstatus_fatal()) return
info=psb_success_
call psb_erractionsave(err_act)
name = 'psb_dinsvi'
name = 'psb_dinsvi_vect_v'
if (.not.desc_a%is_ok()) then
info = psb_err_invalid_cd_state_
@ -399,12 +399,12 @@ subroutine psb_dins_vect_v(m, irw, val, x, desc_a, info, dupl,local)
lval = val%get_vect()
call desc_a%indxmap%g2lip(irl(1:m),info,owned=.true.)
call x%ins(m,irl,lval,dupl_,info)
end if
if (info /= 0) then
call psb_errpush(info,name)
goto 9999
end if
deallocate(irl)
call psb_erractionrestore(err_act)
return

@ -116,7 +116,7 @@ subroutine psb_dspasb(a,desc_a, info, afmt, upd, dupl, mold)
if (a%is_bld()) then
call a%cscnv(info,type=afmt,dupl=dupl, mold=mold)
else if (a%is_upd()) then
call a%asb()
call a%asb(mold=mold)
else
info = psb_err_invalid_mat_state_
call psb_errpush(info,name)

Loading…
Cancel
Save