Maintenance fixes

maint-3.9.0
sfilippone 1 week ago
parent 4e7d899fb6
commit 914dd33e4a

@ -129,6 +129,7 @@ module psb_c_base_vect_mod
procedure, pass(x) :: set_vect => c_base_set_vect
generic, public :: set => set_vect, set_scal
procedure, pass(x) :: get_entry=> c_base_get_entry
procedure, pass(x) :: set_entry=> c_base_set_entry
!
! Gather/scatter. These are needed for MPI interfacing.
! May have to be reworked.
@ -903,14 +904,30 @@ contains
!
function c_base_get_entry(x, index) result(res)
implicit none
class(psb_c_base_vect_type), intent(in) :: x
class(psb_c_base_vect_type), intent(inout) :: x
integer(psb_ipk_), intent(in) :: index
complex(psb_spk_) :: res
res = 0
if (allocated(x%v)) res = x%v(index)
res = czero
if (allocated(x%v)) then
if (x%is_dev()) call x%sync()
res = x%v(index)
end if
end function c_base_get_entry
subroutine c_base_set_entry(x, index, val)
implicit none
class(psb_c_base_vect_type), intent(inout) :: x
integer(psb_ipk_), intent(in) :: index
complex(psb_spk_) :: val
if (allocated(x%v)) then
if (x%is_dev()) call x%sync()
x%v(index) = val
call x%set_host()
end if
end subroutine c_base_set_entry
!
! Overwrite with absolute value

@ -93,6 +93,7 @@ module psb_c_vect_mod
procedure, pass(x) :: set_sync => c_vect_set_sync
procedure, pass(x) :: get_entry => c_vect_get_entry
procedure, pass(x) :: set_entry => c_vect_set_entry
procedure, pass(x) :: dot_v => c_vect_dot_v
procedure, pass(x) :: dot_a => c_vect_dot_a
@ -680,13 +681,22 @@ contains
function c_vect_get_entry(x,index) result(res)
implicit none
class(psb_c_vect_type), intent(in) :: x
class(psb_c_vect_type), intent(inout) :: x
integer(psb_ipk_), intent(in) :: index
complex(psb_spk_) :: res
res = 0
res = czero
if (allocated(x%v)) res = x%v%get_entry(index)
end function c_vect_get_entry
subroutine c_vect_set_entry(x,index,val)
implicit none
class(psb_c_vect_type), intent(inout) :: x
integer(psb_ipk_), intent(in) :: index
complex(psb_spk_) :: val
if (allocated(x%v)) call x%v%set_entry(index,val)
end subroutine c_vect_set_entry
function c_vect_dot_v(n,x,y) result(res)
implicit none
class(psb_c_vect_type), intent(inout) :: x, y

@ -129,6 +129,7 @@ module psb_d_base_vect_mod
procedure, pass(x) :: set_vect => d_base_set_vect
generic, public :: set => set_vect, set_scal
procedure, pass(x) :: get_entry=> d_base_get_entry
procedure, pass(x) :: set_entry=> d_base_set_entry
!
! Gather/scatter. These are needed for MPI interfacing.
! May have to be reworked.
@ -910,14 +911,30 @@ contains
!
function d_base_get_entry(x, index) result(res)
implicit none
class(psb_d_base_vect_type), intent(in) :: x
class(psb_d_base_vect_type), intent(inout) :: x
integer(psb_ipk_), intent(in) :: index
real(psb_dpk_) :: res
res = 0
if (allocated(x%v)) res = x%v(index)
res = dzero
if (allocated(x%v)) then
if (x%is_dev()) call x%sync()
res = x%v(index)
end if
end function d_base_get_entry
subroutine d_base_set_entry(x, index, val)
implicit none
class(psb_d_base_vect_type), intent(inout) :: x
integer(psb_ipk_), intent(in) :: index
real(psb_dpk_) :: val
if (allocated(x%v)) then
if (x%is_dev()) call x%sync()
x%v(index) = val
call x%set_host()
end if
end subroutine d_base_set_entry
!
! Overwrite with absolute value
@ -1810,8 +1827,8 @@ contains
integer(psb_ipk_) :: i
if (x%is_dev()) call x%sync()
#if defined(PSB_OPENMP)
res = HUGE(done)
#if defined(PSB_OPENMP)
!$omp parallel do private(i) reduction(min: res)
do i=1, n
res = min(res,abs(x%v(i)))

@ -93,6 +93,7 @@ module psb_d_vect_mod
procedure, pass(x) :: set_sync => d_vect_set_sync
procedure, pass(x) :: get_entry => d_vect_get_entry
procedure, pass(x) :: set_entry => d_vect_set_entry
procedure, pass(x) :: dot_v => d_vect_dot_v
procedure, pass(x) :: dot_a => d_vect_dot_a
@ -687,13 +688,22 @@ contains
function d_vect_get_entry(x,index) result(res)
implicit none
class(psb_d_vect_type), intent(in) :: x
class(psb_d_vect_type), intent(inout) :: x
integer(psb_ipk_), intent(in) :: index
real(psb_dpk_) :: res
res = 0
res = dzero
if (allocated(x%v)) res = x%v%get_entry(index)
end function d_vect_get_entry
subroutine d_vect_set_entry(x,index,val)
implicit none
class(psb_d_vect_type), intent(inout) :: x
integer(psb_ipk_), intent(in) :: index
real(psb_dpk_) :: val
if (allocated(x%v)) call x%v%set_entry(index,val)
end subroutine d_vect_set_entry
function d_vect_dot_v(n,x,y) result(res)
implicit none
class(psb_d_vect_type), intent(inout) :: x, y
@ -1255,7 +1265,7 @@ contains
if (allocated(x%v)) then
res = x%v%minreal(n)
else
res = dzero
res = HUGE(dzero)
end if
end function d_vect_min

@ -129,6 +129,7 @@ module psb_s_base_vect_mod
procedure, pass(x) :: set_vect => s_base_set_vect
generic, public :: set => set_vect, set_scal
procedure, pass(x) :: get_entry=> s_base_get_entry
procedure, pass(x) :: set_entry=> s_base_set_entry
!
! Gather/scatter. These are needed for MPI interfacing.
! May have to be reworked.
@ -910,14 +911,30 @@ contains
!
function s_base_get_entry(x, index) result(res)
implicit none
class(psb_s_base_vect_type), intent(in) :: x
class(psb_s_base_vect_type), intent(inout) :: x
integer(psb_ipk_), intent(in) :: index
real(psb_spk_) :: res
res = 0
if (allocated(x%v)) res = x%v(index)
res = szero
if (allocated(x%v)) then
if (x%is_dev()) call x%sync()
res = x%v(index)
end if
end function s_base_get_entry
subroutine s_base_set_entry(x, index, val)
implicit none
class(psb_s_base_vect_type), intent(inout) :: x
integer(psb_ipk_), intent(in) :: index
real(psb_spk_) :: val
if (allocated(x%v)) then
if (x%is_dev()) call x%sync()
x%v(index) = val
call x%set_host()
end if
end subroutine s_base_set_entry
!
! Overwrite with absolute value
@ -1810,8 +1827,8 @@ contains
integer(psb_ipk_) :: i
if (x%is_dev()) call x%sync()
#if defined(PSB_OPENMP)
res = HUGE(sone)
#if defined(PSB_OPENMP)
!$omp parallel do private(i) reduction(min: res)
do i=1, n
res = min(res,abs(x%v(i)))

@ -93,6 +93,7 @@ module psb_s_vect_mod
procedure, pass(x) :: set_sync => s_vect_set_sync
procedure, pass(x) :: get_entry => s_vect_get_entry
procedure, pass(x) :: set_entry => s_vect_set_entry
procedure, pass(x) :: dot_v => s_vect_dot_v
procedure, pass(x) :: dot_a => s_vect_dot_a
@ -687,13 +688,22 @@ contains
function s_vect_get_entry(x,index) result(res)
implicit none
class(psb_s_vect_type), intent(in) :: x
class(psb_s_vect_type), intent(inout) :: x
integer(psb_ipk_), intent(in) :: index
real(psb_spk_) :: res
res = 0
res = szero
if (allocated(x%v)) res = x%v%get_entry(index)
end function s_vect_get_entry
subroutine s_vect_set_entry(x,index,val)
implicit none
class(psb_s_vect_type), intent(inout) :: x
integer(psb_ipk_), intent(in) :: index
real(psb_spk_) :: val
if (allocated(x%v)) call x%v%set_entry(index,val)
end subroutine s_vect_set_entry
function s_vect_dot_v(n,x,y) result(res)
implicit none
class(psb_s_vect_type), intent(inout) :: x, y
@ -1255,7 +1265,7 @@ contains
if (allocated(x%v)) then
res = x%v%minreal(n)
else
res = szero
res = HUGE(szero)
end if
end function s_vect_min

@ -129,6 +129,7 @@ module psb_z_base_vect_mod
procedure, pass(x) :: set_vect => z_base_set_vect
generic, public :: set => set_vect, set_scal
procedure, pass(x) :: get_entry=> z_base_get_entry
procedure, pass(x) :: set_entry=> z_base_set_entry
!
! Gather/scatter. These are needed for MPI interfacing.
! May have to be reworked.
@ -903,14 +904,30 @@ contains
!
function z_base_get_entry(x, index) result(res)
implicit none
class(psb_z_base_vect_type), intent(in) :: x
class(psb_z_base_vect_type), intent(inout) :: x
integer(psb_ipk_), intent(in) :: index
complex(psb_dpk_) :: res
res = 0
if (allocated(x%v)) res = x%v(index)
res = zzero
if (allocated(x%v)) then
if (x%is_dev()) call x%sync()
res = x%v(index)
end if
end function z_base_get_entry
subroutine z_base_set_entry(x, index, val)
implicit none
class(psb_z_base_vect_type), intent(inout) :: x
integer(psb_ipk_), intent(in) :: index
complex(psb_dpk_) :: val
if (allocated(x%v)) then
if (x%is_dev()) call x%sync()
x%v(index) = val
call x%set_host()
end if
end subroutine z_base_set_entry
!
! Overwrite with absolute value

@ -93,6 +93,7 @@ module psb_z_vect_mod
procedure, pass(x) :: set_sync => z_vect_set_sync
procedure, pass(x) :: get_entry => z_vect_get_entry
procedure, pass(x) :: set_entry => z_vect_set_entry
procedure, pass(x) :: dot_v => z_vect_dot_v
procedure, pass(x) :: dot_a => z_vect_dot_a
@ -680,13 +681,22 @@ contains
function z_vect_get_entry(x,index) result(res)
implicit none
class(psb_z_vect_type), intent(in) :: x
class(psb_z_vect_type), intent(inout) :: x
integer(psb_ipk_), intent(in) :: index
complex(psb_dpk_) :: res
res = 0
res = zzero
if (allocated(x%v)) res = x%v%get_entry(index)
end function z_vect_get_entry
subroutine z_vect_set_entry(x,index,val)
implicit none
class(psb_z_vect_type), intent(inout) :: x
integer(psb_ipk_), intent(in) :: index
complex(psb_dpk_) :: val
if (allocated(x%v)) call x%v%set_entry(index,val)
end subroutine z_vect_set_entry
function z_vect_dot_v(n,x,y) result(res)
implicit none
class(psb_z_vect_type), intent(inout) :: x, y

@ -443,6 +443,17 @@ Module psb_c_tools_mod
end function
end interface
interface psb_setelem
subroutine psb_c_setelem(index,val,x,desc_a,info)
import
type(psb_c_vect_type), intent(inout) :: x
integer(psb_lpk_), intent(in) :: index
type(psb_desc_type), intent(inout) :: desc_a
integer(psb_ipk_), intent(out) :: info
complex(psb_spk_) ::val
end subroutine psb_c_setelem
end interface
interface psb_remap
subroutine psb_c_remap(np_remap, desc_in, a_in, &
& ipd, isrc, nrsrc, naggr, desc_out, a_out, info)

@ -443,6 +443,17 @@ Module psb_d_tools_mod
end function
end interface
interface psb_setelem
subroutine psb_d_setelem(index,val,x,desc_a,info)
import
type(psb_d_vect_type), intent(inout) :: x
integer(psb_lpk_), intent(in) :: index
type(psb_desc_type), intent(inout) :: desc_a
integer(psb_ipk_), intent(out) :: info
real(psb_dpk_) ::val
end subroutine psb_d_setelem
end interface
interface psb_remap
subroutine psb_d_remap(np_remap, desc_in, a_in, &
& ipd, isrc, nrsrc, naggr, desc_out, a_out, info)

@ -443,6 +443,17 @@ Module psb_s_tools_mod
end function
end interface
interface psb_setelem
subroutine psb_s_setelem(index,val,x,desc_a,info)
import
type(psb_s_vect_type), intent(inout) :: x
integer(psb_lpk_), intent(in) :: index
type(psb_desc_type), intent(inout) :: desc_a
integer(psb_ipk_), intent(out) :: info
real(psb_spk_) ::val
end subroutine psb_s_setelem
end interface
interface psb_remap
subroutine psb_s_remap(np_remap, desc_in, a_in, &
& ipd, isrc, nrsrc, naggr, desc_out, a_out, info)

@ -443,6 +443,17 @@ Module psb_z_tools_mod
end function
end interface
interface psb_setelem
subroutine psb_z_setelem(index,val,x,desc_a,info)
import
type(psb_z_vect_type), intent(inout) :: x
integer(psb_lpk_), intent(in) :: index
type(psb_desc_type), intent(inout) :: desc_a
integer(psb_ipk_), intent(out) :: info
complex(psb_dpk_) ::val
end subroutine psb_z_setelem
end interface
interface psb_remap
subroutine psb_z_remap(np_remap, desc_in, a_in, &
& ipd, isrc, nrsrc, naggr, desc_out, a_out, info)

@ -180,6 +180,30 @@ contains
end function psb_c_cvect_set_scal
function psb_c_cvect_set_scal_bound(x,val,ifirst,ilast) bind(c) result(info)
use psb_base_mod
implicit none
type(psb_c_cvector) :: x
type(psb_c_vect_type), pointer :: xp
integer(psb_c_ipk_) :: info
integer(psb_c_ipk_), value :: ifirst, ilast
complex(c_float_complex) :: val
info = -1;
if (c_associated(x%item)) then
call c_f_pointer(x%item,xp)
else
return
end if
call xp%set(val,first=ifirst,last=ilast)
info = 0
end function psb_c_cvect_set_scal_bound
function psb_c_cvect_set_vect(x,val,n) bind(c) result(info)
use psb_base_mod
implicit none
@ -204,5 +228,49 @@ contains
end function psb_c_cvect_set_vect
function psb_c_cvect_set_entry(x,index,val) bind(c) result(info)
use psb_base_mod
implicit none
type(psb_c_cvector) :: x
type(psb_c_vect_type), pointer :: xp
integer(psb_c_ipk_) :: info
integer(psb_c_ipk_), value :: index
complex(c_float_complex), value :: val
integer(psb_c_ipk_) :: ixb
info = -1;
if (c_associated(x%item)) then
call c_f_pointer(x%item,xp)
else
return
end if
ixb = psb_c_get_index_base()
call xp%set_entry((index+(1-ixb)),val)
info = 0
end function psb_c_cvect_set_entry
function psb_c_cvect_get_entry(x,index) bind(c) result(res)
use psb_base_mod
implicit none
type(psb_c_cvector) :: x
type(psb_c_vect_type), pointer :: xp
integer(psb_c_ipk_), value :: index
complex(c_float_complex) :: res
integer(psb_c_ipk_) :: ixb
if (c_associated(x%item)) then
call c_f_pointer(x%item,xp)
else
return
end if
ixb = psb_c_get_index_base()
res = xp%get_entry((index+(1-ixb)))
end function psb_c_cvect_get_entry
end module psb_c_serial_cbind_mod

@ -452,4 +452,40 @@ contains
end function psb_c_cgetelem
function psb_c_csetelem(index,val,xh,cdh) bind(c) result(res)
implicit none
type(psb_c_cvector) :: xh
integer(psb_c_lpk_), value :: index
type(psb_c_descriptor) :: cdh
complex(c_float_complex), value :: val
integer(psb_c_ipk_) :: res
type(psb_c_vect_type), pointer :: xp
type(psb_desc_type), pointer :: descp
integer(psb_c_ipk_) :: info, ixb
res = -1
if (c_associated(cdh%item)) then
call c_f_pointer(cdh%item,descp)
else
return
end if
if (c_associated(xh%item)) then
call c_f_pointer(xh%item,xp)
else
return
end if
ixb = psb_c_get_index_base()
if (ixb == 1) then
call psb_setelem(index,val,xp,descp,info)
else
call psb_setelem(index+(1-ixb),val,xp,descp,info)
end if
res=info
return
end function psb_c_csetelem
end module psb_c_tools_cbind_mod

@ -180,6 +180,30 @@ contains
end function psb_c_dvect_set_scal
function psb_c_dvect_set_scal_bound(x,val,ifirst,ilast) bind(c) result(info)
use psb_base_mod
implicit none
type(psb_c_dvector) :: x
type(psb_d_vect_type), pointer :: xp
integer(psb_c_ipk_) :: info
integer(psb_c_ipk_), value :: ifirst, ilast
real(c_double) :: val
info = -1;
if (c_associated(x%item)) then
call c_f_pointer(x%item,xp)
else
return
end if
call xp%set(val,first=ifirst,last=ilast)
info = 0
end function psb_c_dvect_set_scal_bound
function psb_c_dvect_set_vect(x,val,n) bind(c) result(info)
use psb_base_mod
implicit none
@ -204,5 +228,49 @@ contains
end function psb_c_dvect_set_vect
function psb_c_dvect_set_entry(x,index,val) bind(c) result(info)
use psb_base_mod
implicit none
type(psb_c_dvector) :: x
type(psb_d_vect_type), pointer :: xp
integer(psb_c_ipk_) :: info
integer(psb_c_ipk_), value :: index
real(c_double), value :: val
integer(psb_c_ipk_) :: ixb
info = -1;
if (c_associated(x%item)) then
call c_f_pointer(x%item,xp)
else
return
end if
ixb = psb_c_get_index_base()
call xp%set_entry((index+(1-ixb)),val)
info = 0
end function psb_c_dvect_set_entry
function psb_c_dvect_get_entry(x,index) bind(c) result(res)
use psb_base_mod
implicit none
type(psb_c_dvector) :: x
type(psb_d_vect_type), pointer :: xp
integer(psb_c_ipk_), value :: index
real(c_double) :: res
integer(psb_c_ipk_) :: ixb
if (c_associated(x%item)) then
call c_f_pointer(x%item,xp)
else
return
end if
ixb = psb_c_get_index_base()
res = xp%get_entry((index+(1-ixb)))
end function psb_c_dvect_get_entry
end module psb_d_serial_cbind_mod

@ -452,4 +452,40 @@ contains
end function psb_c_dgetelem
function psb_c_dsetelem(index,val,xh,cdh) bind(c) result(res)
implicit none
type(psb_c_dvector) :: xh
integer(psb_c_lpk_), value :: index
type(psb_c_descriptor) :: cdh
real(c_double), value :: val
integer(psb_c_ipk_) :: res
type(psb_d_vect_type), pointer :: xp
type(psb_desc_type), pointer :: descp
integer(psb_c_ipk_) :: info, ixb
res = -1
if (c_associated(cdh%item)) then
call c_f_pointer(cdh%item,descp)
else
return
end if
if (c_associated(xh%item)) then
call c_f_pointer(xh%item,xp)
else
return
end if
ixb = psb_c_get_index_base()
if (ixb == 1) then
call psb_setelem(index,val,xp,descp,info)
else
call psb_setelem(index+(1-ixb),val,xp,descp,info)
end if
res=info
return
end function psb_c_dsetelem
end module psb_d_tools_cbind_mod

@ -180,6 +180,30 @@ contains
end function psb_c_svect_set_scal
function psb_c_svect_set_scal_bound(x,val,ifirst,ilast) bind(c) result(info)
use psb_base_mod
implicit none
type(psb_c_svector) :: x
type(psb_s_vect_type), pointer :: xp
integer(psb_c_ipk_) :: info
integer(psb_c_ipk_), value :: ifirst, ilast
real(c_float) :: val
info = -1;
if (c_associated(x%item)) then
call c_f_pointer(x%item,xp)
else
return
end if
call xp%set(val,first=ifirst,last=ilast)
info = 0
end function psb_c_svect_set_scal_bound
function psb_c_svect_set_vect(x,val,n) bind(c) result(info)
use psb_base_mod
implicit none
@ -204,5 +228,49 @@ contains
end function psb_c_svect_set_vect
function psb_c_svect_set_entry(x,index,val) bind(c) result(info)
use psb_base_mod
implicit none
type(psb_c_svector) :: x
type(psb_s_vect_type), pointer :: xp
integer(psb_c_ipk_) :: info
integer(psb_c_ipk_), value :: index
real(c_float), value :: val
integer(psb_c_ipk_) :: ixb
info = -1;
if (c_associated(x%item)) then
call c_f_pointer(x%item,xp)
else
return
end if
ixb = psb_c_get_index_base()
call xp%set_entry((index+(1-ixb)),val)
info = 0
end function psb_c_svect_set_entry
function psb_c_svect_get_entry(x,index) bind(c) result(res)
use psb_base_mod
implicit none
type(psb_c_svector) :: x
type(psb_s_vect_type), pointer :: xp
integer(psb_c_ipk_), value :: index
real(c_float) :: res
integer(psb_c_ipk_) :: ixb
if (c_associated(x%item)) then
call c_f_pointer(x%item,xp)
else
return
end if
ixb = psb_c_get_index_base()
res = xp%get_entry((index+(1-ixb)))
end function psb_c_svect_get_entry
end module psb_s_serial_cbind_mod

@ -452,4 +452,40 @@ contains
end function psb_c_sgetelem
function psb_c_ssetelem(index,val,xh,cdh) bind(c) result(res)
implicit none
type(psb_c_svector) :: xh
integer(psb_c_lpk_), value :: index
type(psb_c_descriptor) :: cdh
real(c_float), value :: val
integer(psb_c_ipk_) :: res
type(psb_s_vect_type), pointer :: xp
type(psb_desc_type), pointer :: descp
integer(psb_c_ipk_) :: info, ixb
res = -1
if (c_associated(cdh%item)) then
call c_f_pointer(cdh%item,descp)
else
return
end if
if (c_associated(xh%item)) then
call c_f_pointer(xh%item,xp)
else
return
end if
ixb = psb_c_get_index_base()
if (ixb == 1) then
call psb_setelem(index,val,xp,descp,info)
else
call psb_setelem(index+(1-ixb),val,xp,descp,info)
end if
res=info
return
end function psb_c_ssetelem
end module psb_s_tools_cbind_mod

@ -180,6 +180,30 @@ contains
end function psb_c_zvect_set_scal
function psb_c_zvect_set_scal_bound(x,val,ifirst,ilast) bind(c) result(info)
use psb_base_mod
implicit none
type(psb_c_zvector) :: x
type(psb_z_vect_type), pointer :: xp
integer(psb_c_ipk_) :: info
integer(psb_c_ipk_), value :: ifirst, ilast
complex(c_double_complex) :: val
info = -1;
if (c_associated(x%item)) then
call c_f_pointer(x%item,xp)
else
return
end if
call xp%set(val,first=ifirst,last=ilast)
info = 0
end function psb_c_zvect_set_scal_bound
function psb_c_zvect_set_vect(x,val,n) bind(c) result(info)
use psb_base_mod
implicit none
@ -204,5 +228,49 @@ contains
end function psb_c_zvect_set_vect
function psb_c_zvect_set_entry(x,index,val) bind(c) result(info)
use psb_base_mod
implicit none
type(psb_c_zvector) :: x
type(psb_z_vect_type), pointer :: xp
integer(psb_c_ipk_) :: info
integer(psb_c_ipk_), value :: index
complex(c_double_complex), value :: val
integer(psb_c_ipk_) :: ixb
info = -1;
if (c_associated(x%item)) then
call c_f_pointer(x%item,xp)
else
return
end if
ixb = psb_c_get_index_base()
call xp%set_entry((index+(1-ixb)),val)
info = 0
end function psb_c_zvect_set_entry
function psb_c_zvect_get_entry(x,index) bind(c) result(res)
use psb_base_mod
implicit none
type(psb_c_zvector) :: x
type(psb_z_vect_type), pointer :: xp
integer(psb_c_ipk_), value :: index
complex(c_double_complex) :: res
integer(psb_c_ipk_) :: ixb
if (c_associated(x%item)) then
call c_f_pointer(x%item,xp)
else
return
end if
ixb = psb_c_get_index_base()
res = xp%get_entry((index+(1-ixb)))
end function psb_c_zvect_get_entry
end module psb_z_serial_cbind_mod

@ -452,4 +452,40 @@ contains
end function psb_c_zgetelem
function psb_c_zsetelem(index,val,xh,cdh) bind(c) result(res)
implicit none
type(psb_c_zvector) :: xh
integer(psb_c_lpk_), value :: index
type(psb_c_descriptor) :: cdh
complex(c_double_complex), value :: val
integer(psb_c_ipk_) :: res
type(psb_z_vect_type), pointer :: xp
type(psb_desc_type), pointer :: descp
integer(psb_c_ipk_) :: info, ixb
res = -1
if (c_associated(cdh%item)) then
call c_f_pointer(cdh%item,descp)
else
return
end if
if (c_associated(xh%item)) then
call c_f_pointer(xh%item,xp)
else
return
end if
ixb = psb_c_get_index_base()
if (ixb == 1) then
call psb_setelem(index,val,xp,descp,info)
else
call psb_setelem(index+(1-ixb),val,xp,descp,info)
end if
res=info
return
end function psb_c_zsetelem
end module psb_z_tools_cbind_mod

@ -230,7 +230,7 @@ contains
f_ => d_null_func_2d
end if
deltah = done/(idim+1)
deltah = done/(idim+2)
sqdeltah = deltah*deltah
deltah2 = (2*done)* deltah
@ -467,8 +467,8 @@ contains
! compute gridpoint coordinates
call idx2ijk(ix,iy,glob_row,idim,idim)
! x, y coordinates
x = (ix-1)*deltah
y = (iy-1)*deltah
x = (ix)*deltah
y = (iy)*deltah
zt(k) = f_(x,y)
! internal point: build discretization

@ -246,7 +246,7 @@ contains
f_ => d_null_func_3d
end if
deltah = done/(idim+1)
deltah = done/(idim+2)
sqdeltah = deltah*deltah
deltah2 = (2*done)* deltah
@ -496,9 +496,9 @@ contains
! compute gridpoint coordinates
call idx2ijk(ix,iy,iz,glob_row,idim,idim,idim)
! x, y, z coordinates
x = (ix-1)*deltah
y = (iy-1)*deltah
z = (iz-1)*deltah
x = (ix)*deltah
y = (iy)*deltah
z = (iz)*deltah
zt(k) = f_(x,y,z)
! internal point: build discretization
!

@ -230,7 +230,7 @@ contains
f_ => s_null_func_2d
end if
deltah = sone/(idim+1)
deltah = sone/(idim+2)
sqdeltah = deltah*deltah
deltah2 = (2*sone)* deltah
@ -467,8 +467,8 @@ contains
! compute gridpoint coordinates
call idx2ijk(ix,iy,glob_row,idim,idim)
! x, y coordinates
x = (ix-1)*deltah
y = (iy-1)*deltah
x = (ix)*deltah
y = (iy)*deltah
zt(k) = f_(x,y)
! internal point: build discretization

@ -246,7 +246,7 @@ contains
f_ => s_null_func_3d
end if
deltah = sone/(idim+1)
deltah = sone/(idim+2)
sqdeltah = deltah*deltah
deltah2 = (2*sone)* deltah
@ -496,9 +496,9 @@ contains
! compute gridpoint coordinates
call idx2ijk(ix,iy,iz,glob_row,idim,idim,idim)
! x, y, z coordinates
x = (ix-1)*deltah
y = (iy-1)*deltah
z = (iz-1)*deltah
x = (ix)*deltah
y = (iy)*deltah
z = (iz)*deltah
zt(k) = f_(x,y,z)
! internal point: build discretization
!

Loading…
Cancel
Save