! 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_c_cuda_elg_allocate_mnnz(m,n,a,nz) use psb_base_mod #ifdef HAVE_SPGPU use elldev_mod use psb_vectordev_mod use psb_c_cuda_elg_mat_mod, psb_protect_name => psb_c_cuda_elg_allocate_mnnz #else use psb_c_cuda_elg_mat_mod #endif implicit none integer(psb_ipk_), intent(in) :: m,n class(psb_c_cuda_elg_sparse_mat), intent(inout) :: a integer(psb_ipk_), intent(in), optional :: nz Integer(Psb_ipk_) :: err_act, info, nz_,ld character(len=20) :: name='allocate_mnz' logical, parameter :: debug=.false. #ifdef HAVE_SPGPU type(elldev_parms) :: gpu_parms #endif call psb_erractionsave(err_act) info = psb_success_ if (m < 0) then info = psb_err_iarg_neg_ call psb_errpush(info,name,i_err=(/ione,izero,izero,izero,izero/)) goto 9999 endif if (n < 0) then info = psb_err_iarg_neg_ call psb_errpush(info,name,i_err=(/2*ione,izero,izero,izero,izero/)) goto 9999 endif if (present(nz)) then nz_ = (max(nz,ione) + m -1 )/m else nz_ = (max(7*m,7*n,ione)+m-1)/m end if if (nz_ < 0) then info = psb_err_iarg_neg_ call psb_errpush(info,name,i_err=(/3*ione,izero,izero,izero,izero/)) goto 9999 endif #ifdef HAVE_SPGPU gpu_parms = FgetEllDeviceParams(m,nz_,nz_*m,n,spgpu_type_complex_float,1) ld = gpu_parms%pitch nz_ = gpu_parms%maxRowSize #else ld = m #endif if (info == psb_success_) call psb_realloc(m,a%irn,info) if (info == psb_success_) call psb_realloc(m,a%idiag,info) if (info == psb_success_) call psb_realloc(ld,nz_,a%ja,info) if (info == psb_success_) call psb_realloc(ld,nz_,a%val,info) if (info == psb_success_) then a%irn = 0 a%idiag = 0 a%nzt = 0 call a%set_nrows(m) call a%set_ncols(n) call a%set_bld() call a%set_triangle(.false.) call a%set_unit(.false.) call a%set_dupl(psb_dupl_def_) end if #ifdef HAVE_SPGPU call a%to_gpu(info,nzrm=nz_) if (info /= 0) goto 9999 #endif call psb_erractionrestore(err_act) return 9999 call psb_error_handler(err_act) return end subroutine psb_c_cuda_elg_allocate_mnnz