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

  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
181
      integer(kind=c_intptr_t), intent(in), value    :: size
182
183
184
185
186
187
188
189
190
191
192
193
194
195
      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

196
      integer(kind=C_intptr_T), value                :: dst
197
      integer(kind=c_intptr_t), intent(in), value    :: dpitch
198
      integer(kind=C_intptr_T), value                :: src
199
200
201
      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
202
203
      integer(kind=C_INT), intent(in), value         :: dir
      integer(kind=C_INT)                            :: istat
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230

    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
231
      integer(kind=c_intptr_t), intent(in), value   :: width_height
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
      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
247
      integer(kind=c_intptr_t), intent(in), value  :: size
248
249
250
251
252
253
254
      integer(kind=C_INT)                        :: istat

    end function cuda_memset_c
  end interface

  ! cuBLAS
  interface
Pavel Kus's avatar
Pavel Kus committed
255
256
    subroutine cublas_dgemm_c(handle, cta, ctb, m, n, k, alpha, a, lda, b, ldb, beta, c, ldc) &
                              bind(C,name='cublasDgemm_elpa_wrapper')
257
258
259
260
261
262
263
264
265

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

268
269
270
    end subroutine cublas_dgemm_c
  end interface

271
  interface
Pavel Kus's avatar
Pavel Kus committed
272
273
    subroutine cublas_sgemm_c(handle, cta, ctb, m, n, k, alpha, a, lda, b, ldb, beta, c, ldc) &
                              bind(C,name='cublasSgemm_elpa_wrapper')
274
275
276
277
278
279
280
281
282

      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
283
284
      integer(kind=C_intptr_T), value         :: handle
      
285
286
287
    end subroutine cublas_sgemm_c
  end interface

288
  interface
Pavel Kus's avatar
Pavel Kus committed
289
290
    subroutine cublas_dtrmm_c(handle, side, uplo, trans, diag, m, n, alpha, a, lda, b, ldb) &
                              bind(C,name='cublasDtrmm_elpa_wrapper')
291
292
293
294
295
296
297
298
299

      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
300
301
      integer(kind=C_intptr_T), value         :: handle
      
302
303
304
    end subroutine cublas_dtrmm_c
  end interface

305
  interface
Pavel Kus's avatar
Pavel Kus committed
306
307
    subroutine cublas_strmm_c(handle, side, uplo, trans, diag, m, n, alpha, a, lda, b, ldb) &
                              bind(C,name='cublasStrmm_elpa_wrapper')
308
309
310
311
312
313
314
315
316

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

319
320
321
    end subroutine cublas_strmm_c
  end interface

322
  interface
Pavel Kus's avatar
Pavel Kus committed
323
324
    subroutine cublas_zgemm_c(handle, cta, ctb, m, n, k, alpha, a, lda, b, ldb, beta, c,ldc) &
                              bind(C,name='cublasZgemm_elpa_wrapper')
325
326
327
328
329
330
331

      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
332
      complex(kind=C_DOUBLE_COMPLEX),value           :: alpha,beta
333
      integer(kind=C_intptr_T), value        :: a, b, c
Pavel Kus's avatar
Pavel Kus committed
334
      integer(kind=C_intptr_T), value         :: handle
335
336
337
338

    end subroutine cublas_zgemm_c
  end interface

339
  interface
Pavel Kus's avatar
Pavel Kus committed
340
341
    subroutine cublas_cgemm_c(handle, cta, ctb, m, n, k, alpha, a, lda, b, ldb, beta, c,ldc) &
                              bind(C,name='cublasCgemm_elpa_wrapper')
342
343
344
345
346
347
348

      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
349
      complex(kind=C_FLOAT_COMPLEX),value            :: alpha,beta
350
      integer(kind=C_intptr_T), value        :: a, b, c
Pavel Kus's avatar
Pavel Kus committed
351
      integer(kind=C_intptr_T), value         :: handle
352
353
354
355

    end subroutine cublas_cgemm_c
  end interface

356
  interface
Pavel Kus's avatar
Pavel Kus committed
357
358
    subroutine cublas_ztrmm_c(handle, side, uplo, trans, diag, m, n, alpha, a, lda, b, ldb) &
                              bind(C,name='cublasZtrmm_elpa_wrapper')
359
360
361
362
363
364
365

      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
366
      complex(kind=C_DOUBLE_COMPLEX), value          :: alpha
367
      integer(kind=C_intptr_T), value        :: a, b
Pavel Kus's avatar
Pavel Kus committed
368
      integer(kind=C_intptr_T), value         :: handle
369
370
371
372

    end subroutine cublas_ztrmm_c
  end interface

373
  interface
Pavel Kus's avatar
Pavel Kus committed
374
375
    subroutine cublas_ctrmm_c(handle, side, uplo, trans, diag, m, n, alpha, a, lda, b, ldb) &
                              bind(C,name='cublasCtrmm_elpa_wrapper')
376
377
378
379
380
381
382

      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
383
      complex(kind=C_FLOAT_COMPLEX), value           :: alpha
384
      integer(kind=C_intptr_T), value        :: a, b
Pavel Kus's avatar
Pavel Kus committed
385
      integer(kind=C_intptr_T), value         :: handle
386
387
388
389

    end subroutine cublas_ctrmm_c
  end interface

390
  interface
Pavel Kus's avatar
Pavel Kus committed
391
392
    subroutine cublas_dgemv_c(handle, cta, m, n, alpha, a, lda, x, incx, beta, y, incy) &
                              bind(C,name='cublasDgemv_elpa_wrapper')
393
394
395
396
397
398
399
400
401

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

404
405
406
407
    end subroutine cublas_dgemv_c
  end interface

  interface
Pavel Kus's avatar
Pavel Kus committed
408
409
    subroutine cublas_sgemv_c(handle, cta, m, n, alpha, a, lda, x, incx, beta, y, incy) &
                              bind(C,name='cublasSgemv_elpa_wrapper')
410
411
412
413
414
415
416
417
418

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

421
422
423
    end subroutine cublas_sgemv_c
  end interface

424
  interface
Pavel Kus's avatar
Pavel Kus committed
425
426
    subroutine cublas_zgemv_c(handle, cta, m, n, alpha, a, lda, x, incx, beta, y, incy) &
                              bind(C,name='cublasZgemv_elpa_wrapper')
427
428
429
430
431
432
433

      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
434
      complex(kind=C_DOUBLE_COMPLEX),value               :: alpha,beta
435
      integer(kind=C_intptr_T), value         :: a, x, y
Pavel Kus's avatar
Pavel Kus committed
436
437
      integer(kind=C_intptr_T), value         :: handle

438
439
440
441
    end subroutine cublas_zgemv_c
  end interface

  interface
Pavel Kus's avatar
Pavel Kus committed
442
443
    subroutine cublas_cgemv_c(handle, cta, m, n, alpha, a, lda, x, incx, beta, y, incy) &
                              bind(C,name='cublasCgemv_elpa_wrapper')
444
445
446
447
448
449
450

      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
451
      complex(kind=C_FLOAT_COMPLEX),value                :: alpha,beta
452
      integer(kind=C_intptr_T), value         :: a, x, y
Pavel Kus's avatar
Pavel Kus committed
453
      integer(kind=C_intptr_T), value         :: handle
Pavel Kus's avatar
Pavel Kus committed
454

Pavel Kus's avatar
Pavel Kus committed
455
    end subroutine cublas_cgemv_c
Pavel Kus's avatar
Pavel Kus committed
456
457
458
  end interface


Pavel Kus's avatar
Pavel Kus committed
459
  contains
Pavel Kus's avatar
Pavel Kus committed
460

Pavel Kus's avatar
Pavel Kus committed
461
    ! functions to set and query the CUDA devices
Pavel Kus's avatar
Pavel Kus committed
462

Pavel Kus's avatar
Pavel Kus committed
463
464
465
   function cublas_create(handle) result(success)
     use iso_c_binding
     implicit none
466

Pavel Kus's avatar
Pavel Kus committed
467
468
469
470
471
472
473
474
     integer(kind=C_intptr_t)                  :: handle
     logical                                   :: success
#ifdef WITH_GPU_VERSION
     success = cublas_create_c(handle) /= 0
#else
     success = .true.
#endif
   end function
475

Pavel Kus's avatar
Pavel Kus committed
476
477
478
   function cublas_destroy(handle) result(success)
     use iso_c_binding
     implicit none
479

Pavel Kus's avatar
Pavel Kus committed
480
481
482
483
484
485
486
487
488
     integer(kind=C_intptr_t)                  :: handle
     logical                                   :: success
#ifdef WITH_GPU_VERSION
     success = cublas_destroy_c(handle) /= 0
#else
     success = .true.
#endif
   end function
    
489
490
491
492
493
    function cuda_setdevice(n) result(success)
      use iso_c_binding

      implicit none

494
495
      integer(kind=ik), intent(in)  :: n
      logical                       :: success
496
497
498
499
500
501
502
503
504
505
506
#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

507
      integer(kind=ik)     :: n
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
      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
528
      success = .true.
529
530
531
532
533
534
535
536
537
538
#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
539
     integer(kind=c_intptr_t), intent(in)        :: width_height
540
541
542
543
     logical                                   :: success
#ifdef WITH_GPU_VERSION
     success = cuda_malloc_c(a, width_height) /= 0
#else
544
     success = .true.
545
546
547
548
549
550
551
552
553
554
555
556
557
#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
558
     success = .true.
559
560
561
562
563
564
565
566
567
568
#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
569
   integer(kind=ik)                        :: val
570
   integer(kind=c_intptr_t), intent(in)      :: size
571
572
573
574
   integer(kind=C_INT)                     :: istat

   logical :: success
#ifdef WITH_GPU_VERSION
575
   success= cuda_memset_c(a, int(val,kind=c_int), int(size,kind=c_intptr_t)) /=0
576
#else
577
   success = .true.
578
579
580
581
582
583
584
585
#endif
 end function cuda_memset

 ! functions to memcopy CUDA memory

 function cuda_memcpyDeviceToDevice() result(flag)
   use iso_c_binding
   implicit none
586
   integer(kind=ik) :: flag
587
588
589
590
591
592
593
594
595
#ifdef WITH_GPU_VERSION
   flag = int(cuda_memcpyDeviceToDevice_c())
#else
   flag = 0
#endif
 end function

 function cuda_memcpyHostToDevice() result(flag)
   use iso_c_binding
596
   use precision
597
   implicit none
598
   integer(kind=ik) :: flag
599
600
601
602
603
604
605
606
607
#ifdef WITH_GPU_VERSION
   flag = int(cuda_memcpyHostToDevice_c())
#else
   flag = 0
#endif
 end function

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

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

 function cuda_hostRegisterMapped() result(flag)
   use iso_c_binding
632
   use precision
633
   implicit none
634
   integer(kind=ik) :: flag
635
636
637
638
639
640
641
642
643
644
645
646
647
648
#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
649
      integer(kind=c_intptr_t), intent(in)    :: size
650
651
      integer(kind=C_INT), intent(in)       :: dir
      logical :: success
652

653
654
655
#ifdef WITH_GPU_VERSION
        success = cuda_memcpy_c(dst, src, size, dir) /= 0
#else
656
        success = .true.
657
658
659
660
661
662
663
664
665
666
#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
667
      integer(kind=c_intptr_t), intent(in) :: dpitch
668
      integer(kind=C_intptr_T)           :: src
669
670
671
      integer(kind=c_intptr_t), intent(in) :: spitch
      integer(kind=c_intptr_t), intent(in) :: width
      integer(kind=c_intptr_t), intent(in) :: height
672
673
674
675
676
      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
677
      success = .true.
678
679
680
681
682
683
684
685
#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
686
      character(1,C_CHAR),value       :: cta, ctb
687
688
689
690
691
      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
692
      call cublas_dgemm_c(cublasHandle, cta, ctb, m, n, k, alpha, a, lda, b, ldb, beta, c, ldc)
693
694
695
#endif
    end subroutine cublas_dgemm

696
697
698
699
    subroutine cublas_sgemm(cta, ctb, m, n, k, alpha, a, lda, b, ldb, beta, c, ldc)
      use iso_c_binding

      implicit none
700
      character(1,C_CHAR),value       :: cta, ctb
701
702
703
704
705
      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
706
      call cublas_sgemm_c(cublasHandle, cta, ctb, m, n, k, alpha, a, lda, b, ldb, beta, c, ldc)
707
708
709
#endif
    end subroutine cublas_sgemm

710
711
712
713
714
    subroutine cublas_dtrmm(side, uplo, trans, diag, m, n, alpha, a, lda, b, ldb)

      use iso_c_binding

      implicit none
715
      character(1,C_CHAR),value       :: side, uplo, trans, diag
716
717
718
719
720
      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
721
      call cublas_dtrmm_c(cublasHandle, side, uplo, trans, diag, m, n, alpha, a, lda, b, ldb)
722
723
724
#endif
    end subroutine cublas_dtrmm

725
726
727
728
729
    subroutine cublas_strmm(side, uplo, trans, diag, m, n, alpha, a, lda, b, ldb)

      use iso_c_binding

      implicit none
730
      character(1,C_CHAR),value       :: side, uplo, trans, diag
731
732
733
734
735
      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
736
      call cublas_strmm_c(cublasHandle, side, uplo, trans, diag, m, n, alpha, a, lda, b, ldb)
737
738
739
#endif
    end subroutine cublas_strmm

740
741
742
743
744
    subroutine cublas_zgemm(cta, ctb, m, n, k, alpha, a, lda, b, ldb, beta, c,ldc)

      use iso_c_binding

      implicit none
745
      character(1,C_CHAR),value       :: cta, ctb
746
747
      integer(kind=C_INT)             :: m,n,k
      integer(kind=C_INT), intent(in) :: lda,ldb,ldc
748
      complex(kind=C_DOUBLE_COMPLEX)          :: alpha,beta
749
750
      integer(kind=C_intptr_T)        :: a, b, c
#ifdef WITH_GPU_VERSION
Pavel Kus's avatar
Pavel Kus committed
751
      call cublas_zgemm_c(cublasHandle, cta, ctb, m, n, k, alpha, a, lda, b, ldb, beta, c,ldc)
752
753
754
#endif
    end subroutine cublas_zgemm

755
756
757
758
759
    subroutine cublas_cgemm(cta, ctb, m, n, k, alpha, a, lda, b, ldb, beta, c,ldc)

      use iso_c_binding

      implicit none
760
      character(1,C_CHAR),value       :: cta, ctb
761
762
      integer(kind=C_INT)             :: m,n,k
      integer(kind=C_INT), intent(in) :: lda,ldb,ldc
763
      complex(kind=C_FLOAT_COMPLEX)           :: alpha,beta
764
765
      integer(kind=C_intptr_T)        :: a, b, c
#ifdef WITH_GPU_VERSION
Pavel Kus's avatar
Pavel Kus committed
766
      call cublas_cgemm_c(cublasHandle, cta, ctb, m, n, k, alpha, a, lda, b, ldb, beta, c,ldc)
767
768
769
#endif
    end subroutine cublas_cgemm

770
771
772
773
774
    subroutine cublas_ztrmm(side, uplo, trans, diag, m, n, alpha, a, lda, b, ldb)

      use iso_c_binding

      implicit none
775
      character(1,C_CHAR),value       :: side, uplo, trans, diag
776
777
      integer(kind=C_INT)             :: m,n
      integer(kind=C_INT), intent(in) :: lda,ldb
778
      complex(kind=C_DOUBLE_COMPLEX)          :: alpha
779
780
      integer(kind=C_intptr_T)        :: a, b
#ifdef WITH_GPU_VERSION
Pavel Kus's avatar
Pavel Kus committed
781
      call cublas_ztrmm_c(cublasHandle, side, uplo, trans, diag, m, n, alpha, a, lda, b, ldb)
782
783
784
#endif
    end subroutine cublas_ztrmm

785
786
787
788
789
    subroutine cublas_ctrmm(side, uplo, trans, diag, m, n, alpha, a, lda, b, ldb)

      use iso_c_binding

      implicit none
790
      character(1,C_CHAR),value       :: side, uplo, trans, diag
791
792
      integer(kind=C_INT)             :: m,n
      integer(kind=C_INT), intent(in) :: lda,ldb
793
      complex(kind=C_FLOAT_COMPLEX)           :: alpha
794
795
      integer(kind=C_intptr_T)        :: a, b
#ifdef WITH_GPU_VERSION
Pavel Kus's avatar
Pavel Kus committed
796
      call cublas_ctrmm_c(cublasHandle, side, uplo, trans, diag, m, n, alpha, a, lda, b, ldb)
797
798
#endif
    end subroutine cublas_ctrmm
799

800
801
802
803
804
805
806
807
808
809
    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
810
      call cublas_dgemv_c(cublasHandle, cta, m, n, alpha, a, lda, x, incx, beta, y, incy)
811
812
813
814
815
816
817
818
819
820
821
822
823
#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
824
      call cublas_sgemv_c(cublasHandle, cta, m, n, alpha, a, lda, x, incx, beta, y, incy)
825
826
827
#endif
    end subroutine cublas_sgemv

828
829
830
831
832
833
834
    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
835
      complex(kind=C_DOUBLE_COMPLEX)             :: alpha,beta
836
837
      integer(kind=C_intptr_T)        :: a, x, y
#ifdef WITH_GPU_VERSION
Pavel Kus's avatar
Pavel Kus committed
838
      call cublas_zgemv_c(cublasHandle, cta, m, n, alpha, a, lda, x, incx, beta, y, incy)
839
840
841
842
843
844
845
846
847
848
#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
849
      complex(kind=C_FLOAT_COMPLEX)              :: alpha,beta
850
851
      integer(kind=C_intptr_T)        :: a, x, y
#ifdef WITH_GPU_VERSION
Pavel Kus's avatar
Pavel Kus committed
852
      call cublas_cgemv_c(cublasHandle, cta, m, n, alpha, a, lda, x, incx, beta, y, incy)
853
854
855
#endif
    end subroutine cublas_cgemv

856

Pavel Kus's avatar
Pavel Kus committed
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
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
!     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
912
913


914
end module cuda_functions