mod_cuda.F90 29.1 KB
Newer Older
1 2 3 4 5 6 7 8 9 10 11
!    This file is part of ELPA.
!
!    The ELPA library was originally created by the ELPA consortium,
!    consisting of the following organizations:
!
!    - Rechenzentrum Garching der Max-Planck-Gesellschaft (RZG),
!    - Bergische Universität Wuppertal, Lehrstuhl für angewandte
!      Informatik,
!    - Technische Universität München, Lehrstuhl für Informatik mit
!      Schwerpunkt Wissenschaftliches Rechnen ,
!    - Fritz-Haber-Institut, Berlin, Abt. Theorie,
12
!    - Max-Plack-Institut für Mathematik in den Naturwissenschaften,
13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43
!      Leipzig, Abt. Komplexe Strukutren in Biologie und Kognition,
!      and
!    - IBM Deutschland GmbH
!
!
!    More information can be found here:
!    http://elpa.rzg.mpg.de/
!
!    ELPA is free software: you can redistribute it and/or modify
!    it under the terms of the version 3 of the license of the
!    GNU Lesser General Public License as published by the Free
!    Software Foundation.
!
!    ELPA is distributed in the hope that it will be useful,
!    but WITHOUT ANY WARRANTY; without even the implied warranty of
!    MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
!    GNU Lesser General Public License for more details.
!
!    You should have received a copy of the GNU Lesser General Public License
!    along with ELPA.  If not, see <http://www.gnu.org/licenses/>
!
!    ELPA reflects a substantial effort on the part of the original
!    ELPA consortium, and we ask you to respect the spirit of the
!    license that we chose: i.e., please contribute any changes you
!    may have back to the original ELPA library distribution, and keep
!    any derivatives of ELPA under the same license that we chose for
!    the original distribution, the GNU Lesser General Public License.
!
! This file was written by A. Marek, MPCDF


44 45 46
#include "config-f90.h"
module cuda_functions
  use iso_c_binding
47
  use precision
48 49 50 51
  implicit none

  public

52 53 54 55 56
  integer(kind=ik) :: cudaMemcpyHostToDevice
  integer(kind=ik) :: cudaMemcpyDeviceToHost
  integer(kind=ik) :: cudaHostRegisterPortable
  integer(kind=ik) :: cudaHostRegisterMapped
  integer(kind=ik) :: cudaMemcpyDeviceToDevice
Andreas Marek's avatar
Andreas Marek committed
57

58 59
  ! TODO global variable, has to be changed
  integer(kind=C_intptr_T) :: cublasHandle = -1
60

61
  integer(kind=c_intptr_t), parameter :: size_of_double_real    = 8_rk8
62
#ifdef WANT_SINGLE_PRECISION_REAL
63
  integer(kind=c_intptr_t), parameter :: size_of_single_real    = 4_rk4
64 65
#endif

66
  integer(kind=c_intptr_t), parameter :: size_of_double_complex = 16_ck8
67
#ifdef WANT_SINGLE_PRECISION_COMPLEX
68
  integer(kind=c_intptr_t), parameter :: size_of_single_complex = 8_ck4
69
#endif
70 71

  ! functions to set and query the CUDA devices
Pavel Kus's avatar
Pavel Kus committed
72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90
  interface
    function cublas_create_c(handle) result(istat) &
             bind(C, name="cublasCreateFromC")
      use iso_c_binding
      implicit none
      integer(kind=C_intptr_T) :: handle
      integer(kind=C_INT)  :: istat
    end function cublas_create_c
  end interface  

  interface
    function cublas_destroy_c(handle) result(istat) &
             bind(C, name="cublasDestroyFromC")
      use iso_c_binding
      implicit none
      integer(kind=C_intptr_T) :: handle
      integer(kind=C_INT)  :: istat
    end function cublas_destroy_c
  end interface  
91

92 93 94 95 96 97 98 99 100
  interface
    function cuda_threadsynchronize_c() result(istat) &
             bind(C,name="cudaThreadSynchronizeFromC")
      use iso_c_binding
      implicit none
      integer(kind=C_INT)  :: istat
    end function cuda_threadsynchronize_c
  end interface

101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189
  interface
    function cuda_setdevice_c(n) result(istat) &
             bind(C, name="cudaSetDeviceFromC")

      use iso_c_binding
      implicit none
      integer(kind=C_INT), value    :: n
      integer(kind=C_INT)           :: istat
    end function cuda_setdevice_c
  end interface

  interface
    function cuda_getdevicecount_c(n) result(istat) &
             bind(C, name="cudaGetDeviceCountFromC")
      use iso_c_binding
      implicit none
      integer(kind=C_INT), intent(out) :: n
      integer(kind=C_INT)              :: istat
    end function cuda_getdevicecount_c
  end interface

  interface
    function cuda_devicesynchronize_c()result(istat) &
             bind(C,name='cudaDeviceSynchronizeFromC')

      use iso_c_binding

      implicit none
      integer(kind=C_INT)                       :: istat

    end function cuda_devicesynchronize_c
  end interface


  ! functions to copy CUDA memory
  interface
    function cuda_memcpyDeviceToDevice_c() result(flag) &
             bind(C, name="cudaMemcpyDeviceToDeviceFromC")
      use iso_c_binding
      implicit none
      integer(kind=c_int) :: flag
    end function
  end interface

  interface
    function cuda_memcpyHostToDevice_c() result(flag) &
             bind(C, name="cudaMemcpyHostToDeviceFromC")
      use iso_c_binding
      implicit none
      integer(kind=c_int) :: flag
    end function
  end interface

  interface
    function cuda_memcpyDeviceToHost_c() result(flag) &
             bind(C, name="cudaMemcpyDeviceToHostFromC")
      use iso_c_binding
      implicit none
      integer(kind=c_int) :: flag
    end function
  end interface

  interface
    function cuda_hostRegisterPortable_c() result(flag) &
             bind(C, name="cudaHostRegisterPortableFromC")
      use iso_c_binding
      implicit none
      integer(kind=c_int) :: flag
    end function
  end interface

  interface
    function cuda_hostRegisterMapped_c() result(flag) &
             bind(C, name="cudaHostRegisterMappedFromC")
      use iso_c_binding
      implicit none
      integer(kind=c_int) :: flag
    end function
  end interface

  interface
    function cuda_memcpy_c(dst, src, size, dir) result(istat) &
             bind(C, name="cudaMemcpyFromC")

      use iso_c_binding

      implicit none
      integer(kind=C_intptr_t), value              :: dst
      integer(kind=C_intptr_t), value              :: src
190
      integer(kind=c_intptr_t), intent(in), value    :: size
191 192 193 194 195 196 197 198 199 200 201 202 203 204
      integer(kind=C_INT), intent(in), value       :: dir
      integer(kind=C_INT)                          :: istat

    end function cuda_memcpy_c
  end interface

  interface
    function cuda_memcpy2d_c(dst, dpitch, src, spitch, width, height , dir) result(istat) &
             bind(C, name="cudaMemcpy2dFromC")

      use iso_c_binding

      implicit none

205
      integer(kind=C_intptr_T), value                :: dst
206
      integer(kind=c_intptr_t), intent(in), value    :: dpitch
207
      integer(kind=C_intptr_T), value                :: src
208 209 210
      integer(kind=c_intptr_t), intent(in), value    :: spitch
      integer(kind=c_intptr_t), intent(in), value    :: width
      integer(kind=c_intptr_t), intent(in), value    :: height
211 212
      integer(kind=C_INT), intent(in), value         :: dir
      integer(kind=C_INT)                            :: istat
213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239

    end function cuda_memcpy2d_c
  end interface

  ! functions to allocate and free CUDA memory

  interface
    function cuda_free_c(a) result(istat) &
             bind(C, name="cudaFreeFromC")

      use iso_c_binding

      implicit none
      integer(kind=C_intptr_T), value  :: a
      integer(kind=C_INT)              :: istat

    end function cuda_free_c
  end interface

  interface
    function cuda_malloc_c(a, width_height) result(istat) &
             bind(C, name="cudaMallocFromC")

      use iso_c_binding
      implicit none

      integer(kind=C_intptr_T)                    :: a
240
      integer(kind=c_intptr_t), intent(in), value   :: width_height
241 242 243 244 245 246 247 248 249 250 251 252 253 254 255
      integer(kind=C_INT)                         :: istat

    end function cuda_malloc_c
  end interface

  interface
    function cuda_memset_c(a, val, size) result(istat) &
             bind(C, name="cudaMemsetFromC")

      use iso_c_binding

      implicit none

      integer(kind=C_intptr_T), value            :: a
      integer(kind=C_INT), value                 :: val
256
      integer(kind=c_intptr_t), intent(in), value  :: size
257 258 259 260 261 262 263
      integer(kind=C_INT)                        :: istat

    end function cuda_memset_c
  end interface

  ! cuBLAS
  interface
Pavel Kus's avatar
Pavel Kus committed
264 265
    subroutine cublas_dgemm_c(handle, cta, ctb, m, n, k, alpha, a, lda, b, ldb, beta, c, ldc) &
                              bind(C,name='cublasDgemm_elpa_wrapper')
266 267 268 269 270 271 272 273 274

      use iso_c_binding

      implicit none
      character(1,C_CHAR),value               :: cta, ctb
      integer(kind=C_INT),value               :: m,n,k
      integer(kind=C_INT), intent(in), value  :: lda,ldb,ldc
      real(kind=C_DOUBLE),value               :: alpha,beta
      integer(kind=C_intptr_T), value         :: a, b, c
Pavel Kus's avatar
Pavel Kus committed
275 276
      integer(kind=C_intptr_T), value         :: handle

277 278 279
    end subroutine cublas_dgemm_c
  end interface

280
  interface
Pavel Kus's avatar
Pavel Kus committed
281 282
    subroutine cublas_sgemm_c(handle, cta, ctb, m, n, k, alpha, a, lda, b, ldb, beta, c, ldc) &
                              bind(C,name='cublasSgemm_elpa_wrapper')
283 284 285 286 287 288 289 290 291

      use iso_c_binding

      implicit none
      character(1,C_CHAR),value               :: cta, ctb
      integer(kind=C_INT),value               :: m,n,k
      integer(kind=C_INT), intent(in), value  :: lda,ldb,ldc
      real(kind=C_FLOAT),value                :: alpha,beta
      integer(kind=C_intptr_T), value         :: a, b, c
Pavel Kus's avatar
Pavel Kus committed
292 293
      integer(kind=C_intptr_T), value         :: handle
      
294 295 296
    end subroutine cublas_sgemm_c
  end interface

297
  interface
Pavel Kus's avatar
Pavel Kus committed
298 299
    subroutine cublas_dtrmm_c(handle, side, uplo, trans, diag, m, n, alpha, a, lda, b, ldb) &
                              bind(C,name='cublasDtrmm_elpa_wrapper')
300 301 302 303 304 305 306 307 308

      use iso_c_binding

      implicit none
      character(1,C_CHAR),value               :: side, uplo, trans, diag
      integer(kind=C_INT),value               :: m,n
      integer(kind=C_INT), intent(in), value  :: lda,ldb
      real(kind=C_DOUBLE), value              :: alpha
      integer(kind=C_intptr_T), value         :: a, b
Pavel Kus's avatar
Pavel Kus committed
309 310
      integer(kind=C_intptr_T), value         :: handle
      
311 312 313
    end subroutine cublas_dtrmm_c
  end interface

314
  interface
Pavel Kus's avatar
Pavel Kus committed
315 316
    subroutine cublas_strmm_c(handle, side, uplo, trans, diag, m, n, alpha, a, lda, b, ldb) &
                              bind(C,name='cublasStrmm_elpa_wrapper')
317 318 319 320 321 322 323 324 325

      use iso_c_binding

      implicit none
      character(1,C_CHAR),value               :: side, uplo, trans, diag
      integer(kind=C_INT),value               :: m,n
      integer(kind=C_INT), intent(in), value  :: lda,ldb
      real(kind=C_FLOAT), value               :: alpha
      integer(kind=C_intptr_T), value         :: a, b
Pavel Kus's avatar
Pavel Kus committed
326 327
      integer(kind=C_intptr_T), value         :: handle

328 329 330
    end subroutine cublas_strmm_c
  end interface

331
  interface
Pavel Kus's avatar
Pavel Kus committed
332 333
    subroutine cublas_zgemm_c(handle, cta, ctb, m, n, k, alpha, a, lda, b, ldb, beta, c,ldc) &
                              bind(C,name='cublasZgemm_elpa_wrapper')
334 335 336 337 338 339 340

      use iso_c_binding

      implicit none
      character(1,C_CHAR),value              :: cta, ctb
      integer(kind=C_INT),value              :: m,n,k
      integer(kind=C_INT), intent(in), value :: lda,ldb,ldc
341
      complex(kind=C_DOUBLE_COMPLEX),value           :: alpha,beta
342
      integer(kind=C_intptr_T), value        :: a, b, c
Pavel Kus's avatar
Pavel Kus committed
343
      integer(kind=C_intptr_T), value         :: handle
344 345 346 347

    end subroutine cublas_zgemm_c
  end interface

348
  interface
Pavel Kus's avatar
Pavel Kus committed
349 350
    subroutine cublas_cgemm_c(handle, cta, ctb, m, n, k, alpha, a, lda, b, ldb, beta, c,ldc) &
                              bind(C,name='cublasCgemm_elpa_wrapper')
351 352 353 354 355 356 357

      use iso_c_binding

      implicit none
      character(1,C_CHAR),value              :: cta, ctb
      integer(kind=C_INT),value              :: m,n,k
      integer(kind=C_INT), intent(in), value :: lda,ldb,ldc
358
      complex(kind=C_FLOAT_COMPLEX),value            :: alpha,beta
359
      integer(kind=C_intptr_T), value        :: a, b, c
Pavel Kus's avatar
Pavel Kus committed
360
      integer(kind=C_intptr_T), value         :: handle
361 362 363 364

    end subroutine cublas_cgemm_c
  end interface

365
  interface
Pavel Kus's avatar
Pavel Kus committed
366 367
    subroutine cublas_ztrmm_c(handle, side, uplo, trans, diag, m, n, alpha, a, lda, b, ldb) &
                              bind(C,name='cublasZtrmm_elpa_wrapper')
368 369 370 371 372 373 374

      use iso_c_binding

      implicit none
      character(1,C_CHAR),value              :: side, uplo, trans, diag
      integer(kind=C_INT),value              :: m,n
      integer(kind=C_INT), intent(in), value :: lda,ldb
375
      complex(kind=C_DOUBLE_COMPLEX), value          :: alpha
376
      integer(kind=C_intptr_T), value        :: a, b
Pavel Kus's avatar
Pavel Kus committed
377
      integer(kind=C_intptr_T), value         :: handle
378 379 380 381

    end subroutine cublas_ztrmm_c
  end interface

382
  interface
Pavel Kus's avatar
Pavel Kus committed
383 384
    subroutine cublas_ctrmm_c(handle, side, uplo, trans, diag, m, n, alpha, a, lda, b, ldb) &
                              bind(C,name='cublasCtrmm_elpa_wrapper')
385 386 387 388 389 390 391

      use iso_c_binding

      implicit none
      character(1,C_CHAR),value              :: side, uplo, trans, diag
      integer(kind=C_INT),value              :: m,n
      integer(kind=C_INT), intent(in), value :: lda,ldb
392
      complex(kind=C_FLOAT_COMPLEX), value           :: alpha
393
      integer(kind=C_intptr_T), value        :: a, b
Pavel Kus's avatar
Pavel Kus committed
394
      integer(kind=C_intptr_T), value         :: handle
395 396 397 398

    end subroutine cublas_ctrmm_c
  end interface

399
  interface
Pavel Kus's avatar
Pavel Kus committed
400 401
    subroutine cublas_dgemv_c(handle, cta, m, n, alpha, a, lda, x, incx, beta, y, incy) &
                              bind(C,name='cublasDgemv_elpa_wrapper')
402 403 404 405 406 407 408 409 410

      use iso_c_binding

      implicit none
      character(1,C_CHAR),value               :: cta
      integer(kind=C_INT),value               :: m,n
      integer(kind=C_INT), intent(in), value  :: lda,incx,incy
      real(kind=C_DOUBLE),value               :: alpha,beta
      integer(kind=C_intptr_T), value         :: a, x, y
Pavel Kus's avatar
Pavel Kus committed
411 412
      integer(kind=C_intptr_T), value         :: handle

413 414 415 416
    end subroutine cublas_dgemv_c
  end interface

  interface
Pavel Kus's avatar
Pavel Kus committed
417 418
    subroutine cublas_sgemv_c(handle, cta, m, n, alpha, a, lda, x, incx, beta, y, incy) &
                              bind(C,name='cublasSgemv_elpa_wrapper')
419 420 421 422 423 424 425 426 427

      use iso_c_binding

      implicit none
      character(1,C_CHAR),value               :: cta
      integer(kind=C_INT),value               :: m,n
      integer(kind=C_INT), intent(in), value  :: lda,incx,incy
      real(kind=C_FLOAT),value                :: alpha,beta
      integer(kind=C_intptr_T), value         :: a, x, y
Pavel Kus's avatar
Pavel Kus committed
428 429
      integer(kind=C_intptr_T), value         :: handle

430 431 432
    end subroutine cublas_sgemv_c
  end interface

433
  interface
Pavel Kus's avatar
Pavel Kus committed
434 435
    subroutine cublas_zgemv_c(handle, cta, m, n, alpha, a, lda, x, incx, beta, y, incy) &
                              bind(C,name='cublasZgemv_elpa_wrapper')
436 437 438 439 440 441 442

      use iso_c_binding

      implicit none
      character(1,C_CHAR),value               :: cta
      integer(kind=C_INT),value               :: m,n
      integer(kind=C_INT), intent(in), value  :: lda,incx,incy
443
      complex(kind=C_DOUBLE_COMPLEX),value               :: alpha,beta
444
      integer(kind=C_intptr_T), value         :: a, x, y
Pavel Kus's avatar
Pavel Kus committed
445 446
      integer(kind=C_intptr_T), value         :: handle

447 448 449 450
    end subroutine cublas_zgemv_c
  end interface

  interface
Pavel Kus's avatar
Pavel Kus committed
451 452
    subroutine cublas_cgemv_c(handle, cta, m, n, alpha, a, lda, x, incx, beta, y, incy) &
                              bind(C,name='cublasCgemv_elpa_wrapper')
453 454 455 456 457 458 459

      use iso_c_binding

      implicit none
      character(1,C_CHAR),value               :: cta
      integer(kind=C_INT),value               :: m,n
      integer(kind=C_INT), intent(in), value  :: lda,incx,incy
460
      complex(kind=C_FLOAT_COMPLEX),value                :: alpha,beta
461
      integer(kind=C_intptr_T), value         :: a, x, y
Pavel Kus's avatar
Pavel Kus committed
462
      integer(kind=C_intptr_T), value         :: handle
Pavel Kus's avatar
Pavel Kus committed
463

Pavel Kus's avatar
Pavel Kus committed
464
    end subroutine cublas_cgemv_c
Pavel Kus's avatar
Pavel Kus committed
465 466 467
  end interface


Pavel Kus's avatar
Pavel Kus committed
468
  contains
Pavel Kus's avatar
Pavel Kus committed
469

Pavel Kus's avatar
Pavel Kus committed
470
    ! functions to set and query the CUDA devices
Pavel Kus's avatar
Pavel Kus committed
471

Pavel Kus's avatar
Pavel Kus committed
472 473 474
   function cublas_create(handle) result(success)
     use iso_c_binding
     implicit none
475

Pavel Kus's avatar
Pavel Kus committed
476 477 478 479 480 481 482 483
     integer(kind=C_intptr_t)                  :: handle
     logical                                   :: success
#ifdef WITH_GPU_VERSION
     success = cublas_create_c(handle) /= 0
#else
     success = .true.
#endif
   end function
484

Pavel Kus's avatar
Pavel Kus committed
485 486 487
   function cublas_destroy(handle) result(success)
     use iso_c_binding
     implicit none
488

Pavel Kus's avatar
Pavel Kus committed
489 490 491 492 493 494 495 496 497
     integer(kind=C_intptr_t)                  :: handle
     logical                                   :: success
#ifdef WITH_GPU_VERSION
     success = cublas_destroy_c(handle) /= 0
#else
     success = .true.
#endif
   end function
    
498 499 500 501 502 503 504 505 506 507 508 509 510
    function cuda_threadsynchronize() result(success)
      use iso_c_binding

      implicit none

      logical :: success
#ifdef WITH_GPU_VERSION
      success = cuda_threadsynchronize_c() /= 0
#else
      success = .true.
#endif
    end function cuda_threadsynchronize

511 512 513 514 515
    function cuda_setdevice(n) result(success)
      use iso_c_binding

      implicit none

516 517
      integer(kind=ik), intent(in)  :: n
      logical                       :: success
518 519 520 521 522 523 524 525 526 527 528
#ifdef WITH_GPU_VERSION
      success = cuda_setdevice_c(int(n,kind=c_int)) /= 0
#else
      success = .true.
#endif
    end function cuda_setdevice

    function cuda_getdevicecount(n) result(success)
      use iso_c_binding
      implicit none

529
      integer(kind=ik)     :: n
530 531 532 533 534 535 536 537 538 539 540 541 542 543 544 545 546 547 548 549
      integer(kind=c_int)  :: nCasted
      logical              :: success
#ifdef WITH_GPU_VERSION
      success = cuda_getdevicecount_c(nCasted) /=0
      n = int(nCasted)
#else
      success = .true.
      n = 0
#endif
    end function cuda_getdevicecount

    function cuda_devicesynchronize()result(success)

      use iso_c_binding

      implicit none
      logical :: success
#ifdef WITH_GPU_VERSION
      success = cuda_devicesynchronize_c() /=0
#else
550
      success = .true.
551 552 553 554 555 556 557 558 559 560
#endif
    end function cuda_devicesynchronize
    ! functions to allocate and free memory

    function cuda_malloc(a, width_height) result(success)

     use iso_c_binding
     implicit none

     integer(kind=C_intptr_t)                  :: a
561
     integer(kind=c_intptr_t), intent(in)        :: width_height
562 563 564 565
     logical                                   :: success
#ifdef WITH_GPU_VERSION
     success = cuda_malloc_c(a, width_height) /= 0
#else
566
     success = .true.
567 568 569 570 571 572 573 574 575 576 577 578 579
#endif
   end function

   function cuda_free(a) result(success)

     use iso_c_binding

     implicit none
     integer(kind=C_intptr_T) :: a
     logical                  :: success
#ifdef WITH_GPU_VERSION
     success = cuda_free_c(a) /= 0
#else
580
     success = .true.
581 582 583 584 585 586 587 588 589 590
#endif
   end function cuda_free

 function cuda_memset(a, val, size) result(success)

   use iso_c_binding

   implicit none

   integer(kind=c_intptr_t)                :: a
591
   integer(kind=ik)                        :: val
592
   integer(kind=c_intptr_t), intent(in)      :: size
593 594 595 596
   integer(kind=C_INT)                     :: istat

   logical :: success
#ifdef WITH_GPU_VERSION
597
   success= cuda_memset_c(a, int(val,kind=c_int), int(size,kind=c_intptr_t)) /=0
598
#else
599
   success = .true.
600 601 602 603 604 605 606 607
#endif
 end function cuda_memset

 ! functions to memcopy CUDA memory

 function cuda_memcpyDeviceToDevice() result(flag)
   use iso_c_binding
   implicit none
608
   integer(kind=ik) :: flag
609 610 611 612 613 614 615 616 617
#ifdef WITH_GPU_VERSION
   flag = int(cuda_memcpyDeviceToDevice_c())
#else
   flag = 0
#endif
 end function

 function cuda_memcpyHostToDevice() result(flag)
   use iso_c_binding
618
   use precision
619
   implicit none
620
   integer(kind=ik) :: flag
621 622 623 624 625 626 627 628 629
#ifdef WITH_GPU_VERSION
   flag = int(cuda_memcpyHostToDevice_c())
#else
   flag = 0
#endif
 end function

 function cuda_memcpyDeviceToHost() result(flag)
   use iso_c_binding
630
   use precision
631
   implicit none
632
   integer(kind=ik) :: flag
633 634 635 636 637 638 639 640 641
#ifdef WITH_GPU_VERSION
   flag = int( cuda_memcpyDeviceToHost_c())
#else
   flag = 0
#endif
 end function

 function cuda_hostRegisterPortable() result(flag)
   use iso_c_binding
642
   use precision
643
   implicit none
644
   integer(kind=ik) :: flag
645 646 647 648 649 650 651 652 653
#ifdef WITH_GPU_VERSION
   flag = int(cuda_hostRegisterPortable_c())
#else
   flag = 0
#endif
 end function

 function cuda_hostRegisterMapped() result(flag)
   use iso_c_binding
654
   use precision
655
   implicit none
656
   integer(kind=ik) :: flag
657 658 659 660 661 662 663 664 665 666 667 668 669 670
#ifdef WITH_GPU_VERSION
   flag = int(cuda_hostRegisterMapped_c())
#else
   flag = 0
#endif
 end function

 function cuda_memcpy(dst, src, size, dir) result(success)

      use iso_c_binding

      implicit none
      integer(kind=C_intptr_t)              :: dst
      integer(kind=C_intptr_t)              :: src
671
      integer(kind=c_intptr_t), intent(in)    :: size
672 673
      integer(kind=C_INT), intent(in)       :: dir
      logical :: success
674

675 676 677
#ifdef WITH_GPU_VERSION
        success = cuda_memcpy_c(dst, src, size, dir) /= 0
#else
678
        success = .true.
679 680 681 682 683 684 685 686 687 688
#endif
    end function

    function cuda_memcpy2d(dst, dpitch, src, spitch, width, height , dir) result(success)

      use iso_c_binding

      implicit none

      integer(kind=C_intptr_T)           :: dst
689
      integer(kind=c_intptr_t), intent(in) :: dpitch
690
      integer(kind=C_intptr_T)           :: src
691 692 693
      integer(kind=c_intptr_t), intent(in) :: spitch
      integer(kind=c_intptr_t), intent(in) :: width
      integer(kind=c_intptr_t), intent(in) :: height
694 695 696 697 698
      integer(kind=C_INT), intent(in)    :: dir
      logical                            :: success
#ifdef WITH_GPU_VERSION
      success = cuda_memcpy2d_c(dst, dpitch, src, spitch, width, height , dir) /= 0
#else
699
      success = .true.
700 701 702 703 704 705 706 707
#endif
    end function cuda_memcpy2d

    ! cuBLAS
    subroutine cublas_dgemm(cta, ctb, m, n, k, alpha, a, lda, b, ldb, beta, c, ldc)
      use iso_c_binding

      implicit none
708
      character(1,C_CHAR),value       :: cta, ctb
709 710 711 712 713
      integer(kind=C_INT)             :: m,n,k
      integer(kind=C_INT), intent(in) :: lda,ldb,ldc
      real(kind=C_DOUBLE)             :: alpha,beta
      integer(kind=C_intptr_T)        :: a, b, c
#ifdef WITH_GPU_VERSION
Pavel Kus's avatar
Pavel Kus committed
714
      call cublas_dgemm_c(cublasHandle, cta, ctb, m, n, k, alpha, a, lda, b, ldb, beta, c, ldc)
715 716 717
#endif
    end subroutine cublas_dgemm

718 719 720 721
    subroutine cublas_sgemm(cta, ctb, m, n, k, alpha, a, lda, b, ldb, beta, c, ldc)
      use iso_c_binding

      implicit none
722
      character(1,C_CHAR),value       :: cta, ctb
723 724 725 726 727
      integer(kind=C_INT)             :: m,n,k
      integer(kind=C_INT), intent(in) :: lda,ldb,ldc
      real(kind=C_FLOAT)              :: alpha,beta
      integer(kind=C_intptr_T)        :: a, b, c
#ifdef WITH_GPU_VERSION
Pavel Kus's avatar
Pavel Kus committed
728
      call cublas_sgemm_c(cublasHandle, cta, ctb, m, n, k, alpha, a, lda, b, ldb, beta, c, ldc)
729 730 731
#endif
    end subroutine cublas_sgemm

732 733 734 735 736
    subroutine cublas_dtrmm(side, uplo, trans, diag, m, n, alpha, a, lda, b, ldb)

      use iso_c_binding

      implicit none
737
      character(1,C_CHAR),value       :: side, uplo, trans, diag
738 739 740 741 742
      integer(kind=C_INT)             :: m,n
      integer(kind=C_INT), intent(in) :: lda,ldb
      real(kind=C_DOUBLE)             :: alpha
      integer(kind=C_intptr_T)        :: a, b
#ifdef WITH_GPU_VERSION
Pavel Kus's avatar
Pavel Kus committed
743
      call cublas_dtrmm_c(cublasHandle, side, uplo, trans, diag, m, n, alpha, a, lda, b, ldb)
744 745 746
#endif
    end subroutine cublas_dtrmm

747 748 749 750 751
    subroutine cublas_strmm(side, uplo, trans, diag, m, n, alpha, a, lda, b, ldb)

      use iso_c_binding

      implicit none
752
      character(1,C_CHAR),value       :: side, uplo, trans, diag
753 754 755 756 757
      integer(kind=C_INT)             :: m,n
      integer(kind=C_INT), intent(in) :: lda,ldb
      real(kind=C_FLOAT)              :: alpha
      integer(kind=C_intptr_T)        :: a, b
#ifdef WITH_GPU_VERSION
Pavel Kus's avatar
Pavel Kus committed
758
      call cublas_strmm_c(cublasHandle, side, uplo, trans, diag, m, n, alpha, a, lda, b, ldb)
759 760 761
#endif
    end subroutine cublas_strmm

762 763 764 765 766
    subroutine cublas_zgemm(cta, ctb, m, n, k, alpha, a, lda, b, ldb, beta, c,ldc)

      use iso_c_binding

      implicit none
767
      character(1,C_CHAR),value       :: cta, ctb
768 769
      integer(kind=C_INT)             :: m,n,k
      integer(kind=C_INT), intent(in) :: lda,ldb,ldc
770
      complex(kind=C_DOUBLE_COMPLEX)          :: alpha,beta
771 772
      integer(kind=C_intptr_T)        :: a, b, c
#ifdef WITH_GPU_VERSION
Pavel Kus's avatar
Pavel Kus committed
773
      call cublas_zgemm_c(cublasHandle, cta, ctb, m, n, k, alpha, a, lda, b, ldb, beta, c,ldc)
774 775 776
#endif
    end subroutine cublas_zgemm

777 778 779 780 781
    subroutine cublas_cgemm(cta, ctb, m, n, k, alpha, a, lda, b, ldb, beta, c,ldc)

      use iso_c_binding

      implicit none
782
      character(1,C_CHAR),value       :: cta, ctb
783 784
      integer(kind=C_INT)             :: m,n,k
      integer(kind=C_INT), intent(in) :: lda,ldb,ldc
785
      complex(kind=C_FLOAT_COMPLEX)           :: alpha,beta
786 787
      integer(kind=C_intptr_T)        :: a, b, c
#ifdef WITH_GPU_VERSION
Pavel Kus's avatar
Pavel Kus committed
788
      call cublas_cgemm_c(cublasHandle, cta, ctb, m, n, k, alpha, a, lda, b, ldb, beta, c,ldc)
789 790 791
#endif
    end subroutine cublas_cgemm

792 793 794 795 796
    subroutine cublas_ztrmm(side, uplo, trans, diag, m, n, alpha, a, lda, b, ldb)

      use iso_c_binding

      implicit none
797
      character(1,C_CHAR),value       :: side, uplo, trans, diag
798 799
      integer(kind=C_INT)             :: m,n
      integer(kind=C_INT), intent(in) :: lda,ldb
800
      complex(kind=C_DOUBLE_COMPLEX)          :: alpha
801 802
      integer(kind=C_intptr_T)        :: a, b
#ifdef WITH_GPU_VERSION
Pavel Kus's avatar
Pavel Kus committed
803
      call cublas_ztrmm_c(cublasHandle, side, uplo, trans, diag, m, n, alpha, a, lda, b, ldb)
804 805 806
#endif
    end subroutine cublas_ztrmm

807 808 809 810 811
    subroutine cublas_ctrmm(side, uplo, trans, diag, m, n, alpha, a, lda, b, ldb)

      use iso_c_binding

      implicit none
812
      character(1,C_CHAR),value       :: side, uplo, trans, diag
813 814
      integer(kind=C_INT)             :: m,n
      integer(kind=C_INT), intent(in) :: lda,ldb
815
      complex(kind=C_FLOAT_COMPLEX)           :: alpha
816 817
      integer(kind=C_intptr_T)        :: a, b
#ifdef WITH_GPU_VERSION
Pavel Kus's avatar
Pavel Kus committed
818
      call cublas_ctrmm_c(cublasHandle, side, uplo, trans, diag, m, n, alpha, a, lda, b, ldb)
819 820
#endif
    end subroutine cublas_ctrmm
821

822 823 824 825 826 827 828 829 830 831
    subroutine cublas_dgemv(cta, m, n, alpha, a, lda, x, incx, beta, y, incy)
      use iso_c_binding

      implicit none
      character(1,C_CHAR),value       :: cta
      integer(kind=C_INT)             :: m,n
      integer(kind=C_INT), intent(in) :: lda,incx,incy
      real(kind=C_DOUBLE)             :: alpha,beta
      integer(kind=C_intptr_T)        :: a, x, y
#ifdef WITH_GPU_VERSION
Pavel Kus's avatar
Pavel Kus committed
832
      call cublas_dgemv_c(cublasHandle, cta, m, n, alpha, a, lda, x, incx, beta, y, incy)
833 834 835 836 837 838 839 840 841 842 843 844 845
#endif
    end subroutine cublas_dgemv

    subroutine cublas_sgemv(cta, m, n, alpha, a, lda, x, incx, beta, y, incy)
      use iso_c_binding

      implicit none
      character(1,C_CHAR),value       :: cta
      integer(kind=C_INT)             :: m,n
      integer(kind=C_INT), intent(in) :: lda,incx,incy
      real(kind=C_FLOAT)              :: alpha,beta
      integer(kind=C_intptr_T)        :: a, x, y
#ifdef WITH_GPU_VERSION
Pavel Kus's avatar
Pavel Kus committed
846
      call cublas_sgemv_c(cublasHandle, cta, m, n, alpha, a, lda, x, incx, beta, y, incy)
847 848 849
#endif
    end subroutine cublas_sgemv

850 851 852 853 854 855 856
    subroutine cublas_zgemv(cta, m, n, alpha, a, lda, x, incx, beta, y, incy)
      use iso_c_binding

      implicit none
      character(1,C_CHAR),value       :: cta
      integer(kind=C_INT)             :: m,n
      integer(kind=C_INT), intent(in) :: lda,incx,incy
857
      complex(kind=C_DOUBLE_COMPLEX)             :: alpha,beta
858 859
      integer(kind=C_intptr_T)        :: a, x, y
#ifdef WITH_GPU_VERSION
Pavel Kus's avatar
Pavel Kus committed
860
      call cublas_zgemv_c(cublasHandle, cta, m, n, alpha, a, lda, x, incx, beta, y, incy)
861 862 863 864 865 866 867 868 869 870
#endif
    end subroutine cublas_zgemv

    subroutine cublas_cgemv(cta, m, n, alpha, a, lda, x, incx, beta, y, incy)
      use iso_c_binding

      implicit none
      character(1,C_CHAR),value       :: cta
      integer(kind=C_INT)             :: m,n
      integer(kind=C_INT), intent(in) :: lda,incx,incy
871
      complex(kind=C_FLOAT_COMPLEX)              :: alpha,beta
872 873
      integer(kind=C_intptr_T)        :: a, x, y
#ifdef WITH_GPU_VERSION
Pavel Kus's avatar
Pavel Kus committed
874
      call cublas_cgemv_c(cublasHandle, cta, m, n, alpha, a, lda, x, incx, beta, y, incy)
875 876 877
#endif
    end subroutine cublas_cgemv

878

Pavel Kus's avatar
Pavel Kus committed
879 880 881 882 883 884 885 886 887 888 889 890 891 892 893 894 895 896 897 898 899 900 901 902 903 904 905 906 907 908 909 910 911 912 913 914 915 916 917 918 919 920 921 922 923 924 925 926 927 928 929 930 931 932 933
!     subroutine cublas_dsymv(cta, n, alpha, a, lda, x, incx, beta, y, incy)
!       use iso_c_binding
! 
!       implicit none
!       character(1,C_CHAR),value       :: cta
!       integer(kind=C_INT)             :: n
!       integer(kind=C_INT), intent(in) :: lda,incx,incy
!       real(kind=C_DOUBLE)             :: alpha,beta
!       integer(kind=C_intptr_T)        :: a, x, y
! #ifdef WITH_GPU_VERSION
!       call cublas_dsymv_c(cta, n, alpha, a, lda, x, incx, beta, y, incy)
! #endif
!     end subroutine cublas_dsymv
! 
!     subroutine cublas_ssymv(cta, n, alpha, a, lda, x, incx, beta, y, incy)
!       use iso_c_binding
! 
!       implicit none
!       character(1,C_CHAR),value       :: cta
!       integer(kind=C_INT)             :: n
!       integer(kind=C_INT), intent(in) :: lda,incx,incy
!       real(kind=C_FLOAT)              :: alpha,beta
!       integer(kind=C_intptr_T)        :: a, x, y
! #ifdef WITH_GPU_VERSION
!       call cublas_ssymv_c(cta, n, alpha, a, lda, x, incx, beta, y, incy)
! #endif
!     end subroutine cublas_ssymv
! 
!     subroutine cublas_zsymv(cta, n, alpha, a, lda, x, incx, beta, y, incy)
!       use iso_c_binding
! 
!       implicit none
!       character(1,C_CHAR),value       :: cta
!       integer(kind=C_INT)             :: n
!       integer(kind=C_INT), intent(in) :: lda,incx,incy
!       complex(kind=C_DOUBLE_COMPLEX)             :: alpha,beta
!       integer(kind=C_intptr_T)        :: a, x, y
! #ifdef WITH_GPU_VERSION
! !       call cublas_zsymv_c(cta, n, alpha, a, lda, x, incx, beta, y, incy)
! #endif
!     end subroutine cublas_zsymv
! 
!     subroutine cublas_csymv(cta, n, alpha, a, lda, x, incx, beta, y, incy)
!       use iso_c_binding
! 
!       implicit none
!       character(1,C_CHAR),value       :: cta
!       integer(kind=C_INT)             :: n
!       integer(kind=C_INT), intent(in) :: lda,incx,incy
!       complex(kind=C_FLOAT_COMPLEX)              :: alpha,beta
!       integer(kind=C_intptr_T)        :: a, x, y
! #ifdef WITH_GPU_VERSION
! !       call cublas_csymv_c(cta, n, alpha, a, lda, x, incx, beta, y, incy)
! #endif
!     end subroutine cublas_csymv
Pavel Kus's avatar
Pavel Kus committed
934 935


936
end module cuda_functions