psblas3-mcbind:

cbind/Makefile
 cbind/base/Makefile
 cbind/base/psb_c_serial_cbind_mod.F90
 cbind/base/psb_d_serial_cbind_mod.F90
 cbind/base/psb_s_serial_cbind_mod.F90
 cbind/base/psb_z_serial_cbind_mod.F90
 cbind/test/pargen/ppdec.c

Added print interface.
psblas3-mcbind
Salvatore Filippone 8 years ago
parent 49cddbb6fb
commit 855f557673

@ -13,7 +13,7 @@ lib: based precd
based:
cd base && $(MAKE) lib LIBNAME=$(LIBNAME)
precd:
precd: based
cd prec && $(MAKE) lib LIBNAME=$(LIBNAME)

@ -52,7 +52,7 @@ psb_base_cbind_mod.o: psb_cpenv_mod.o psb_objhandle_mod.o psb_base_tools_cbind_m
psb_s_comm_cbind_mod.o psb_d_comm_cbind_mod.o \
psb_c_comm_cbind_mod.o psb_z_comm_cbind_mod.o
psb_base_tools_cbind_mod.o: psb_objhandle_mod.o psb_base_string_cbind_mod.o
psb_base_tools_cbind_mod.o: psb_cpenv_mod.o psb_objhandle_mod.o psb_base_string_cbind_mod.o
psb_s_tools_cbind_mod.o psb_s_serial_cbind_mod.o \
psb_d_tools_cbind_mod.o psb_d_serial_cbind_mod.o \

@ -113,7 +113,32 @@ contains
end function psb_c_cmat_get_ncols
function psb_c_cmat_name_print(mh,name) bind(c) result(res)
use psb_base_mod
use psb_objhandle_mod
use psb_base_string_cbind_mod
implicit none
integer(psb_c_int) :: res
character(c_char) :: name(*)
type(psb_c_cspmat) :: mh
type(psb_cspmat_type), pointer :: ap
integer :: info
character(1024) :: fname
res = 0
if (c_associated(mh%item)) then
call c_f_pointer(mh%item,ap)
else
return
end if
call stringc2f(name,fname)
call ap%print(fname,head='PSBLAS Cbinding Interface')
end function psb_c_cmat_name_print
end module psb_c_serial_cbind_mod

@ -113,7 +113,32 @@ contains
end function psb_c_dmat_get_ncols
function psb_c_dmat_name_print(mh,name) bind(c) result(res)
use psb_base_mod
use psb_objhandle_mod
use psb_base_string_cbind_mod
implicit none
integer(psb_c_int) :: res
character(c_char) :: name(*)
type(psb_c_dspmat) :: mh
type(psb_dspmat_type), pointer :: ap
integer :: info
character(1024) :: fname
res = 0
if (c_associated(mh%item)) then
call c_f_pointer(mh%item,ap)
else
return
end if
call stringc2f(name,fname)
call ap%print(fname,head='PSBLAS Cbinding Interface')
end function psb_c_dmat_name_print
end module psb_d_serial_cbind_mod

@ -113,7 +113,32 @@ contains
end function psb_c_smat_get_ncols
function psb_c_smat_name_print(mh,name) bind(c) result(res)
use psb_base_mod
use psb_objhandle_mod
use psb_base_string_cbind_mod
implicit none
integer(psb_c_int) :: res
character(c_char) :: name(*)
type(psb_c_sspmat) :: mh
type(psb_sspmat_type), pointer :: ap
integer :: info
character(1024) :: fname
res = 0
if (c_associated(mh%item)) then
call c_f_pointer(mh%item,ap)
else
return
end if
call stringc2f(name,fname)
call ap%print(fname,head='PSBLAS Cbinding Interface')
end function psb_c_smat_name_print
end module psb_s_serial_cbind_mod

@ -113,7 +113,32 @@ contains
end function psb_c_zmat_get_ncols
function psb_c_zmat_name_print(mh,name) bind(c) result(res)
use psb_base_mod
use psb_objhandle_mod
use psb_base_string_cbind_mod
implicit none
integer(psb_c_int) :: res
character(c_char) :: name(*)
type(psb_c_zspmat) :: mh
type(psb_zspmat_type), pointer :: ap
integer :: info
character(1024) :: fname
res = 0
if (c_associated(mh%item)) then
call c_f_pointer(mh%item,ap)
else
return
end if
call stringc2f(name,fname)
call ap%print(fname,head='PSBLAS Cbinding Interface')
end function psb_c_zmat_name_print
end module psb_z_serial_cbind_mod

@ -214,7 +214,10 @@ int matgen(int ictxt, int ng,int idim,int vg[],psb_c_dspmat *ah,psb_c_descriptor
icol[el]=(x)*idim*idim+(y-1)*idim+(z);
el=el+1;
}
for (i=0; i<el; i++) irow[i]=glob_row;
for (i=0; i<el; i++) {
irow[i]=glob_row-1;
icol[i]--;
}
if ((ret=psb_c_dspins(el,irow,icol,val,ah,cdh))!=0)
fprintf(stderr,"From psb_c_dspins: %d\n",ret);
psb_c_dgeins(1,irow,zt,bh,cdh);
@ -329,7 +332,7 @@ int main(int argc, char *argv[])
fprintf(stderr,"Error during matrix build loop\n");
psb_c_abort(ictxt);
}
psb_c_dmat_name_print(ah,"cbindmat.mtx");
psb_c_barrier(ictxt);
/* Set up the preconditioner */
ph = psb_c_new_dprec();

Loading…
Cancel
Save