! Parallel Sparse BLAS GPU plugin ! (C) Copyright 2013 ! ! Salvatore Filippone ! Alessandro Fanfarillo ! ! Redistribution and use in source and binary forms, with or without ! modification, are permitted provided that the following conditions ! are met: ! 1. Redistributions of source code must retain the above copyright ! notice, this list of conditions and the following disclaimer. ! 2. Redistributions in binary form must reproduce the above copyright ! notice, this list of conditions, and the following disclaimer in the ! documentation and/or other materials provided with the distribution. ! 3. The name of the PSBLAS group or the names of its contributors may ! not be used to endorse or promote products derived from this ! software without specific written permission. ! ! THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS ! ``AS IS'' AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED ! TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR ! PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE PSBLAS GROUP OR ITS CONTRIBUTORS ! BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR ! CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF ! SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS ! INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN ! CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ! ARISING IN ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE ! POSSIBILITY OF SUCH DAMAGE. ! subroutine psb_s_cuda_diag_to_gpu(a,info,nzrm) use psb_base_mod use diagdev_mod use psb_vectordev_mod use psb_s_cuda_diag_mat_mod, psb_protect_name => psb_s_cuda_diag_to_gpu use iso_c_binding implicit none class(psb_s_cuda_diag_sparse_mat), intent(inout) :: a integer(psb_ipk_), intent(out) :: info integer(psb_ipk_), intent(in), optional :: nzrm integer(psb_ipk_) :: m, nzm, n, c,pitch,maxrowsize,d type(diagdev_parms) :: gpu_parms info = 0 if ((.not.allocated(a%data)).or.(.not.allocated(a%offset))) return n = size(a%data,1) d = size(a%data,2) c = a%get_ncols() !allocsize = a%get_size() !write(*,*) 'Create the DIAG matrix' gpu_parms = FgetDiagDeviceParams(n,c,d,spgpu_type_float) if (c_associated(a%deviceMat)) then call freeDiagDevice(a%deviceMat) endif info = FallocDiagDevice(a%deviceMat,n,c,d,spgpu_type_float) if (info == 0) info = & & writeDiagDevice(a%deviceMat,a%data,a%offset,n) ! if (info /= 0) goto 9999 end subroutine psb_s_cuda_diag_to_gpu