From d94b3dae5822a8e4cf027d30afd15ecd2db609fc Mon Sep 17 00:00:00 2001 From: Salvatore Filippone Date: Fri, 15 Mar 2024 06:43:16 -0400 Subject: [PATCH] Fix tracking of vector allocations. --- cuda/psb_c_cuda_vect_mod.F90 | 6 +++--- cuda/psb_d_cuda_vect_mod.F90 | 6 +++--- cuda/psb_i_cuda_vect_mod.F90 | 6 +++--- cuda/psb_s_cuda_vect_mod.F90 | 6 +++--- cuda/psb_z_cuda_vect_mod.F90 | 6 +++--- 5 files changed, 15 insertions(+), 15 deletions(-) diff --git a/cuda/psb_c_cuda_vect_mod.F90 b/cuda/psb_c_cuda_vect_mod.F90 index fb0b8dc7..bf6d2b87 100644 --- a/cuda/psb_c_cuda_vect_mod.F90 +++ b/cuda/psb_c_cuda_vect_mod.F90 @@ -729,7 +729,7 @@ contains if (c_associated(x%deviceVect)) then nd = getMultiVecDeviceSize(x%deviceVect) if (nd < nh ) then - call trackCudaFree(' c_vect_cuda ',x%sizeof()) + call trackCudaFree(' c_vect_cuda 1',x%sizeof()) call freeMultiVecDevice(x%deviceVect) x%deviceVect=c_null_ptr end if @@ -790,13 +790,13 @@ contains integer(psb_ipk_), intent(out) :: info info = 0 - if (allocated(x%v)) deallocate(x%v, stat=info) if (c_associated(x%deviceVect)) then !!$ write(0,*)'d_cuda_free Calling freeMultiVecDevice' - call trackCudaFree(' c_vect_cuda ',x%sizeof()) + call trackCudaFree(' c_vect_cuda 2',x%sizeof()) call freeMultiVecDevice(x%deviceVect) x%deviceVect=c_null_ptr end if + if (allocated(x%v)) deallocate(x%v, stat=info) call x%free_buffer(info) call x%set_sync() end subroutine c_cuda_free diff --git a/cuda/psb_d_cuda_vect_mod.F90 b/cuda/psb_d_cuda_vect_mod.F90 index b1a8d3d4..311527db 100644 --- a/cuda/psb_d_cuda_vect_mod.F90 +++ b/cuda/psb_d_cuda_vect_mod.F90 @@ -729,7 +729,7 @@ contains if (c_associated(x%deviceVect)) then nd = getMultiVecDeviceSize(x%deviceVect) if (nd < nh ) then - call trackCudaFree(' d_vect_cuda ',x%sizeof()) + call trackCudaFree(' d_vect_cuda 1',x%sizeof()) call freeMultiVecDevice(x%deviceVect) x%deviceVect=c_null_ptr end if @@ -790,13 +790,13 @@ contains integer(psb_ipk_), intent(out) :: info info = 0 - if (allocated(x%v)) deallocate(x%v, stat=info) if (c_associated(x%deviceVect)) then !!$ write(0,*)'d_cuda_free Calling freeMultiVecDevice' - call trackCudaFree(' d_vect_cuda ',x%sizeof()) + call trackCudaFree(' d_vect_cuda 2',x%sizeof()) call freeMultiVecDevice(x%deviceVect) x%deviceVect=c_null_ptr end if + if (allocated(x%v)) deallocate(x%v, stat=info) call x%free_buffer(info) call x%set_sync() end subroutine d_cuda_free diff --git a/cuda/psb_i_cuda_vect_mod.F90 b/cuda/psb_i_cuda_vect_mod.F90 index 91fb8a91..573f14be 100644 --- a/cuda/psb_i_cuda_vect_mod.F90 +++ b/cuda/psb_i_cuda_vect_mod.F90 @@ -711,7 +711,7 @@ contains if (c_associated(x%deviceVect)) then nd = getMultiVecDeviceSize(x%deviceVect) if (nd < nh ) then - call trackCudaFree(' i_vect_cuda ',x%sizeof()) + call trackCudaFree(' i_vect_cuda 1',x%sizeof()) call freeMultiVecDevice(x%deviceVect) x%deviceVect=c_null_ptr end if @@ -772,13 +772,13 @@ contains integer(psb_ipk_), intent(out) :: info info = 0 - if (allocated(x%v)) deallocate(x%v, stat=info) if (c_associated(x%deviceVect)) then !!$ write(0,*)'d_cuda_free Calling freeMultiVecDevice' - call trackCudaFree(' i_vect_cuda ',x%sizeof()) + call trackCudaFree(' i_vect_cuda 2',x%sizeof()) call freeMultiVecDevice(x%deviceVect) x%deviceVect=c_null_ptr end if + if (allocated(x%v)) deallocate(x%v, stat=info) call x%free_buffer(info) call x%set_sync() end subroutine i_cuda_free diff --git a/cuda/psb_s_cuda_vect_mod.F90 b/cuda/psb_s_cuda_vect_mod.F90 index 9aec58f9..e0e229d7 100644 --- a/cuda/psb_s_cuda_vect_mod.F90 +++ b/cuda/psb_s_cuda_vect_mod.F90 @@ -729,7 +729,7 @@ contains if (c_associated(x%deviceVect)) then nd = getMultiVecDeviceSize(x%deviceVect) if (nd < nh ) then - call trackCudaFree(' s_vect_cuda ',x%sizeof()) + call trackCudaFree(' s_vect_cuda 1',x%sizeof()) call freeMultiVecDevice(x%deviceVect) x%deviceVect=c_null_ptr end if @@ -790,13 +790,13 @@ contains integer(psb_ipk_), intent(out) :: info info = 0 - if (allocated(x%v)) deallocate(x%v, stat=info) if (c_associated(x%deviceVect)) then !!$ write(0,*)'d_cuda_free Calling freeMultiVecDevice' - call trackCudaFree(' s_vect_cuda ',x%sizeof()) + call trackCudaFree(' s_vect_cuda 2',x%sizeof()) call freeMultiVecDevice(x%deviceVect) x%deviceVect=c_null_ptr end if + if (allocated(x%v)) deallocate(x%v, stat=info) call x%free_buffer(info) call x%set_sync() end subroutine s_cuda_free diff --git a/cuda/psb_z_cuda_vect_mod.F90 b/cuda/psb_z_cuda_vect_mod.F90 index 579fdd41..b8a53de1 100644 --- a/cuda/psb_z_cuda_vect_mod.F90 +++ b/cuda/psb_z_cuda_vect_mod.F90 @@ -729,7 +729,7 @@ contains if (c_associated(x%deviceVect)) then nd = getMultiVecDeviceSize(x%deviceVect) if (nd < nh ) then - call trackCudaFree(' z_vect_cuda ',x%sizeof()) + call trackCudaFree(' z_vect_cuda 1',x%sizeof()) call freeMultiVecDevice(x%deviceVect) x%deviceVect=c_null_ptr end if @@ -790,13 +790,13 @@ contains integer(psb_ipk_), intent(out) :: info info = 0 - if (allocated(x%v)) deallocate(x%v, stat=info) if (c_associated(x%deviceVect)) then !!$ write(0,*)'d_cuda_free Calling freeMultiVecDevice' - call trackCudaFree(' z_vect_cuda ',x%sizeof()) + call trackCudaFree(' z_vect_cuda 2',x%sizeof()) call freeMultiVecDevice(x%deviceVect) x%deviceVect=c_null_ptr end if + if (allocated(x%v)) deallocate(x%v, stat=info) call x%free_buffer(info) call x%set_sync() end subroutine z_cuda_free