The current version of test/omp seems to be working. To be completed

omp-threadsafe
sfilippone 2 years ago
parent c05b32c202
commit 0f1603a2e9

@ -1036,7 +1036,12 @@ contains
#endif
else if (.not.use_openmp) then
#ifdef OPENMP
! $ omp parallel
! $ omp critical
!write(0,*) 'In cnv: ',omp_get_num_threads()
#endif
isLoopValid = .true.
if (idxmap%is_bld()) then
if (present(lidx)) then
@ -1067,7 +1072,7 @@ contains
if (info /= psb_success_) then
call psb_errpush(psb_err_from_subroutine_ai_,name,&
&a_err='psb_ensure_size',i_err=(/info/))
goto 9999
isLoopValid = .false.
end if
idxmap%loc_to_glob(nxt) = ip
call idxmap%set_lc(ncol)
@ -1076,7 +1081,7 @@ contains
else
call psb_errpush(psb_err_from_subroutine_ai_,name,&
& a_err='SearchInsKeyVal',i_err=(/info/))
goto 9999
isLoopValid = .false.
end if
end if
idx(i) = lip
@ -1114,7 +1119,7 @@ contains
info=1
call psb_errpush(psb_err_from_subroutine_ai_,name,&
&a_err='psb_ensure_size',i_err=(/info/))
goto 9999
isLoopValid = .false.
end if
idxmap%loc_to_glob(nxt) = ip
call idxmap%set_lc(ncol)
@ -1123,7 +1128,7 @@ contains
else
call psb_errpush(psb_err_from_subroutine_ai_,name,&
& a_err='SearchInsKeyVal',i_err=(/info/))
goto 9999
isLoopValid = .false.
end if
end if
idx(i) = lip
@ -1160,7 +1165,7 @@ contains
info=1
call psb_errpush(psb_err_from_subroutine_ai_,name,&
& a_err='psb_ensure_size',i_err=(/info/))
goto 9999
isLoopValid = .false.
end if
idxmap%loc_to_glob(nxt) = ip
call idxmap%set_lc(ncol)
@ -1169,7 +1174,7 @@ contains
else
call psb_errpush(psb_err_from_subroutine_ai_,name,&
& a_err='SearchInsKeyVal',i_err=(/info/))
goto 9999
isLoopValid = .false.
end if
idx(i) = lip
info = psb_success_
@ -1205,7 +1210,8 @@ contains
ch_err='psb_ensure_size'
call psb_errpush(psb_err_from_subroutine_ai_,name,&
&a_err=ch_err,i_err=(/info,izero,izero,izero,izero/))
goto 9999
isLoopValid = .false.
end if
idxmap%loc_to_glob(nxt) = ip
call idxmap%set_lc(ncol)
@ -1215,7 +1221,7 @@ contains
ch_err='SearchInsKeyVal'
call psb_errpush(psb_err_from_subroutine_ai_,name,&
& a_err=ch_err,i_err=(/info,izero,izero,izero,izero/))
goto 9999
isLoopValid = .false.
end if
idx(i) = lip
info = psb_success_
@ -1229,6 +1235,12 @@ contains
idx = -1
info = -1
end if
#ifdef OPENMP
! $ omp end critical
! $ omp end parallel
#endif
if (.not. isLoopValid) goto 9999
end if
call psb_erractionrestore(err_act)
return

@ -2818,6 +2818,9 @@ subroutine psb_d_coo_csput_a(nz,ia,ja,val,a,imin,imax,jmin,jmax,info)
use psb_realloc_mod
use psb_sort_mod
use psb_d_base_mat_mod, psb_protect_name => psb_d_coo_csput_a
#if defined(OPENMP)
use omp_lib
#endif
implicit none
class(psb_d_coo_sparse_mat), intent(inout) :: a
@ -2829,7 +2832,7 @@ subroutine psb_d_coo_csput_a(nz,ia,ja,val,a,imin,imax,jmin,jmax,info)
integer(psb_ipk_) :: err_act
character(len=20) :: name='d_coo_csput_a_impl'
logical, parameter :: debug=.false.
integer(psb_ipk_) :: nza, i,j,k, nzl, isza, debug_level, debug_unit
integer(psb_ipk_) :: nza, i,j,k, nzl, isza, debug_level, debug_unit, nzaold
info = psb_success_
debug_unit = psb_get_debug_unit()
@ -2862,9 +2865,11 @@ subroutine psb_d_coo_csput_a(nz,ia,ja,val,a,imin,imax,jmin,jmax,info)
if (nz == 0) return
if (a%is_bld()) then
!$omp critical
nza = a%get_nzeros()
isza = a%get_size()
if (a%is_bld()) then
! Build phase. Must handle reallocations in a sensible way.
if (isza < (nza+nz)) then
call a%reallocate(max(nza+nz,int(1.5*isza)))
@ -2872,16 +2877,23 @@ subroutine psb_d_coo_csput_a(nz,ia,ja,val,a,imin,imax,jmin,jmax,info)
isza = a%get_size()
if (isza < (nza+nz)) then
info = psb_err_alloc_dealloc_; call psb_errpush(info,name)
goto 9999
else
nzaold = nza
nza = nza + nz
call a%set_nzeros(nza)
#if defined(OPENMP)
!write(0,*) 'From thread ',omp_get_thread_num(),nzaold,nz,nza,a%get_nzeros()
#endif
end if
call psb_inner_ins(nz,ia,ja,val,nza,a%ia,a%ja,a%val,isza,&
!$omp end critical
if (info /= 0) goto 9999
call psb_inner_ins(nz,ia,ja,val,nzaold,a%ia,a%ja,a%val,isza,&
& imin,imax,jmin,jmax,info)
call a%set_nzeros(nza)
call a%set_sorted(.false.)
else if (a%is_upd()) then
nza = a%get_nzeros()
isza = a%get_size()
if (a%is_dev()) call a%sync()
@ -2951,7 +2963,7 @@ contains
end do
! $OMP END PARALLEL DO
nza = nza + nz
!nza = nza + nz
#else
do i=1, nz
ir = ia(i)

@ -51,6 +51,9 @@
subroutine psb_dspins(nz,ia,ja,val,a,desc_a,info,rebuild,local)
use psb_base_mod, psb_protect_name => psb_dspins
use psi_mod
#if defined(OPENMP)
use omp_lib
#endif
implicit none
!....parameters...
@ -82,7 +85,9 @@ subroutine psb_dspins(nz,ia,ja,val,a,desc_a,info,rebuild,local)
ctxt = desc_a%get_context()
call psb_info(ctxt, me, np)
#if defined(OPENMP)
!write(0,*) name,omp_get_num_threads(),omp_get_thread_num()
#endif
if (nz < 0) then
info = 1111
call psb_errpush(info,name)
@ -131,15 +136,26 @@ subroutine psb_dspins(nz,ia,ja,val,a,desc_a,info,rebuild,local)
& a_err='allocate',i_err=(/info/))
goto 9999
end if
#if defined(OPENMP)
!$omp parallel private(ila,jla,nrow,ncol)
#endif
call desc_a%indxmap%g2l(ia(1:nz),ila(1:nz),info,owned=.true.)
#if defined(OPENMP)
!$omp critical
#endif
if (info == 0) call desc_a%indxmap%g2l_ins(ja(1:nz),jla(1:nz),info,&
& mask=(ila(1:nz)>0))
!!$ write(0,*) omp_get_thread_num(),'Check g2l: ',ila(1:nz),':',&
!!$ & jla(1:nz),':',ja(1:nz)
#if defined(OPENMP)
!$omp end critical
#endif
if (info /= psb_success_) then
call psb_errpush(psb_err_from_subroutine_ai_,name,&
& a_err='psb_cdins',i_err=(/info/))
goto 9999
!goto 9999
end if
nrow = desc_a%get_local_rows()
ncol = desc_a%get_local_cols()
@ -149,7 +165,7 @@ subroutine psb_dspins(nz,ia,ja,val,a,desc_a,info,rebuild,local)
if (info /= psb_success_) then
info=psb_err_from_subroutine_
call psb_errpush(info,name,a_err='a%csput')
goto 9999
!goto 9999
end if
if (a%is_remote_build()) then
@ -175,8 +191,13 @@ subroutine psb_dspins(nz,ia,ja,val,a,desc_a,info,rebuild,local)
else
info = psb_err_invalid_a_and_cd_state_
call psb_errpush(info,name)
goto 9999
!goto 9999
end if
#if defined(OPENMP)
!$omp end parallel
#endif
endif
else if (desc_a%is_asb()) then

@ -0,0 +1,56 @@
INSTALLDIR=../..
INCDIR=$(INSTALLDIR)/include
MODDIR=$(INSTALLDIR)/modules/
include $(INCDIR)/Make.inc.psblas
#
# Libraries used
LIBDIR=$(INSTALLDIR)/lib
PSBLAS_LIB= -L$(LIBDIR) -lpsb_util -lpsb_krylov -lpsb_prec -lpsb_base
LDLIBS=$(PSBLDLIBS)
#
# Compilers and such
#
CCOPT= -g
FINCLUDES=$(FMFLAG)$(MODDIR) $(FMFLAG).
EXEDIR=./runs
all: runsd psb_tomp #psb_d_pde3d psb_s_pde3d psb_d_pde2d psb_s_pde2d
runsd:
(if test ! -d runs ; then mkdir runs; fi)
psb_tomp: psb_tomp.o
$(FLINK) psb_tomp.o -o psb_tomp $(PSBLAS_LIB) $(LDLIBS)
/bin/mv psb_tomp $(EXEDIR)
psb_d_pde3d: psb_d_pde3d.o
$(FLINK) psb_d_pde3d.o -o psb_d_pde3d $(PSBLAS_LIB) $(LDLIBS)
/bin/mv psb_d_pde3d $(EXEDIR)
psb_s_pde3d: psb_s_pde3d.o
$(FLINK) psb_s_pde3d.o -o psb_s_pde3d $(PSBLAS_LIB) $(LDLIBS)
/bin/mv psb_s_pde3d $(EXEDIR)
psb_d_pde2d: psb_d_pde2d.o
$(FLINK) psb_d_pde2d.o -o psb_d_pde2d $(PSBLAS_LIB) $(LDLIBS)
/bin/mv psb_d_pde2d $(EXEDIR)
psb_s_pde2d: psb_s_pde2d.o
$(FLINK) psb_s_pde2d.o -o psb_s_pde2d $(PSBLAS_LIB) $(LDLIBS)
/bin/mv psb_s_pde2d $(EXEDIR)
clean:
/bin/rm -f psb_tomp.o psb_d_pde3d.o psb_s_pde3d.o psb_d_pde2d.o psb_s_pde2d.o *$(.mod) \
$(EXEDIR)/psb_d_pde3d $(EXEDIR)/psb_s_pde3d $(EXEDIR)/psb_d_pde2d $(EXEDIR)/psb_s_pde2d
verycleanlib:
(cd ../..; make veryclean)
lib:
(cd ../../; make library)

File diff suppressed because it is too large Load Diff

@ -799,6 +799,7 @@ program psb_d_pde3d
if(psb_errstatus_fatal()) goto 9999
name='pde3d90'
call psb_set_errverbosity(itwo)
call psb_cd_set_large_threshold(2000_psb_ipk_)
!
! Hello world
!

Loading…
Cancel
Save