Infrastructure for non-local matrix build

remotebuild
Salvatore Filippone 3 years ago
parent 7c3852109f
commit 1337009f91

@ -746,7 +746,8 @@ contains
if (info >=0) then
if (nxt == lip) then
ncol = max(nxt,ncol)
call psb_ensure_size(ncol,idxmap%loc_to_glob,info,pad=-1_psb_lpk_,addsz=laddsz)
call psb_ensure_size(ncol,idxmap%loc_to_glob,info,&
& pad=-1_psb_lpk_,addsz=laddsz)
if (info /= psb_success_) then
info=1
call psb_errpush(psb_err_from_subroutine_ai_,name,&

@ -204,6 +204,9 @@ module psb_const_mod
integer(psb_ipk_), parameter :: psb_spmat_null_=0, psb_spmat_bld_=1
integer(psb_ipk_), parameter :: psb_spmat_asb_=2, psb_spmat_upd_=4
integer(psb_ipk_), parameter :: psb_matbld_noremote_=0, psb_matbld_remote_=1
integer(psb_ipk_), parameter :: psb_ireg_flgs_=10, psb_ip2_=0
integer(psb_ipk_), parameter :: psb_iflag_=2, psb_ichk_=3
integer(psb_ipk_), parameter :: psb_nnzt_=4, psb_zero_=5,psb_ipc_=6

@ -85,6 +85,8 @@ module psb_c_mat_mod
type :: psb_cspmat_type
class(psb_c_base_sparse_mat), allocatable :: a
integer(psb_ipk_) :: remote_build=psb_matbld_noremote_
class(psb_lc_base_sparse_mat), allocatable :: rmta
contains
! Getters

@ -85,6 +85,8 @@ module psb_d_mat_mod
type :: psb_dspmat_type
class(psb_d_base_sparse_mat), allocatable :: a
integer(psb_ipk_) :: remote_build=psb_matbld_noremote_
class(psb_ld_base_sparse_mat), allocatable :: rmta
contains
! Getters

@ -85,6 +85,8 @@ module psb_s_mat_mod
type :: psb_sspmat_type
class(psb_s_base_sparse_mat), allocatable :: a
integer(psb_ipk_) :: remote_build=psb_matbld_noremote_
class(psb_ls_base_sparse_mat), allocatable :: rmta
contains
! Getters

@ -85,6 +85,8 @@ module psb_z_mat_mod
type :: psb_zspmat_type
class(psb_z_base_sparse_mat), allocatable :: a
integer(psb_ipk_) :: remote_build=psb_matbld_noremote_
class(psb_lz_base_sparse_mat), allocatable :: rmta
contains
! Getters

@ -239,13 +239,13 @@ Module psb_c_tools_mod
interface psb_spall
subroutine psb_cspalloc(a, desc_a, info, nnz)
subroutine psb_cspalloc(a, desc_a, info, nnz, bldmode)
import
implicit none
type(psb_desc_type), intent(in) :: desc_a
type(psb_cspmat_type), intent(inout) :: a
integer(psb_ipk_), intent(out) :: info
integer(psb_ipk_), optional, intent(in) :: nnz
integer(psb_ipk_), optional, intent(in) :: nnz, bldmode
end subroutine psb_cspalloc
end interface

@ -239,13 +239,13 @@ Module psb_d_tools_mod
interface psb_spall
subroutine psb_dspalloc(a, desc_a, info, nnz)
subroutine psb_dspalloc(a, desc_a, info, nnz, bldmode)
import
implicit none
type(psb_desc_type), intent(in) :: desc_a
type(psb_dspmat_type), intent(inout) :: a
integer(psb_ipk_), intent(out) :: info
integer(psb_ipk_), optional, intent(in) :: nnz
integer(psb_ipk_), optional, intent(in) :: nnz, bldmode
end subroutine psb_dspalloc
end interface

@ -239,13 +239,13 @@ Module psb_s_tools_mod
interface psb_spall
subroutine psb_sspalloc(a, desc_a, info, nnz)
subroutine psb_sspalloc(a, desc_a, info, nnz, bldmode)
import
implicit none
type(psb_desc_type), intent(in) :: desc_a
type(psb_sspmat_type), intent(inout) :: a
integer(psb_ipk_), intent(out) :: info
integer(psb_ipk_), optional, intent(in) :: nnz
integer(psb_ipk_), optional, intent(in) :: nnz, bldmode
end subroutine psb_sspalloc
end interface

@ -239,13 +239,13 @@ Module psb_z_tools_mod
interface psb_spall
subroutine psb_zspalloc(a, desc_a, info, nnz)
subroutine psb_zspalloc(a, desc_a, info, nnz, bldmode)
import
implicit none
type(psb_desc_type), intent(in) :: desc_a
type(psb_zspmat_type), intent(inout) :: a
integer(psb_ipk_), intent(out) :: info
integer(psb_ipk_), optional, intent(in) :: nnz
integer(psb_ipk_), optional, intent(in) :: nnz, bldmode
end subroutine psb_zspalloc
end interface

@ -6121,7 +6121,6 @@ subroutine psb_lc_coo_csput_a(nz,ia,ja,val,a,imin,imax,jmin,jmax,info)
if (nz == 0) return
nza = a%get_nzeros()
isza = a%get_size()
if (a%is_bld()) then

@ -675,7 +675,12 @@ subroutine psb_c_free(a)
call a%a%free()
deallocate(a%a)
endif
if (allocated(a%rmta)) then
call a%rmta%free()
deallocate(a%rmta)
end if
a%remote_build = psb_matbld_noremote_
end subroutine psb_c_free

@ -6121,7 +6121,6 @@ subroutine psb_ld_coo_csput_a(nz,ia,ja,val,a,imin,imax,jmin,jmax,info)
if (nz == 0) return
nza = a%get_nzeros()
isza = a%get_size()
if (a%is_bld()) then

@ -675,7 +675,12 @@ subroutine psb_d_free(a)
call a%a%free()
deallocate(a%a)
endif
if (allocated(a%rmta)) then
call a%rmta%free()
deallocate(a%rmta)
end if
a%remote_build = psb_matbld_noremote_
end subroutine psb_d_free

@ -6121,7 +6121,6 @@ subroutine psb_ls_coo_csput_a(nz,ia,ja,val,a,imin,imax,jmin,jmax,info)
if (nz == 0) return
nza = a%get_nzeros()
isza = a%get_size()
if (a%is_bld()) then

@ -675,7 +675,12 @@ subroutine psb_s_free(a)
call a%a%free()
deallocate(a%a)
endif
if (allocated(a%rmta)) then
call a%rmta%free()
deallocate(a%rmta)
end if
a%remote_build = psb_matbld_noremote_
end subroutine psb_s_free

@ -6121,7 +6121,6 @@ subroutine psb_lz_coo_csput_a(nz,ia,ja,val,a,imin,imax,jmin,jmax,info)
if (nz == 0) return
nza = a%get_nzeros()
isza = a%get_size()
if (a%is_bld()) then

@ -675,7 +675,12 @@ subroutine psb_z_free(a)
call a%a%free()
deallocate(a%a)
endif
if (allocated(a%rmta)) then
call a%rmta%free()
deallocate(a%rmta)
end if
a%remote_build = psb_matbld_noremote_
end subroutine psb_z_free

@ -41,7 +41,7 @@
! nnz - integer(optional). The number of nonzeroes in the matrix.
! (local, user estimate)
!
subroutine psb_cspalloc(a, desc_a, info, nnz)
subroutine psb_cspalloc(a, desc_a, info, nnz, bldmode)
use psb_base_mod, psb_protect_name => psb_cspalloc
implicit none
@ -49,13 +49,13 @@ subroutine psb_cspalloc(a, desc_a, info, nnz)
type(psb_desc_type), intent(in) :: desc_a
type(psb_cspmat_type), intent(inout) :: a
integer(psb_ipk_), intent(out) :: info
integer(psb_ipk_), optional, intent(in) :: nnz
integer(psb_ipk_), optional, intent(in) :: nnz, bldmode
!locals
type(psb_ctxt_type) :: ctxt
integer(psb_ipk_) :: np, me, err_act
integer(psb_ipk_) :: loc_row,loc_col, nnz_, dectype
integer(psb_lpk_) :: m, n
integer(psb_ipk_) :: loc_row,loc_col, nnz_, dectype, bldmode_
integer(psb_lpk_) :: m, n, nnzrmt_
integer(psb_ipk_) :: debug_level, debug_unit
character(len=20) :: name
@ -96,7 +96,12 @@ subroutine psb_cspalloc(a, desc_a, info, nnz)
else
nnz_ = max(1,5*loc_row)
endif
if (present(bldmode)) then
bldmode_ = bldmode
else
bldmode_ = psb_matbld_noremote_
end if
if (debug_level >= psb_debug_ext_) &
& write(debug_unit,*) me,' ',trim(name), &
& ':allocating size:',loc_row,loc_col,nnz_
@ -109,6 +114,25 @@ subroutine psb_cspalloc(a, desc_a, info, nnz)
goto 9999
end if
write(0,*) name,'Setting a%remote_build ',&
& bldmode_,psb_matbld_noremote_,psb_matbld_remote_
a%remote_build = bldmode_
select case(a%remote_build)
case (psb_matbld_noremote_)
! nothing needed
write(0,*) name,' matbld_noremote_ nothing needed'
case (psb_matbld_remote_)
write(0,*) name,' matbld_remote_ start '
allocate(psb_lc_coo_sparse_mat :: a%rmta)
nnzrmt_ = max(100,(nnz_/100))
call a%rmta%allocate(m,n,nnzrmt_)
case default
write(0,*) name,'Invalid value for remote_build '
a%remote_build = psb_matbld_noremote_
end select
if (debug_level >= psb_debug_ext_) &
& write(debug_unit,*) me,' ',trim(name),': ', &
& desc_a%get_dectype(),psb_desc_bld_

@ -111,7 +111,15 @@ subroutine psb_cspasb(a,desc_a, info, afmt, upd, dupl, mold)
call a%set_ncols(n_col)
end if
if (a%is_bld()) then
if (a%is_bld()) then
select case(a%remote_build)
case (psb_matbld_noremote_)
! nothing needed
case (psb_matbld_remote_)
write(0,*) 'Need to implement data movement! '
end select
call a%cscnv(info,type=afmt,dupl=dupl, mold=mold)
else if (a%is_upd()) then
call a%asb(mold=mold)

@ -70,6 +70,10 @@ subroutine psb_cspins(nz,ia,ja,val,a,desc_a,info,rebuild,local)
integer(psb_ipk_), parameter :: relocsz=200
logical :: rebuild_, local_
integer(psb_ipk_), allocatable :: ila(:),jla(:)
integer(psb_ipk_) :: i,k
integer(psb_lpk_) :: nnl
integer(psb_lpk_), allocatable :: lila(:),ljla(:)
complex(psb_spk_), allocatable :: lval(:)
character(len=20) :: name
info = psb_success_
@ -147,6 +151,30 @@ subroutine psb_cspins(nz,ia,ja,val,a,desc_a,info,rebuild,local)
call psb_errpush(info,name,a_err='a%csput')
goto 9999
end if
select case(a%remote_build)
case (psb_matbld_noremote_)
! Do nothing
case (psb_matbld_remote_)
nnl = count(ila(1:nz)<0)
allocate(lila(nnl),ljla(nnl),lval(nnl))
k = 0
do i=1,nz
if (ila(i)<0) then
k=k+1
lila(k) = ia(k)
ljla(k) = ja(k)
lval(k) = val(k)
end if
end do
if (k /= nnl) write(0,*) name,' Wrong conversion?',k,nnl
call a%rmta%csput(nnl,lila,ljla,lval,1_psb_lpk_,desc_a%get_global_rows(),&
& 1_psb_lpk_,desc_a%get_global_rows(),info)
case default
write(0,*) name,' Ignoring wrong value for %remote_build'
end select
else
info = psb_err_invalid_a_and_cd_state_
call psb_errpush(info,name)

@ -41,7 +41,7 @@
! nnz - integer(optional). The number of nonzeroes in the matrix.
! (local, user estimate)
!
subroutine psb_dspalloc(a, desc_a, info, nnz)
subroutine psb_dspalloc(a, desc_a, info, nnz, bldmode)
use psb_base_mod, psb_protect_name => psb_dspalloc
implicit none
@ -49,13 +49,13 @@ subroutine psb_dspalloc(a, desc_a, info, nnz)
type(psb_desc_type), intent(in) :: desc_a
type(psb_dspmat_type), intent(inout) :: a
integer(psb_ipk_), intent(out) :: info
integer(psb_ipk_), optional, intent(in) :: nnz
integer(psb_ipk_), optional, intent(in) :: nnz, bldmode
!locals
type(psb_ctxt_type) :: ctxt
integer(psb_ipk_) :: np, me, err_act
integer(psb_ipk_) :: loc_row,loc_col, nnz_, dectype
integer(psb_lpk_) :: m, n
integer(psb_ipk_) :: loc_row,loc_col, nnz_, dectype, bldmode_
integer(psb_lpk_) :: m, n, nnzrmt_
integer(psb_ipk_) :: debug_level, debug_unit
character(len=20) :: name
@ -96,7 +96,12 @@ subroutine psb_dspalloc(a, desc_a, info, nnz)
else
nnz_ = max(1,5*loc_row)
endif
if (present(bldmode)) then
bldmode_ = bldmode
else
bldmode_ = psb_matbld_noremote_
end if
if (debug_level >= psb_debug_ext_) &
& write(debug_unit,*) me,' ',trim(name), &
& ':allocating size:',loc_row,loc_col,nnz_
@ -109,6 +114,25 @@ subroutine psb_dspalloc(a, desc_a, info, nnz)
goto 9999
end if
write(0,*) name,'Setting a%remote_build ',&
& bldmode_,psb_matbld_noremote_,psb_matbld_remote_
a%remote_build = bldmode_
select case(a%remote_build)
case (psb_matbld_noremote_)
! nothing needed
write(0,*) name,' matbld_noremote_ nothing needed'
case (psb_matbld_remote_)
write(0,*) name,' matbld_remote_ start '
allocate(psb_ld_coo_sparse_mat :: a%rmta)
nnzrmt_ = max(100,(nnz_/100))
call a%rmta%allocate(m,n,nnzrmt_)
case default
write(0,*) name,'Invalid value for remote_build '
a%remote_build = psb_matbld_noremote_
end select
if (debug_level >= psb_debug_ext_) &
& write(debug_unit,*) me,' ',trim(name),': ', &
& desc_a%get_dectype(),psb_desc_bld_

@ -111,7 +111,15 @@ subroutine psb_dspasb(a,desc_a, info, afmt, upd, dupl, mold)
call a%set_ncols(n_col)
end if
if (a%is_bld()) then
if (a%is_bld()) then
select case(a%remote_build)
case (psb_matbld_noremote_)
! nothing needed
case (psb_matbld_remote_)
write(0,*) 'Need to implement data movement! '
end select
call a%cscnv(info,type=afmt,dupl=dupl, mold=mold)
else if (a%is_upd()) then
call a%asb(mold=mold)

@ -70,6 +70,10 @@ subroutine psb_dspins(nz,ia,ja,val,a,desc_a,info,rebuild,local)
integer(psb_ipk_), parameter :: relocsz=200
logical :: rebuild_, local_
integer(psb_ipk_), allocatable :: ila(:),jla(:)
integer(psb_ipk_) :: i,k
integer(psb_lpk_) :: nnl
integer(psb_lpk_), allocatable :: lila(:),ljla(:)
real(psb_dpk_), allocatable :: lval(:)
character(len=20) :: name
info = psb_success_
@ -147,6 +151,30 @@ subroutine psb_dspins(nz,ia,ja,val,a,desc_a,info,rebuild,local)
call psb_errpush(info,name,a_err='a%csput')
goto 9999
end if
select case(a%remote_build)
case (psb_matbld_noremote_)
! Do nothing
case (psb_matbld_remote_)
nnl = count(ila(1:nz)<0)
allocate(lila(nnl),ljla(nnl),lval(nnl))
k = 0
do i=1,nz
if (ila(i)<0) then
k=k+1
lila(k) = ia(k)
ljla(k) = ja(k)
lval(k) = val(k)
end if
end do
if (k /= nnl) write(0,*) name,' Wrong conversion?',k,nnl
call a%rmta%csput(nnl,lila,ljla,lval,1_psb_lpk_,desc_a%get_global_rows(),&
& 1_psb_lpk_,desc_a%get_global_rows(),info)
case default
write(0,*) name,' Ignoring wrong value for %remote_build'
end select
else
info = psb_err_invalid_a_and_cd_state_
call psb_errpush(info,name)

@ -41,7 +41,7 @@
! nnz - integer(optional). The number of nonzeroes in the matrix.
! (local, user estimate)
!
subroutine psb_sspalloc(a, desc_a, info, nnz)
subroutine psb_sspalloc(a, desc_a, info, nnz, bldmode)
use psb_base_mod, psb_protect_name => psb_sspalloc
implicit none
@ -49,13 +49,13 @@ subroutine psb_sspalloc(a, desc_a, info, nnz)
type(psb_desc_type), intent(in) :: desc_a
type(psb_sspmat_type), intent(inout) :: a
integer(psb_ipk_), intent(out) :: info
integer(psb_ipk_), optional, intent(in) :: nnz
integer(psb_ipk_), optional, intent(in) :: nnz, bldmode
!locals
type(psb_ctxt_type) :: ctxt
integer(psb_ipk_) :: np, me, err_act
integer(psb_ipk_) :: loc_row,loc_col, nnz_, dectype
integer(psb_lpk_) :: m, n
integer(psb_ipk_) :: loc_row,loc_col, nnz_, dectype, bldmode_
integer(psb_lpk_) :: m, n, nnzrmt_
integer(psb_ipk_) :: debug_level, debug_unit
character(len=20) :: name
@ -96,7 +96,12 @@ subroutine psb_sspalloc(a, desc_a, info, nnz)
else
nnz_ = max(1,5*loc_row)
endif
if (present(bldmode)) then
bldmode_ = bldmode
else
bldmode_ = psb_matbld_noremote_
end if
if (debug_level >= psb_debug_ext_) &
& write(debug_unit,*) me,' ',trim(name), &
& ':allocating size:',loc_row,loc_col,nnz_
@ -109,6 +114,25 @@ subroutine psb_sspalloc(a, desc_a, info, nnz)
goto 9999
end if
write(0,*) name,'Setting a%remote_build ',&
& bldmode_,psb_matbld_noremote_,psb_matbld_remote_
a%remote_build = bldmode_
select case(a%remote_build)
case (psb_matbld_noremote_)
! nothing needed
write(0,*) name,' matbld_noremote_ nothing needed'
case (psb_matbld_remote_)
write(0,*) name,' matbld_remote_ start '
allocate(psb_ls_coo_sparse_mat :: a%rmta)
nnzrmt_ = max(100,(nnz_/100))
call a%rmta%allocate(m,n,nnzrmt_)
case default
write(0,*) name,'Invalid value for remote_build '
a%remote_build = psb_matbld_noremote_
end select
if (debug_level >= psb_debug_ext_) &
& write(debug_unit,*) me,' ',trim(name),': ', &
& desc_a%get_dectype(),psb_desc_bld_

@ -111,7 +111,15 @@ subroutine psb_sspasb(a,desc_a, info, afmt, upd, dupl, mold)
call a%set_ncols(n_col)
end if
if (a%is_bld()) then
if (a%is_bld()) then
select case(a%remote_build)
case (psb_matbld_noremote_)
! nothing needed
case (psb_matbld_remote_)
write(0,*) 'Need to implement data movement! '
end select
call a%cscnv(info,type=afmt,dupl=dupl, mold=mold)
else if (a%is_upd()) then
call a%asb(mold=mold)

@ -70,6 +70,10 @@ subroutine psb_sspins(nz,ia,ja,val,a,desc_a,info,rebuild,local)
integer(psb_ipk_), parameter :: relocsz=200
logical :: rebuild_, local_
integer(psb_ipk_), allocatable :: ila(:),jla(:)
integer(psb_ipk_) :: i,k
integer(psb_lpk_) :: nnl
integer(psb_lpk_), allocatable :: lila(:),ljla(:)
real(psb_spk_), allocatable :: lval(:)
character(len=20) :: name
info = psb_success_
@ -147,6 +151,30 @@ subroutine psb_sspins(nz,ia,ja,val,a,desc_a,info,rebuild,local)
call psb_errpush(info,name,a_err='a%csput')
goto 9999
end if
select case(a%remote_build)
case (psb_matbld_noremote_)
! Do nothing
case (psb_matbld_remote_)
nnl = count(ila(1:nz)<0)
allocate(lila(nnl),ljla(nnl),lval(nnl))
k = 0
do i=1,nz
if (ila(i)<0) then
k=k+1
lila(k) = ia(k)
ljla(k) = ja(k)
lval(k) = val(k)
end if
end do
if (k /= nnl) write(0,*) name,' Wrong conversion?',k,nnl
call a%rmta%csput(nnl,lila,ljla,lval,1_psb_lpk_,desc_a%get_global_rows(),&
& 1_psb_lpk_,desc_a%get_global_rows(),info)
case default
write(0,*) name,' Ignoring wrong value for %remote_build'
end select
else
info = psb_err_invalid_a_and_cd_state_
call psb_errpush(info,name)

@ -41,7 +41,7 @@
! nnz - integer(optional). The number of nonzeroes in the matrix.
! (local, user estimate)
!
subroutine psb_zspalloc(a, desc_a, info, nnz)
subroutine psb_zspalloc(a, desc_a, info, nnz, bldmode)
use psb_base_mod, psb_protect_name => psb_zspalloc
implicit none
@ -49,13 +49,13 @@ subroutine psb_zspalloc(a, desc_a, info, nnz)
type(psb_desc_type), intent(in) :: desc_a
type(psb_zspmat_type), intent(inout) :: a
integer(psb_ipk_), intent(out) :: info
integer(psb_ipk_), optional, intent(in) :: nnz
integer(psb_ipk_), optional, intent(in) :: nnz, bldmode
!locals
type(psb_ctxt_type) :: ctxt
integer(psb_ipk_) :: np, me, err_act
integer(psb_ipk_) :: loc_row,loc_col, nnz_, dectype
integer(psb_lpk_) :: m, n
integer(psb_ipk_) :: loc_row,loc_col, nnz_, dectype, bldmode_
integer(psb_lpk_) :: m, n, nnzrmt_
integer(psb_ipk_) :: debug_level, debug_unit
character(len=20) :: name
@ -96,7 +96,12 @@ subroutine psb_zspalloc(a, desc_a, info, nnz)
else
nnz_ = max(1,5*loc_row)
endif
if (present(bldmode)) then
bldmode_ = bldmode
else
bldmode_ = psb_matbld_noremote_
end if
if (debug_level >= psb_debug_ext_) &
& write(debug_unit,*) me,' ',trim(name), &
& ':allocating size:',loc_row,loc_col,nnz_
@ -109,6 +114,25 @@ subroutine psb_zspalloc(a, desc_a, info, nnz)
goto 9999
end if
write(0,*) name,'Setting a%remote_build ',&
& bldmode_,psb_matbld_noremote_,psb_matbld_remote_
a%remote_build = bldmode_
select case(a%remote_build)
case (psb_matbld_noremote_)
! nothing needed
write(0,*) name,' matbld_noremote_ nothing needed'
case (psb_matbld_remote_)
write(0,*) name,' matbld_remote_ start '
allocate(psb_lz_coo_sparse_mat :: a%rmta)
nnzrmt_ = max(100,(nnz_/100))
call a%rmta%allocate(m,n,nnzrmt_)
case default
write(0,*) name,'Invalid value for remote_build '
a%remote_build = psb_matbld_noremote_
end select
if (debug_level >= psb_debug_ext_) &
& write(debug_unit,*) me,' ',trim(name),': ', &
& desc_a%get_dectype(),psb_desc_bld_

@ -111,7 +111,15 @@ subroutine psb_zspasb(a,desc_a, info, afmt, upd, dupl, mold)
call a%set_ncols(n_col)
end if
if (a%is_bld()) then
if (a%is_bld()) then
select case(a%remote_build)
case (psb_matbld_noremote_)
! nothing needed
case (psb_matbld_remote_)
write(0,*) 'Need to implement data movement! '
end select
call a%cscnv(info,type=afmt,dupl=dupl, mold=mold)
else if (a%is_upd()) then
call a%asb(mold=mold)

@ -70,6 +70,10 @@ subroutine psb_zspins(nz,ia,ja,val,a,desc_a,info,rebuild,local)
integer(psb_ipk_), parameter :: relocsz=200
logical :: rebuild_, local_
integer(psb_ipk_), allocatable :: ila(:),jla(:)
integer(psb_ipk_) :: i,k
integer(psb_lpk_) :: nnl
integer(psb_lpk_), allocatable :: lila(:),ljla(:)
complex(psb_dpk_), allocatable :: lval(:)
character(len=20) :: name
info = psb_success_
@ -147,6 +151,30 @@ subroutine psb_zspins(nz,ia,ja,val,a,desc_a,info,rebuild,local)
call psb_errpush(info,name,a_err='a%csput')
goto 9999
end if
select case(a%remote_build)
case (psb_matbld_noremote_)
! Do nothing
case (psb_matbld_remote_)
nnl = count(ila(1:nz)<0)
allocate(lila(nnl),ljla(nnl),lval(nnl))
k = 0
do i=1,nz
if (ila(i)<0) then
k=k+1
lila(k) = ia(k)
ljla(k) = ja(k)
lval(k) = val(k)
end if
end do
if (k /= nnl) write(0,*) name,' Wrong conversion?',k,nnl
call a%rmta%csput(nnl,lila,ljla,lval,1_psb_lpk_,desc_a%get_global_rows(),&
& 1_psb_lpk_,desc_a%get_global_rows(),info)
case default
write(0,*) name,' Ignoring wrong value for %remote_build'
end select
else
info = psb_err_invalid_a_and_cd_state_
call psb_errpush(info,name)

@ -204,38 +204,5 @@ contains
end function psb_c_cvect_set_vect
function psb_c_g2l(cdh,gindex,cowned) bind(c) result(lindex)
use psb_base_mod
implicit none
integer(psb_c_lpk_), value :: gindex
logical(c_bool), value :: cowned
type(psb_c_descriptor) :: cdh
integer(psb_c_ipk_) :: lindex
type(psb_desc_type), pointer :: descp
integer(psb_ipk_) :: info, localindex, ixb, iam, np
logical :: owned
ixb = psb_c_get_index_base()
owned = cowned
lindex = -1
if (c_associated(cdh%item)) then
call c_f_pointer(cdh%item,descp)
else
return
end if
call psb_info(descp%get_context(),iam,np)
if (ixb == 1) then
call descp%indxmap%g2l(gindex,localindex,info,owned=owned)
lindex = localindex
else
call descp%indxmap%g2l(gindex+(1-ixb),localindex,info,owned=owned)
lindex = localindex-(1-ixb)
endif
end function psb_c_g2l
end module psb_c_serial_cbind_mod

@ -429,7 +429,8 @@ contains
end select
if (info == psb_success_) call psb_spall(a,desc_a,info,nnz=nnz)
if (info == psb_success_) call psb_spall(a,desc_a,info,nnz=nnz, &
& bldmode=psb_matbld_remote_)
! define rhs from boundary conditions; also build initial guess
if (info == psb_success_) call psb_geall(xv,desc_a,info)
if (info == psb_success_) call psb_geall(bv,desc_a,info)

@ -429,7 +429,8 @@ contains
end select
if (info == psb_success_) call psb_spall(a,desc_a,info,nnz=nnz)
if (info == psb_success_) call psb_spall(a,desc_a,info,nnz=nnz, &
& bldmode=psb_matbld_remote_)
! define rhs from boundary conditions; also build initial guess
if (info == psb_success_) call psb_geall(xv,desc_a,info)
if (info == psb_success_) call psb_geall(bv,desc_a,info)

Loading…
Cancel
Save