mod_cuda.F90 31.2 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

Pavel Kus's avatar
Pavel Kus committed
331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364
  !TODO so far only double real
  interface
    subroutine cublas_dsyrk_c(handle, uplo, trans, n, k, alpha, a, lda, beta, c, ldc) &
                              bind(C,name='cublasDsyrk_elpa_wrapper')
      use iso_c_binding

      implicit none
      character(1, C_CHAR), value             :: uplo, trans
      integer(kind=C_INT), value              :: n, k
      integer(kind=C_INT), intent(in), value  :: lda, ldc
      real(kind=C_DOUBLE), value              :: alpha, beta
      integer(kind=C_intptr_T), value         :: a, c
      integer(kind=C_intptr_T), value         :: handle

    end subroutine cublas_dsyrk_c
  end interface

  !TODO so far only double real
  interface
    subroutine cublas_dscal_c(handle, n, alpha, x, incx) &
                              bind(C,name='cublasDscal_elpa_wrapper')
      use iso_c_binding

      implicit none
      integer(kind=C_INT), value              :: n
      integer(kind=C_INT), intent(in), value  :: incx
      real(kind=C_DOUBLE), value              :: alpha
      integer(kind=C_intptr_T), value         :: x
      integer(kind=C_intptr_T), value         :: handle

    end subroutine cublas_dscal_c
  end interface


365
  interface
Pavel Kus's avatar
Pavel Kus committed
366 367
    subroutine cublas_zgemm_c(handle, cta, ctb, m, n, k, alpha, a, lda, b, ldb, beta, c,ldc) &
                              bind(C,name='cublasZgemm_elpa_wrapper')
368 369 370 371 372 373 374

      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
375
      complex(kind=C_DOUBLE_COMPLEX),value           :: alpha,beta
376
      integer(kind=C_intptr_T), value        :: a, b, c
Pavel Kus's avatar
Pavel Kus committed
377
      integer(kind=C_intptr_T), value         :: handle
378 379 380 381

    end subroutine cublas_zgemm_c
  end interface

382
  interface
Pavel Kus's avatar
Pavel Kus committed
383 384
    subroutine cublas_cgemm_c(handle, cta, ctb, m, n, k, alpha, a, lda, b, ldb, beta, c,ldc) &
                              bind(C,name='cublasCgemm_elpa_wrapper')
385 386 387 388 389 390 391

      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
392
      complex(kind=C_FLOAT_COMPLEX),value            :: alpha,beta
393
      integer(kind=C_intptr_T), value        :: a, b, c
Pavel Kus's avatar
Pavel Kus committed
394
      integer(kind=C_intptr_T), value         :: handle
395 396 397 398

    end subroutine cublas_cgemm_c
  end interface

399
  interface
Pavel Kus's avatar
Pavel Kus committed
400 401
    subroutine cublas_ztrmm_c(handle, side, uplo, trans, diag, m, n, alpha, a, lda, b, ldb) &
                              bind(C,name='cublasZtrmm_elpa_wrapper')
402 403 404 405 406 407 408

      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
409
      complex(kind=C_DOUBLE_COMPLEX), value          :: alpha
410
      integer(kind=C_intptr_T), value        :: a, b
Pavel Kus's avatar
Pavel Kus committed
411
      integer(kind=C_intptr_T), value         :: handle
412 413 414 415

    end subroutine cublas_ztrmm_c
  end interface

416
  interface
Pavel Kus's avatar
Pavel Kus committed
417 418
    subroutine cublas_ctrmm_c(handle, side, uplo, trans, diag, m, n, alpha, a, lda, b, ldb) &
                              bind(C,name='cublasCtrmm_elpa_wrapper')
419 420 421 422 423 424 425

      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
426
      complex(kind=C_FLOAT_COMPLEX), value           :: alpha
427
      integer(kind=C_intptr_T), value        :: a, b
Pavel Kus's avatar
Pavel Kus committed
428
      integer(kind=C_intptr_T), value         :: handle
429 430 431 432

    end subroutine cublas_ctrmm_c
  end interface

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

      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
445 446
      integer(kind=C_intptr_T), value         :: handle

447 448 449 450
    end subroutine cublas_dgemv_c
  end interface

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

      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
462 463
      integer(kind=C_intptr_T), value         :: handle

464 465 466
    end subroutine cublas_sgemv_c
  end interface

467
  interface
Pavel Kus's avatar
Pavel Kus committed
468 469
    subroutine cublas_zgemv_c(handle, cta, m, n, alpha, a, lda, x, incx, beta, y, incy) &
                              bind(C,name='cublasZgemv_elpa_wrapper')
470 471 472 473 474 475 476

      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
477
      complex(kind=C_DOUBLE_COMPLEX),value               :: alpha,beta
478
      integer(kind=C_intptr_T), value         :: a, x, y
Pavel Kus's avatar
Pavel Kus committed
479 480
      integer(kind=C_intptr_T), value         :: handle

481 482 483 484
    end subroutine cublas_zgemv_c
  end interface

  interface
Pavel Kus's avatar
Pavel Kus committed
485 486
    subroutine cublas_cgemv_c(handle, cta, m, n, alpha, a, lda, x, incx, beta, y, incy) &
                              bind(C,name='cublasCgemv_elpa_wrapper')
487 488 489 490 491 492 493

      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
494
      complex(kind=C_FLOAT_COMPLEX),value                :: alpha,beta
495
      integer(kind=C_intptr_T), value         :: a, x, y
Pavel Kus's avatar
Pavel Kus committed
496
      integer(kind=C_intptr_T), value         :: handle
Pavel Kus's avatar
Pavel Kus committed
497

Pavel Kus's avatar
Pavel Kus committed
498
    end subroutine cublas_cgemv_c
Pavel Kus's avatar
Pavel Kus committed
499 500 501
  end interface


Pavel Kus's avatar
Pavel Kus committed
502
  contains
Pavel Kus's avatar
Pavel Kus committed
503

Pavel Kus's avatar
Pavel Kus committed
504
    ! functions to set and query the CUDA devices
Pavel Kus's avatar
Pavel Kus committed
505

Pavel Kus's avatar
Pavel Kus committed
506 507 508
   function cublas_create(handle) result(success)
     use iso_c_binding
     implicit none
509

Pavel Kus's avatar
Pavel Kus committed
510 511 512 513 514 515 516 517
     integer(kind=C_intptr_t)                  :: handle
     logical                                   :: success
#ifdef WITH_GPU_VERSION
     success = cublas_create_c(handle) /= 0
#else
     success = .true.
#endif
   end function
518

Pavel Kus's avatar
Pavel Kus committed
519 520 521
   function cublas_destroy(handle) result(success)
     use iso_c_binding
     implicit none
522

Pavel Kus's avatar
Pavel Kus committed
523 524 525 526 527 528 529 530 531
     integer(kind=C_intptr_t)                  :: handle
     logical                                   :: success
#ifdef WITH_GPU_VERSION
     success = cublas_destroy_c(handle) /= 0
#else
     success = .true.
#endif
   end function
    
532 533 534 535 536 537 538 539 540 541 542 543 544
    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

545 546 547 548 549
    function cuda_setdevice(n) result(success)
      use iso_c_binding

      implicit none

550 551
      integer(kind=ik), intent(in)  :: n
      logical                       :: success
552 553 554 555 556 557 558 559 560 561 562
#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

563
      integer(kind=ik)     :: n
564 565 566 567 568 569 570 571 572 573 574 575 576 577 578 579 580 581 582 583
      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
584
      success = .true.
585 586 587 588 589 590 591 592 593 594
#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
595
     integer(kind=c_intptr_t), intent(in)        :: width_height
596 597 598 599
     logical                                   :: success
#ifdef WITH_GPU_VERSION
     success = cuda_malloc_c(a, width_height) /= 0
#else
600
     success = .true.
601 602 603 604 605 606 607 608 609 610 611 612 613
#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
614
     success = .true.
615 616 617 618 619 620 621 622 623 624
#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
625
   integer(kind=ik)                        :: val
626
   integer(kind=c_intptr_t), intent(in)      :: size
627 628 629 630
   integer(kind=C_INT)                     :: istat

   logical :: success
#ifdef WITH_GPU_VERSION
631
   success= cuda_memset_c(a, int(val,kind=c_int), int(size,kind=c_intptr_t)) /=0
632
#else
633
   success = .true.
634 635 636 637 638 639 640 641
#endif
 end function cuda_memset

 ! functions to memcopy CUDA memory

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

 function cuda_memcpyHostToDevice() result(flag)
   use iso_c_binding
652
   use precision
653
   implicit none
654
   integer(kind=ik) :: flag
655 656 657 658 659 660 661 662 663
#ifdef WITH_GPU_VERSION
   flag = int(cuda_memcpyHostToDevice_c())
#else
   flag = 0
#endif
 end function

 function cuda_memcpyDeviceToHost() result(flag)
   use iso_c_binding
664
   use precision
665
   implicit none
666
   integer(kind=ik) :: flag
667 668 669 670 671 672 673 674 675
#ifdef WITH_GPU_VERSION
   flag = int( cuda_memcpyDeviceToHost_c())
#else
   flag = 0
#endif
 end function

 function cuda_hostRegisterPortable() result(flag)
   use iso_c_binding
676
   use precision
677
   implicit none
678
   integer(kind=ik) :: flag
679 680 681 682 683 684 685 686 687
#ifdef WITH_GPU_VERSION
   flag = int(cuda_hostRegisterPortable_c())
#else
   flag = 0
#endif
 end function

 function cuda_hostRegisterMapped() result(flag)
   use iso_c_binding
688
   use precision
689
   implicit none
690
   integer(kind=ik) :: flag
691 692 693 694 695 696 697 698 699 700 701 702 703 704
#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
705
      integer(kind=c_intptr_t), intent(in)    :: size
706 707
      integer(kind=C_INT), intent(in)       :: dir
      logical :: success
708

709 710 711
#ifdef WITH_GPU_VERSION
        success = cuda_memcpy_c(dst, src, size, dir) /= 0
#else
712
        success = .true.
713 714 715 716 717 718 719 720 721 722
#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
723
      integer(kind=c_intptr_t), intent(in) :: dpitch
724
      integer(kind=C_intptr_T)           :: src
725 726 727
      integer(kind=c_intptr_t), intent(in) :: spitch
      integer(kind=c_intptr_t), intent(in) :: width
      integer(kind=c_intptr_t), intent(in) :: height
728 729 730 731 732
      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
733
      success = .true.
734 735 736 737 738 739 740 741
#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
742
      character(1,C_CHAR),value       :: cta, ctb
743 744 745 746 747
      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
748
      call cublas_dgemm_c(cublasHandle, cta, ctb, m, n, k, alpha, a, lda, b, ldb, beta, c, ldc)
749 750 751
#endif
    end subroutine cublas_dgemm

752 753 754 755
    subroutine cublas_sgemm(cta, ctb, m, n, k, alpha, a, lda, b, ldb, beta, c, ldc)
      use iso_c_binding

      implicit none
756
      character(1,C_CHAR),value       :: cta, ctb
757 758 759 760 761
      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
762
      call cublas_sgemm_c(cublasHandle, cta, ctb, m, n, k, alpha, a, lda, b, ldb, beta, c, ldc)
763 764 765
#endif
    end subroutine cublas_sgemm

766 767 768 769 770
    subroutine cublas_dtrmm(side, uplo, trans, diag, m, n, alpha, a, lda, b, ldb)

      use iso_c_binding

      implicit none
771
      character(1,C_CHAR),value       :: side, uplo, trans, diag
772 773 774 775 776
      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
777
      call cublas_dtrmm_c(cublasHandle, side, uplo, trans, diag, m, n, alpha, a, lda, b, ldb)
778 779 780
#endif
    end subroutine cublas_dtrmm

781 782 783 784 785
    subroutine cublas_strmm(side, uplo, trans, diag, m, n, alpha, a, lda, b, ldb)

      use iso_c_binding

      implicit none
786
      character(1,C_CHAR),value       :: side, uplo, trans, diag
787 788 789 790 791
      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
792
      call cublas_strmm_c(cublasHandle, side, uplo, trans, diag, m, n, alpha, a, lda, b, ldb)
793 794 795
#endif
    end subroutine cublas_strmm

Pavel Kus's avatar
Pavel Kus committed
796 797 798 799 800 801 802 803 804 805 806 807 808 809 810 811 812 813 814 815 816 817 818 819 820 821 822 823 824 825 826 827
  !TODO so far only double real
    subroutine cublas_dsyrk(uplo, trans, n, k, alpha, a, lda, beta, c, ldc)

      use iso_c_binding

      implicit none
      character(1, C_CHAR), value      :: uplo, trans
      integer(kind=C_INT)              :: n, k
      integer(kind=C_INT), intent(in)  :: lda, ldc
      real(kind=C_DOUBLE)              :: alpha, beta
      integer(kind=C_intptr_T)         :: a, c
#ifdef WITH_GPU_VERSION
      call cublas_dsyrk_c(cublasHandle, uplo, trans, n, k, alpha, a, lda, beta, c, ldc)
#endif
    end subroutine cublas_dsyrk

  !TODO so far only double real
    subroutine cublas_dscal(n, alpha, x, incx)
      use iso_c_binding

      implicit none
      integer(kind=C_INT)              :: n
      integer(kind=C_INT), intent(in)  :: incx
      real(kind=C_DOUBLE)              :: alpha
      integer(kind=C_intptr_T)         :: x
#ifdef WITH_GPU_VERSION
      call cublas_dscal_c(cublasHandle, n, alpha, x, incx)
#endif

    end subroutine cublas_dscal


828 829 830 831 832
    subroutine cublas_zgemm(cta, ctb, m, n, k, alpha, a, lda, b, ldb, beta, c,ldc)

      use iso_c_binding

      implicit none
833
      character(1,C_CHAR),value       :: cta, ctb
834 835
      integer(kind=C_INT)             :: m,n,k
      integer(kind=C_INT), intent(in) :: lda,ldb,ldc
836
      complex(kind=C_DOUBLE_COMPLEX)          :: alpha,beta
837 838
      integer(kind=C_intptr_T)        :: a, b, c
#ifdef WITH_GPU_VERSION
Pavel Kus's avatar
Pavel Kus committed
839
      call cublas_zgemm_c(cublasHandle, cta, ctb, m, n, k, alpha, a, lda, b, ldb, beta, c,ldc)
840 841 842
#endif
    end subroutine cublas_zgemm

843 844 845 846 847
    subroutine cublas_cgemm(cta, ctb, m, n, k, alpha, a, lda, b, ldb, beta, c,ldc)

      use iso_c_binding

      implicit none
848
      character(1,C_CHAR),value       :: cta, ctb
849 850
      integer(kind=C_INT)             :: m,n,k
      integer(kind=C_INT), intent(in) :: lda,ldb,ldc
851
      complex(kind=C_FLOAT_COMPLEX)           :: alpha,beta
852 853
      integer(kind=C_intptr_T)        :: a, b, c
#ifdef WITH_GPU_VERSION
Pavel Kus's avatar
Pavel Kus committed
854
      call cublas_cgemm_c(cublasHandle, cta, ctb, m, n, k, alpha, a, lda, b, ldb, beta, c,ldc)
855 856 857
#endif
    end subroutine cublas_cgemm

858 859 860 861 862
    subroutine cublas_ztrmm(side, uplo, trans, diag, m, n, alpha, a, lda, b, ldb)

      use iso_c_binding

      implicit none
863
      character(1,C_CHAR),value       :: side, uplo, trans, diag
864 865
      integer(kind=C_INT)             :: m,n
      integer(kind=C_INT), intent(in) :: lda,ldb
866
      complex(kind=C_DOUBLE_COMPLEX)          :: alpha
867 868
      integer(kind=C_intptr_T)        :: a, b
#ifdef WITH_GPU_VERSION
Pavel Kus's avatar
Pavel Kus committed
869
      call cublas_ztrmm_c(cublasHandle, side, uplo, trans, diag, m, n, alpha, a, lda, b, ldb)
870 871 872
#endif
    end subroutine cublas_ztrmm

873 874 875 876 877
    subroutine cublas_ctrmm(side, uplo, trans, diag, m, n, alpha, a, lda, b, ldb)

      use iso_c_binding

      implicit none
878
      character(1,C_CHAR),value       :: side, uplo, trans, diag
879 880
      integer(kind=C_INT)             :: m,n
      integer(kind=C_INT), intent(in) :: lda,ldb
881
      complex(kind=C_FLOAT_COMPLEX)           :: alpha
882 883
      integer(kind=C_intptr_T)        :: a, b
#ifdef WITH_GPU_VERSION
Pavel Kus's avatar
Pavel Kus committed
884
      call cublas_ctrmm_c(cublasHandle, side, uplo, trans, diag, m, n, alpha, a, lda, b, ldb)
885 886
#endif
    end subroutine cublas_ctrmm
887

888 889 890 891 892 893 894 895 896 897
    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
898
      call cublas_dgemv_c(cublasHandle, cta, m, n, alpha, a, lda, x, incx, beta, y, incy)
899 900 901 902 903 904 905 906 907 908 909 910 911
#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
912
      call cublas_sgemv_c(cublasHandle, cta, m, n, alpha, a, lda, x, incx, beta, y, incy)
913 914 915
#endif
    end subroutine cublas_sgemv

916 917 918 919 920 921 922
    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
923
      complex(kind=C_DOUBLE_COMPLEX)             :: alpha,beta
924 925
      integer(kind=C_intptr_T)        :: a, x, y
#ifdef WITH_GPU_VERSION
Pavel Kus's avatar
Pavel Kus committed
926
      call cublas_zgemv_c(cublasHandle, cta, m, n, alpha, a, lda, x, incx, beta, y, incy)
927 928 929 930 931 932 933 934 935 936
#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
937
      complex(kind=C_FLOAT_COMPLEX)              :: alpha,beta
938 939
      integer(kind=C_intptr_T)        :: a, x, y
#ifdef WITH_GPU_VERSION
Pavel Kus's avatar
Pavel Kus committed
940
      call cublas_cgemv_c(cublasHandle, cta, m, n, alpha, a, lda, x, incx, beta, y, incy)
941 942 943
#endif
    end subroutine cublas_cgemv

944

Pavel Kus's avatar
Pavel Kus committed
945 946 947 948 949 950 951 952 953 954 955 956 957 958 959 960 961 962 963 964 965 966 967 968 969 970 971 972 973 974 975 976 977 978 979 980 981 982 983 984 985 986 987 988 989 990 991 992 993 994 995 996 997 998 999
!     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
1000 1001


1002
end module cuda_functions