-
Notifications
You must be signed in to change notification settings - Fork 15
/
tensor_algebra_gpu_nvidia.cu
7196 lines (7062 loc) · 356 KB
/
tensor_algebra_gpu_nvidia.cu
1
2
3
4
5
6
7
8
9
10
11
12
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
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
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
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
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
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
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
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
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
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
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
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
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
1000
/** Tensor Algebra Library for NVidia GPU: NV-TAL (CUDA based).
AUTHOR: Dmitry I. Lyakh (Liakh): [email protected], [email protected]
REVISION: 2021/12/20
Copyright (C) 2014-2022 Dmitry I. Lyakh (Liakh)
Copyright (C) 2014-2022 Oak Ridge National Laboratory (UT-Battelle)
LICENSE: BSD 3-Clause
-------------------------------------------------------------------
OPTIONS:
# -D CUDA_ARCH=350: target GPU compute capability (default is 130);
# -D NO_GPU: disables GPU usage;
# -D NO_BLAS: disables cuBLAS calls, they will be replaced by in-house routines (slower);
# -D USE_CUTT: enables an optimized tensor transpose via the cuTT library;
# -D DEBUG_GPU: collection of debugging information will be activated;
NOTES:
# Minimal required compute capability is 1.1 (1.3 for double precision).
# cuBLAS.v2 is required when BLAS is enabled.
# Non-blocking tensor algebra functions carry an additional output argument <cuda_task> (task handle).
# Non-blocking tensor algebra functions carry an additional input argument <coherence_ctrl>
which controls the tensor data consistency synchronization accross different devices
after the tensor operation has completed successfully.
FOR DEVELOPERS ONLY:
# Currently used device resources:
- Global memory pointer (any device);
- Argument buffer entry handle (any device);
- Multi-index entry * (Host pinned memory, entry length = MAX_TENSOR_RANK);
- GPU constant-memory entry handle (Nvidia GPU);
- CUDA stream handle (Nvidia GPU);
- CUDA event handle (Nvidia GPU).
# A life cycle of a C object (for example, tensBlck_t):
a) Allocate memory for the object itself, if needed: Suffix _alloc or _create (includes cleaning);
b) Clean (initialize to null) an allocated (empty) object: Suffix _clean (normally included in _create);
c) Construct (define or redefine) an existing object (resources will be acquired/released): Suffix _construct;
d) Destruct a defined object (resources will be released, the object will be reset to clean): Suffix _destruct;
e) Free the memory occupied by an object: Suffix _free or _destroy (may include _destruct, if needed).
Thus, as a rule, the device resource acquisition/release occurs solely in _construct and _destruct functions.
# A state of a C object:
a) Undefined: After the memory allocation (either dynamic or static);
b) Defined-empty (clean): After cleaning or destruction (dynamic object creation produces a clean object);
c) Defined to a value (value-defined): After construction;
d) Dead: After memory deallocation (if it was allocated dynamically).
# Resource acquisition/release:
- Tensor block constructor/destructor acquires/releases global memory resources, including
both pointers and buffer entries, as well as multi-index bank entries (pinned Host memory).
- CUDA task constructor/destructor acquires/releases CUDA resources (stream, events).
- Tensor operation scheduling functions acquire GPU global memory resources,
GPU constant memory resources, Host pinned multi-index entries.
- CUDA task completion/error query functions release GPU global memory resources,
GPU constant memory resources, Host pinned multi-index entries.
- Coherence control is only applied to successfully finished CUDA tasks.
# Functions which construct tensor blocks or perform asynchronous operations on them
allocate resources (global/constant memory, etc). In case the corresponding resource
allocator returns TRY_LATER or DEVICE_UNABLE (or an error), the corresponding function
must clean the partially created tensor block or the CUDA task before returning:
The corresponding object will be kept in its initial state if no SUCCESS.
# Some CUDA kernels operating on two or more arguments assume no aliases
for GPU pointers (__restrict__). Check each specific operation to see whether
it is ok for the two tensor arguments to refer to the same tensor body.
TO BE FIXED:
# In tensor operation scheduling functions, if a scheduling error occurs after
the data transfer or CUDA kernel has been scheduled, the CUDA task finalization
must not begin until the partially scheduled CUDA task has completed on GPU.
Insert cudaStreamSynchronize in the finalization procedure.
# Invoke cudaDeviceCanAccessPeer() in tensor operations to check whether
two devices of the same kind can access each other's memory.
# Account for implicit data transfers to/from peer GPUs in their statistics.
# User-provided Alpha factors for gpu_tensor_block_contract() and
gpu_tensor_block_add() reside on Host, thus requiring a slab in GPU
memory (either global or constant) as a temporary for BLAS references.
**/
#include "device_algebra.h"
#include "tensor_algebra.h"
#include "mem_manager.h"
#include "talsh_complex.h"
#include <cstdio>
#include <cstdlib>
#include <ctime>
#ifndef NO_GPU
//PARAMETERS:
#define GPU_DEBUG_DUMP_SIZE 128 //size of the GPU debug dump (int array)
#endif /*NO_GPU*/
//----------------------------------------------------------------------
//FUNCTION PROTOTYPES:
// LOCAL (PRIVATE):
static int prmn_convert(int n, const int *o2n, int *n2o);
static int non_trivial_prmn(int n, const int *prm);
#ifndef NO_GPU
static int cuda_stream_get(int gpu_num, int * cuda_stream_handle);
static int cuda_stream_release(int gpu_num, int cuda_stream_handle);
static cudaStream_t * cuda_stream_ptr(int gpu_num, int cuda_stream_handle);
static int cuda_event_get(int gpu_num, int * cuda_event_handle);
static int cuda_event_release(int gpu_num, int cuda_event_handle);
static cudaEvent_t * cuda_event_ptr(int gpu_num, int cuda_event_handle);
static void limit_cuda_blocks2d(int max_blocks, int *bx, int *by);
static int tens_op_best_gpu(const tensBlck_t *tens0 = NULL, const tensBlck_t *tens1 = NULL, const tensBlck_t *tens2 = NULL);
static int cuda_task_set_arg(cudaTask_t *cuda_task, unsigned int arg_num, tensBlck_t *tens_p);
static int cuda_task_set_prefactor(cudaTask_t *cuda_task, talshComplex4 prefactor);
static int cuda_task_set_prefactor(cudaTask_t *cuda_task, talshComplex8 prefactor);
static int cuda_task_record(cudaTask_t *cuda_task, unsigned int coh_ctrl, unsigned int err_code = 0);
static int cuda_task_finalize(cudaTask_t *cuda_task);
// CUDA KERNELS:
template <typename T>
__global__ void gpu_array_norm2__(size_t tsize, const T * __restrict__ arr, volatile double * bnorm2);
template <typename T>
__global__ void gpu_array_init__(size_t tsize, T * arr, T val);
template <typename T>
__global__ void gpu_scalar_multiply__(const T * left_arg, const T * right_arg, T * dest_arg, T alpha,
int left_conj = 0, int right_conj = 0);
template <typename T>
__global__ void gpu_array_scale__(size_t tsize, T * arr, T alpha);
template <typename T>
__global__ void gpu_array_add__(size_t tsize, T * __restrict__ arr0, const T * __restrict__ arr1,
T alpha, int left_conj = 0);
template <typename T>
__global__ void gpu_array_add__(size_t tsize, T * __restrict__ arr0, const T * __restrict__ arr1, const T * __restrict__ scalar,
T alpha, int left_conj = 0);
template <typename T>
__global__ void gpu_array_dot_product__(size_t tsize, const T * arr1, const T * arr2, volatile T * dprod,
T alpha, int left_conj = 0, int right_conj = 0);
template <typename T>
__global__ void gpu_array_product__(size_t tsize1, const T * arr1, size_t tsize2, const T * arr2, T * arr0,
T alpha, int left_conj = 0, int right_conj = 0);
template <typename T>
__global__ void gpu_tensor_block_add_dlf__(int dmo, int drc, int dim_num, int const_args_pos,
const T * __restrict__ tens_in, T * __restrict__ tens_out);
template <typename T>
__global__ void gpu_tensor_block_copy_dlf__(int dmo, int drc, int dim_num, int const_args_pos,
const T * __restrict__ tens_in, T * __restrict__ tens_out);
template <typename T>
__global__ void gpu_tensor_block_copy_cmplx_split_in_dlf__(int dmo, int drc, int dim_num, int const_args_pos,
const T * __restrict__ tens_in, T * __restrict__ tens_out);
template <typename T>
__global__ void gpu_tensor_block_copy_cmplx_split_out_dlf__(int dmo, int drc, int dim_num, int const_args_pos,
const T * __restrict__ tens_in, T * __restrict__ tens_out);
template <typename T>
__global__ void gpu_tensor_block_copy_scatter_dlf__(int dmo, int drc, int dim_num, int const_args_pos,
const T * __restrict__ tens_in, T * __restrict__ tens_out);
template <typename T>
__global__ void gpu_matrix_multiply_tn__(size_t ll, size_t lr, size_t lc, const T * arg1, const T * arg2, T * arg0, T alpha);
template <typename T>
__global__ void gpu_matrix_multiply_nt__(size_t ll, size_t lr, size_t lc, const T * arg1, const T * arg2, T * arg0, T alpha);
template <typename T>
__global__ void gpu_matrix_multiply_nn__(size_t ll, size_t lr, size_t lc, const T * arg1, const T * arg2, T * arg0, T alpha);
template <typename T>
__global__ void gpu_matrix_multiply_tt__(size_t ll, size_t lr, size_t lc, const T * arg1, const T * arg2, T * arg0, T alpha);
#endif /*NO_GPU*/
//---------------------------------------------------------------------------------------------------------------------------
//PARAMETERS:
static int VERBOSE=1; //verbosity for error messages
static int DEBUG=0; //debugging mode
#ifndef NO_GPU
//GLOBAL DATA:
// GPU control on the current MPI process:
int gpu_up[MAX_GPUS_PER_NODE]={GPU_OFF}; //GPU_OFF(0): GPU is disabled; GPU_MINE(1): GPU is enabled; GPU_MINE_CUBLAS(2): GPU is BLAS enabled
cudaDeviceProp gpu_prop[MAX_GPUS_PER_NODE]; //properties of all GPUs present on the node
talsh_stats_t gpu_stats[MAX_GPUS_PER_NODE]; //runtime statistics for all GPUs present on the node
#ifndef NO_BLAS
// Infrastructure for CUBLAS:
cublasHandle_t cublas_handle[MAX_GPUS_PER_NODE]; //each GPU present on a node obtains its own cuBLAS context handle
#endif /*NO_BLAS*/
#ifdef USE_CUTENSOR
// Infrastructure for cuTensor:
cutensorHandle_t cutensor_handle[MAX_GPUS_PER_NODE]; //each GPU present on a node obtains its own cuTensor context handle
void * cutensor_workspace[MAX_GPUS_PER_NODE] = {NULL}; //cuTensor workspace (in GPU memory)
size_t cutensor_worksize[MAX_GPUS_PER_NODE] = {0}; //cuTensor workspace size
const size_t CUTENSOR_WORKSPACE_SIZE = 128 * 1048576; //default cuTensor workspace size
#endif /*USE_CUTENSOR*/
// Slabs for the GPU asynchronous resources:
// CUDA stream handles:
cudaStream_t CUDAStreamBank[MAX_GPUS_PER_NODE][MAX_CUDA_TASKS]; //pre-allocated CUDA stream handles (for each CUDA device)
int CUDAStreamFreeHandle[MAX_GPUS_PER_NODE][MAX_CUDA_TASKS]; //free CUDA stream handles
int CUDAStreamFFE[MAX_GPUS_PER_NODE]; //number of free handles left in CUDAStreamFreeHandle
// CUDA event handles:
cudaEvent_t CUDAEventBank[MAX_GPUS_PER_NODE][MAX_CUDA_EVENTS]; //pre-allocated CUDA event handles (for each CUDA device)
int CUDAEventFreeHandle[MAX_GPUS_PER_NODE][MAX_CUDA_EVENTS]; //free CUDA event handles
int CUDAEventFFE[MAX_GPUS_PER_NODE]; //number of free handles left in CUDAEventFreeHandle
// Mapped slab of tensor operation prefactors for GPU usage:
slab_t prefactors; //mapped slab of prefactors
void * gpu_prefs_base_ptr; //mapped device pointer of the slab base
// Slab of GPU constant memory arguments for each GPU (managed by "mem_manager.cpp"):
__constant__ int const_args_dims[MAX_GPU_ARGS][MAX_TENSOR_RANK]; //storage for device constant memory arguments: dimension extents
__constant__ int const_args_prmn[MAX_GPU_ARGS][MAX_TENSOR_RANK]; //storage for device constant memory arguments: permutation
// GPU error control and debugging for each GPU:
__device__ int gpu_error_count=0; //total number of CUDA errors registered on device till the current moment
__device__ int gpu_debug_dump[GPU_DEBUG_DUMP_SIZE]; //debug dump
// Global CUDA event recording policy:
static int PRINT_TIMING=1; //non-zero value enables time printing statements
// Infrastructure for function <gpu_tensor_block_copy_dlf> (blocking and non-blocking):
#ifdef USE_CUTT
static int TRANS_SHMEM=EFF_TRN_ON_CUTT; //switch between shared-memory tensor transpose and scatter tensor transpose
#else
static int TRANS_SHMEM=EFF_TRN_ON; //switch between shared-memory tensor transpose and scatter tensor transpose
#endif /*USE_CUTT*/
// Infrastructure for <gpu_tensor_block_contract_dlf> (non-blocking):
#ifndef NO_BLAS
static int DISABLE_BLAS=0; //non-zero value will disable cuBLAS usage (if it had been cuBLAS compiled/linked)
#else
static int DISABLE_BLAS=1; //non-zero value will disable cuBLAS usage (if it had been cuBLAS compiled/linked)
#endif /*NO_BLAS*/
cudaTask_t * LastTask[MAX_GPUS_PER_NODE]; //last CUDA task successfully scheduled on each GPU
float h_sgemm_beta_one=1.0f;
float h_sgemm_beta_zero=0.0f;
double h_dgemm_beta_one=1.0;
double h_dgemm_beta_zero=0.0;
cuComplex h_cgemm_beta_one={1.0f,0.0f};
cuComplex h_cgemm_beta_zero={0.0f,0.0f};
cuDoubleComplex h_zgemm_beta_one={1.0,0.0};
cuDoubleComplex h_zgemm_beta_zero={0.0,0.0};
__constant__ float sgemm_alpha_plus=1.0f; //default alpha constant for SGEMM
__constant__ float sgemm_alpha_minus=-1.0f; //default alpha constant for SGEMM
__constant__ float sgemm_beta_one=1.0f; //default beta constant SGEMM
__constant__ float sgemm_beta_zero=0.0f; //zero beta constant SGEMM
__constant__ double dgemm_alpha_plus=1.0; //default alpha constant for DGEMM
__constant__ double dgemm_alpha_minus=-1.0; //default alpha constant for DGEMM
__constant__ double dgemm_beta_one=1.0; //default beta constant DGEMM
__constant__ double dgemm_beta_zero=0.0; //zero beta constant DGEMM
__constant__ cuComplex cgemm_alpha_plus={1.0f,0.0f}; //default alpha constant CGEMM
__constant__ cuComplex cgemm_alpha_minus={-1.0f,0.0f}; //default alpha constant CGEMM
__constant__ cuComplex cgemm_beta_one={1.0f,0.0f}; //default beta constant CGEMM
__constant__ cuComplex cgemm_beta_zero={0.0f,0.0f}; //zero beta constant CGEMM
__constant__ cuDoubleComplex zgemm_alpha_plus={1.0,0.0}; //default alpha constant ZGEMM
__constant__ cuDoubleComplex zgemm_alpha_minus={-1.0,0.0}; //default alpha constant ZGEMM
__constant__ cuDoubleComplex zgemm_beta_one={1.0,0.0}; //default beta constant ZGEMM
__constant__ cuDoubleComplex zgemm_beta_zero={0.0,0.0}; //zero beta constant ZGEMM
// Infrastructure for kernels <gpu_array_norm2__>:
__device__ int norm2_wr_lock=0; //write lock shared by all <gpu_array_norm2__> running on GPU
// Infrastructure for kernels <gpu_array_dot_product__>:
__device__ int dot_product_wr_lock=0; //write lock shared by all <gpu_array_dot_product__> running on GPU
#endif /*NO_GPU*/
//--------------------------------------------------------------------------------------------------------------
#ifndef NO_GPU
//CUDA KERNELS:
// SUM OF THE SQUARES OF ABSOLUTE VALUES OF ALL ARRAY ELEMENTS:
// REAL:
template <typename T>
__global__ void gpu_array_norm2__(size_t tsize, const T * __restrict__ arr, volatile double * bnorm2)
/** Computes the squared 2-norm of array arr(0:tsize-1)
INPUT:
# tsize - size of the array;
# arr(0:tsize-1) - array;
# bnorm2 - must be zero on entrance (resides on device as well);
OUTPUT:
# bnorm2 - squared 2-norm of the array (resides on device as well);
**/
{
size_t i,n;
double _thread_norm2;
extern __shared__ double thread_norms2[]; //size = blockDim.x*sizeof(double) Bytes per thread block
n=gridDim.x*blockDim.x; _thread_norm2=0.0;
for(i=blockIdx.x*blockDim.x+threadIdx.x;i<tsize;i+=n) _thread_norm2+=arr[i]*arr[i];
thread_norms2[threadIdx.x]=_thread_norm2;
__syncthreads();
if(threadIdx.x == 0){ //global reduction among thread blocks
_thread_norm2=thread_norms2[0];
for(i=1;i<blockDim.x;i++) _thread_norm2+=thread_norms2[i];
i=1; while(i == 1){i=atomicMax(&norm2_wr_lock,1);} //waiting for the lock to unlock, then lock
__threadfence();
*bnorm2+=_thread_norm2;
__threadfence();
i=atomicExch(&norm2_wr_lock,0); //unlock
}
__syncthreads();
return;
}
// COMPLEX4:
template <>
__global__ void gpu_array_norm2__<talshComplex4>(size_t tsize, const talshComplex4 * __restrict__ arr, volatile double * bnorm2)
/** Computes the squared 2-norm of array arr(0:tsize-1)
INPUT:
# tsize - size of the array;
# arr(0:tsize-1) - array;
# bnorm2 - must be zero on entrance (resides on device as well);
OUTPUT:
# bnorm2 - squared 2-norm of the array (resides on device as well);
**/
{
size_t i,n;
double _thread_norm2;
extern __shared__ double thread_norms2[]; //size = blockDim.x*sizeof(double) Bytes per thread block
n=gridDim.x*blockDim.x; _thread_norm2=0.0;
for(i=blockIdx.x*blockDim.x+threadIdx.x;i<tsize;i+=n) _thread_norm2+=talshComplex4Asq(arr[i]);
thread_norms2[threadIdx.x]=_thread_norm2;
__syncthreads();
if(threadIdx.x == 0){ //global reduction among thread blocks (one thread per block)
_thread_norm2=thread_norms2[0];
for(i=1;i<blockDim.x;i++) _thread_norm2+=thread_norms2[i];
i=1; while(i == 1){i=atomicMax(&norm2_wr_lock,1);} //waiting for the lock to unlock, then lock
__threadfence();
*bnorm2+=_thread_norm2;
__threadfence();
i=atomicExch(&norm2_wr_lock,0); //unlock
}
__syncthreads();
return;
}
// COMPLEX8:
template <>
__global__ void gpu_array_norm2__<talshComplex8>(size_t tsize, const talshComplex8 * __restrict__ arr, volatile double * bnorm2)
/** Computes the squared 2-norm of array arr(0:tsize-1)
INPUT:
# tsize - size of the array;
# arr(0:tsize-1) - array;
# bnorm2 - must be zero on entrance (resides on device as well);
OUTPUT:
# bnorm2 - squared 2-norm of the array (resides on device as well);
**/
{
size_t i,n;
double _thread_norm2;
extern __shared__ double thread_norms2[]; //size = blockDim.x*sizeof(double) Bytes per thread block
n=gridDim.x*blockDim.x; _thread_norm2=0.0;
for(i=blockIdx.x*blockDim.x+threadIdx.x;i<tsize;i+=n) _thread_norm2+=talshComplex8Asq(arr[i]);
thread_norms2[threadIdx.x]=_thread_norm2;
__syncthreads();
if(threadIdx.x == 0){ //global reduction among thread blocks (one thread per block)
_thread_norm2=thread_norms2[0];
for(i=1;i<blockDim.x;i++) _thread_norm2+=thread_norms2[i];
i=1; while(i == 1){i=atomicMax(&norm2_wr_lock,1);} //waiting for the lock to unlock, then lock
__threadfence();
*bnorm2+=_thread_norm2;
__threadfence();
i=atomicExch(&norm2_wr_lock,0); //unlock
}
__syncthreads();
return;
}
//------------------------------------------------------------
// ARRAY INITIALIZATION:
template <typename T>
__global__ void gpu_array_init__(size_t tsize, T * arr, T val)
/** arr(0:tsize-1)=val **/
{
size_t _ti = blockIdx.x*blockDim.x + threadIdx.x;
size_t _gd = gridDim.x*blockDim.x;
for(size_t l = _ti; l < tsize; l += _gd) arr[l] = val;
return;
}
//---------------------------------------------------------------------------------------------------
// SCALAR MULTIPLICATION:
// REAL:
template <typename T>
__global__ void gpu_scalar_multiply__(const T * left_arg, const T * right_arg, T * dest_arg, T alpha,
int left_conj, int right_conj)
/** Scalar += Scalar * Scalar * Alpha **/
{
if(blockIdx.x == 0 && threadIdx.x == 0){
*dest_arg+=(*left_arg)*(*right_arg)*alpha;
}
return;
}
// COMPLEX4:
template <>
__global__ void gpu_scalar_multiply__<talshComplex4>(const talshComplex4 * left_arg, const talshComplex4 * right_arg,
talshComplex4 * dest_arg, talshComplex4 alpha,
int left_conj, int right_conj)
/** Scalar += Scalar * Scalar * Alpha **/
{
if(blockIdx.x == 0 && threadIdx.x == 0){
if(left_conj != 0){
if(right_conj != 0){
*dest_arg=talshComplex4Add(*dest_arg,talshComplex4Mul(talshComplex4Mul(talshComplex4Conjg(*left_arg),talshComplex4Conjg(*right_arg)),alpha));
}else{
*dest_arg=talshComplex4Add(*dest_arg,talshComplex4Mul(talshComplex4Mul(talshComplex4Conjg(*left_arg),*right_arg),alpha));
}
}else{
if(right_conj != 0){
*dest_arg=talshComplex4Add(*dest_arg,talshComplex4Mul(talshComplex4Mul(*left_arg,talshComplex4Conjg(*right_arg)),alpha));
}else{
*dest_arg=talshComplex4Add(*dest_arg,talshComplex4Mul(talshComplex4Mul(*left_arg,*right_arg),alpha));
}
}
}
return;
}
// COMPLEX8:
template <>
__global__ void gpu_scalar_multiply__<talshComplex8>(const talshComplex8 * left_arg, const talshComplex8 * right_arg,
talshComplex8 * dest_arg, talshComplex8 alpha,
int left_conj, int right_conj)
/** Scalar += Scalar * Scalar * Alpha **/
{
if(blockIdx.x == 0 && threadIdx.x == 0){
if(left_conj != 0){
if(right_conj != 0){
*dest_arg=talshComplex8Add(*dest_arg,talshComplex8Mul(talshComplex8Mul(talshComplex8Conjg(*left_arg),talshComplex8Conjg(*right_arg)),alpha));
}else{
*dest_arg=talshComplex8Add(*dest_arg,talshComplex8Mul(talshComplex8Mul(talshComplex8Conjg(*left_arg),*right_arg),alpha));
}
}else{
if(right_conj != 0){
*dest_arg=talshComplex8Add(*dest_arg,talshComplex8Mul(talshComplex8Mul(*left_arg,talshComplex8Conjg(*right_arg)),alpha));
}else{
*dest_arg=talshComplex8Add(*dest_arg,talshComplex8Mul(talshComplex8Mul(*left_arg,*right_arg),alpha));
}
}
}
return;
}
//---------------------------------------------------------------
// ARRAY RESCALING:
// REAL:
template <typename T>
__global__ void gpu_array_scale__(size_t tsize, T * arr, T alpha)
/** arr(0:tsize-1)*=alpha **/
{
size_t _ti = blockIdx.x*blockDim.x + threadIdx.x;
size_t _gd = gridDim.x*blockDim.x;
for(size_t l = _ti; l < tsize; l += _gd) arr[l]*=alpha;
return;
}
// COMPLEX4:
template <>
__global__ void gpu_array_scale__<talshComplex4>(size_t tsize, talshComplex4 * arr, talshComplex4 alpha)
/** arr(0:tsize-1)*=alpha **/
{
size_t _ti = blockIdx.x*blockDim.x + threadIdx.x;
size_t _gd = gridDim.x*blockDim.x;
for(size_t l = _ti; l < tsize; l += _gd) arr[l]=talshComplex4Mul(arr[l],alpha);
return;
}
// COMPLEX8:
template <>
__global__ void gpu_array_scale__<talshComplex8>(size_t tsize, talshComplex8 * arr, talshComplex8 alpha)
/** arr(0:tsize-1)*=alpha **/
{
size_t _ti = blockIdx.x*blockDim.x + threadIdx.x;
size_t _gd = gridDim.x*blockDim.x;
for(size_t l = _ti; l < tsize; l += _gd) arr[l]=talshComplex8Mul(arr[l],alpha);
return;
}
//-----------------------------------------------------------------------------------------------------------------------
// ARRAY ADDITION:
// REAL:
template <typename T>
__global__ void gpu_array_add__(size_t tsize, T * __restrict__ arr0, const T * __restrict__ arr1, T alpha, int left_conj)
/** arr0(0:tsize-1)+=arr1(0:tsize-1)*alpha **/
{
size_t _ti = blockIdx.x*blockDim.x + threadIdx.x;
size_t _gd = gridDim.x*blockDim.x;
for(size_t l = _ti; l < tsize; l += _gd) arr0[l]+=(arr1[l]*alpha);
return;
}
// COMPLEX4:
template <>
__global__ void gpu_array_add__<talshComplex4>(size_t tsize, talshComplex4 * __restrict__ arr0, const talshComplex4 * __restrict__ arr1,
talshComplex4 alpha, int left_conj)
/** arr0(0:tsize-1)+=arr1(0:tsize-1)*alpha **/
{
size_t _ti = blockIdx.x*blockDim.x + threadIdx.x;
size_t _gd = gridDim.x*blockDim.x;
if(left_conj != 0){
for(size_t l = _ti; l < tsize; l += _gd) arr0[l]=talshComplex4Add(arr0[l],talshComplex4Mul(talshComplex4Conjg(arr1[l]),alpha));
}else{
for(size_t l = _ti; l < tsize; l += _gd) arr0[l]=talshComplex4Add(arr0[l],talshComplex4Mul(arr1[l],alpha));
}
return;
}
// COMPLEX8:
template <>
__global__ void gpu_array_add__<talshComplex8>(size_t tsize, talshComplex8 * __restrict__ arr0, const talshComplex8 * __restrict__ arr1,
talshComplex8 alpha, int left_conj)
/** arr0(0:tsize-1)+=arr1(0:tsize-1)*alpha **/
{
size_t _ti = blockIdx.x*blockDim.x + threadIdx.x;
size_t _gd = gridDim.x*blockDim.x;
if(left_conj != 0){
for(size_t l = _ti; l < tsize; l += _gd) arr0[l]=talshComplex8Add(arr0[l],talshComplex8Mul(talshComplex8Conjg(arr1[l]),alpha));
}else{
for(size_t l = _ti; l < tsize; l += _gd) arr0[l]=talshComplex8Add(arr0[l],talshComplex8Mul(arr1[l],alpha));
}
return;
}
//------------------------------------------------------------------------------------------------------------------------------
// ARRAY ADDITION AND SCALING:
// REAL:
template <typename T>
__global__ void gpu_array_add__(size_t tsize, T * __restrict__ arr0, const T * __restrict__ arr1, const T * __restrict__ scalar,
T alpha, int left_conj)
/** arr0(0:tsize-1)+=arr1(0:tsize-1)*scalar*alpha **/
{
size_t _ti = blockIdx.x*blockDim.x + threadIdx.x;
size_t _gd = gridDim.x*blockDim.x;
T pref = (*scalar) * alpha;
for(size_t l = _ti; l < tsize; l += _gd) arr0[l]+=(arr1[l]*pref);
return;
}
// COMPLEX4:
template <>
__global__ void gpu_array_add__<talshComplex4>(size_t tsize, talshComplex4 * __restrict__ arr0, const talshComplex4 * __restrict__ arr1,
const talshComplex4 * __restrict__ scalar, talshComplex4 alpha, int left_conj)
/** arr0(0:tsize-1)+=arr1(0:tsize-1)*scalar*alpha **/
{
size_t _ti = blockIdx.x*blockDim.x + threadIdx.x;
size_t _gd = gridDim.x*blockDim.x;
talshComplex4 pref = talshComplex4Mul(*scalar,alpha);
if(left_conj != 0){
for(size_t l = _ti; l < tsize; l += _gd) arr0[l]=talshComplex4Add(arr0[l],talshComplex4Mul(talshComplex4Conjg(arr1[l]),pref));
}else{
for(size_t l = _ti; l < tsize; l += _gd) arr0[l]=talshComplex4Add(arr0[l],talshComplex4Mul(arr1[l],pref));
}
return;
}
// COMPLEX8:
template <>
__global__ void gpu_array_add__<talshComplex8>(size_t tsize, talshComplex8 * __restrict__ arr0, const talshComplex8 * __restrict__ arr1,
const talshComplex8 * __restrict__ scalar, talshComplex8 alpha, int left_conj)
/** arr0(0:tsize-1)+=arr1(0:tsize-1)*scalar*alpha **/
{
size_t _ti = blockIdx.x*blockDim.x + threadIdx.x;
size_t _gd = gridDim.x*blockDim.x;
talshComplex8 pref = talshComplex8Mul(*scalar,alpha);
if(left_conj != 0){
for(size_t l = _ti; l < tsize; l += _gd) arr0[l]=talshComplex8Add(arr0[l],talshComplex8Mul(talshComplex8Conjg(arr1[l]),pref));
}else{
for(size_t l = _ti; l < tsize; l += _gd) arr0[l]=talshComplex8Add(arr0[l],talshComplex8Mul(arr1[l],pref));
}
return;
}
//-------------------------------------------------------------------------------------------------------
// ARRAY DOT-PRODUCT:
// REAL:
template <typename T>
__global__ void gpu_array_dot_product__(size_t tsize, const T * arr1, const T * arr2, volatile T * dprod,
T alpha, int left_conj, int right_conj)
/** Scalar (GPU) += arr1(0:tsize-1) * arr2(0:tsize-1) * alpha **/
{
extern __shared__ char sh_buf[]; //size = blockDim.x * sizeof(T) Bytes per thread block
T *dprs;
T dpr;
size_t l;
unsigned int j,s;
int i;
dprs=(T*)(&sh_buf[0]); //dynamic shared memory buffer
dpr=static_cast<T>(0.0);
for(l = blockIdx.x*blockDim.x+threadIdx.x; l < tsize; l += gridDim.x*blockDim.x) dpr+=arr1[l]*arr2[l];
dprs[threadIdx.x]=dpr*alpha;
__syncthreads();
s = blockDim.x;
while(s > 1){
j = (s+1U)>>1;
if(threadIdx.x + j < s) dprs[threadIdx.x] += dprs[threadIdx.x+j];
__syncthreads();
s=j;
}
if(threadIdx.x == 0){
i=1; while(i != 0){i=atomicMax(&dot_product_wr_lock,1); if(i == 0) *dprod+=dprs[0];}
__threadfence();
i=atomicExch(&dot_product_wr_lock,0); //unlock
}
__syncthreads();
return;
}
// COMPLEX4:
template <>
__global__ void gpu_array_dot_product__<talshComplex4>(size_t tsize, const talshComplex4 * arr1, const talshComplex4 * arr2,
volatile talshComplex4 * dprod, talshComplex4 alpha,
int left_conj, int right_conj)
/** Scalar (GPU) += arr1(0:tsize-1) * arr2(0:tsize-1) * alpha **/
{
extern __shared__ char sh_buf[]; //size = blockDim.x * sizeof(T) Bytes per thread block
talshComplex4 *dprs;
talshComplex4 dpr;
size_t l;
unsigned int j,s;
int i;
dprs=(talshComplex4*)(&sh_buf[0]); //dynamic shared memory buffer
dpr=talshComplex4Set(0.0f,0.0f);
if(left_conj != 0){
if(right_conj != 0){
for(l = blockIdx.x*blockDim.x+threadIdx.x; l < tsize; l += gridDim.x*blockDim.x){
dpr=talshComplex4Add(dpr,talshComplex4Mul(talshComplex4Conjg(arr1[l]),talshComplex4Conjg(arr2[l])));
}
}else{
for(l = blockIdx.x*blockDim.x+threadIdx.x; l < tsize; l += gridDim.x*blockDim.x){
dpr=talshComplex4Add(dpr,talshComplex4Mul(talshComplex4Conjg(arr1[l]),arr2[l]));
}
}
}else{
if(right_conj != 0){
for(l = blockIdx.x*blockDim.x+threadIdx.x; l < tsize; l += gridDim.x*blockDim.x){
dpr=talshComplex4Add(dpr,talshComplex4Mul(arr1[l],talshComplex4Conjg(arr2[l])));
}
}else{
for(l = blockIdx.x*blockDim.x+threadIdx.x; l < tsize; l += gridDim.x*blockDim.x){
dpr=talshComplex4Add(dpr,talshComplex4Mul(arr1[l],arr2[l]));
}
}
}
dprs[threadIdx.x]=talshComplex4Mul(dpr,alpha);
__syncthreads();
s = blockDim.x;
while(s > 1){
j = (s+1U)>>1;
if(threadIdx.x + j < s) dprs[threadIdx.x] = talshComplex4Add(dprs[threadIdx.x],dprs[threadIdx.x+j]);
__syncthreads();
s=j;
}
if(threadIdx.x == 0){
i=1; while(i != 0){
i=atomicMax(&dot_product_wr_lock,1);
if(i == 0){
dprod->x += talshComplex4Real(dprs[0]);
dprod->y += talshComplex4Imag(dprs[0]);
}
}
__threadfence();
i=atomicExch(&dot_product_wr_lock,0); //unlock
}
__syncthreads();
return;
}
// COMPLEX8:
template <>
__global__ void gpu_array_dot_product__<talshComplex8>(size_t tsize, const talshComplex8 * arr1, const talshComplex8 * arr2,
volatile talshComplex8 * dprod, talshComplex8 alpha,
int left_conj, int right_conj)
/** Scalar (GPU) += arr1(0:tsize-1) * arr2(0:tsize-1) * alpha **/
{
extern __shared__ char sh_buf[]; //size = blockDim.x * sizeof(T) Bytes per thread block
talshComplex8 *dprs;
talshComplex8 dpr;
size_t l;
unsigned int j,s;
int i;
dprs=(talshComplex8*)(&sh_buf[0]); //dynamic shared memory buffer
dpr=talshComplex8Set(0.0,0.0);
if(left_conj != 0){
if(right_conj != 0){
for(l = blockIdx.x*blockDim.x+threadIdx.x; l < tsize; l += gridDim.x*blockDim.x){
dpr=talshComplex8Add(dpr,talshComplex8Mul(talshComplex8Conjg(arr1[l]),talshComplex8Conjg(arr2[l])));
}
}else{
for(l = blockIdx.x*blockDim.x+threadIdx.x; l < tsize; l += gridDim.x*blockDim.x){
dpr=talshComplex8Add(dpr,talshComplex8Mul(talshComplex8Conjg(arr1[l]),arr2[l]));
}
}
}else{
if(right_conj != 0){
for(l = blockIdx.x*blockDim.x+threadIdx.x; l < tsize; l += gridDim.x*blockDim.x){
dpr=talshComplex8Add(dpr,talshComplex8Mul(arr1[l],talshComplex8Conjg(arr2[l])));
}
}else{
for(l = blockIdx.x*blockDim.x+threadIdx.x; l < tsize; l += gridDim.x*blockDim.x){
dpr=talshComplex8Add(dpr,talshComplex8Mul(arr1[l],arr2[l]));
}
}
}
dprs[threadIdx.x]=talshComplex8Mul(dpr,alpha);
__syncthreads();
s = blockDim.x;
while(s > 1){
j = (s+1U)>>1;
if(threadIdx.x + j < s) dprs[threadIdx.x] = talshComplex8Add(dprs[threadIdx.x],dprs[threadIdx.x+j]);
__syncthreads();
s=j;
}
if(threadIdx.x == 0){
i=1; while(i != 0){
i=atomicMax(&dot_product_wr_lock,1);
if(i == 0){
dprod->x += talshComplex8Real(dprs[0]);
dprod->y += talshComplex8Imag(dprs[0]);
}
}
__threadfence();
i=atomicExch(&dot_product_wr_lock,0); //unlock
}
__syncthreads();
return;
}
//---------------------------------------------------------------------------------------------------------
// ARRAY DIRECT PRODUCT:
// REAL:
template <typename T>
__global__ void gpu_array_product__(size_t tsize1, const T * arr1, size_t tsize2, const T * arr2, T * arr0,
T alpha, int left_conj, int right_conj)
/** arr0[0:tsize2-1][0:tsize1-1]+=arr1[0:tsize1-1]*arr2[0:tsize2-1]*alpha **/
{
__shared__ T lbuf[THRDS_ARRAY_PRODUCT+1],rbuf[THRDS_ARRAY_PRODUCT];
size_t _ib,_in,_jb,_jn,_tx,_jc;
_tx=(size_t)threadIdx.x;
for(_jb = blockIdx.y*THRDS_ARRAY_PRODUCT; _jb < tsize2; _jb += gridDim.y*THRDS_ARRAY_PRODUCT){
if(_jb+THRDS_ARRAY_PRODUCT > tsize2){_jn=tsize2-_jb;}else{_jn=THRDS_ARRAY_PRODUCT;}
if(_tx < _jn) rbuf[_tx]=arr2[_jb+_tx]*alpha;
for(_ib = blockIdx.x*THRDS_ARRAY_PRODUCT; _ib < tsize1; _ib += gridDim.x*THRDS_ARRAY_PRODUCT){
if(_ib+THRDS_ARRAY_PRODUCT > tsize1){_in=tsize1-_ib;}else{_in=THRDS_ARRAY_PRODUCT;}
if(_tx < _in) lbuf[_tx]=arr1[_ib+_tx];
__syncthreads();
for(_jc = 0; _jc < _jn; _jc++){if(_tx < _in) arr0[(_jb+_jc)*tsize1+_ib+_tx]+=lbuf[_tx]*rbuf[_jc];}
__syncthreads();
}
}
return;
}
// COMPLEX4:
template <>
__global__ void gpu_array_product__<talshComplex4>(size_t tsize1, const talshComplex4 * arr1,
size_t tsize2, const talshComplex4 * arr2,
talshComplex4 * arr0, talshComplex4 alpha,
int left_conj, int right_conj)
/** arr0[0:tsize2-1][0:tsize1-1]+=arr1[0:tsize1-1]*arr2[0:tsize2-1]*alpha **/
{
__shared__ talshComplex4 lbuf[THRDS_ARRAY_PRODUCT+1],rbuf[THRDS_ARRAY_PRODUCT];
size_t _ib,_in,_jb,_jn,_tx,_jc,_ja;
_tx=(size_t)threadIdx.x;
for(_jb = blockIdx.y*THRDS_ARRAY_PRODUCT; _jb < tsize2; _jb += gridDim.y*THRDS_ARRAY_PRODUCT){
if(_jb+THRDS_ARRAY_PRODUCT > tsize2){_jn=tsize2-_jb;}else{_jn=THRDS_ARRAY_PRODUCT;}
if(right_conj != 0){
if(_tx < _jn) rbuf[_tx]=talshComplex4Mul(talshComplex4Conjg(arr2[_jb+_tx]),alpha);
}else{
if(_tx < _jn) rbuf[_tx]=talshComplex4Mul(arr2[_jb+_tx],alpha);
}
for(_ib = blockIdx.x*THRDS_ARRAY_PRODUCT; _ib < tsize1; _ib += gridDim.x*THRDS_ARRAY_PRODUCT){
if(_ib+THRDS_ARRAY_PRODUCT > tsize1){_in=tsize1-_ib;}else{_in=THRDS_ARRAY_PRODUCT;}
if(left_conj != 0){
if(_tx < _in) lbuf[_tx]=talshComplex4Conjg(arr1[_ib+_tx]);
}else{
if(_tx < _in) lbuf[_tx]=arr1[_ib+_tx];
}
__syncthreads();
for(_jc = 0; _jc < _jn; _jc++){
if(_tx < _in){
_ja = (_jb+_jc)*tsize1 + (_ib+_tx);
arr0[_ja]=talshComplex4Add(arr0[_ja],talshComplex4Mul(lbuf[_tx],rbuf[_jc]));
}
}
__syncthreads();
}
}
return;
}
// COMPLEX8:
template <>
__global__ void gpu_array_product__<talshComplex8>(size_t tsize1, const talshComplex8 * arr1,
size_t tsize2, const talshComplex8 * arr2,
talshComplex8 * arr0, talshComplex8 alpha,
int left_conj, int right_conj)
/** arr0[0:tsize2-1][0:tsize1-1]+=arr1[0:tsize1-1]*arr2[0:tsize2-1]*alpha **/
{
__shared__ talshComplex8 lbuf[THRDS_ARRAY_PRODUCT+1],rbuf[THRDS_ARRAY_PRODUCT];
size_t _ib,_in,_jb,_jn,_tx,_jc,_ja;
_tx=(size_t)threadIdx.x;
for(_jb = blockIdx.y*THRDS_ARRAY_PRODUCT; _jb < tsize2; _jb += gridDim.y*THRDS_ARRAY_PRODUCT){
if(_jb+THRDS_ARRAY_PRODUCT > tsize2){_jn=tsize2-_jb;}else{_jn=THRDS_ARRAY_PRODUCT;}
if(right_conj != 0){
if(_tx < _jn) rbuf[_tx]=talshComplex8Mul(talshComplex8Conjg(arr2[_jb+_tx]),alpha);
}else{
if(_tx < _jn) rbuf[_tx]=talshComplex8Mul(arr2[_jb+_tx],alpha);
}
for(_ib = blockIdx.x*THRDS_ARRAY_PRODUCT; _ib < tsize1; _ib += gridDim.x*THRDS_ARRAY_PRODUCT){
if(_ib+THRDS_ARRAY_PRODUCT > tsize1){_in=tsize1-_ib;}else{_in=THRDS_ARRAY_PRODUCT;}
if(left_conj != 0){
if(_tx < _in) lbuf[_tx]=talshComplex8Conjg(arr1[_ib+_tx]);
}else{
if(_tx < _in) lbuf[_tx]=arr1[_ib+_tx];
}
__syncthreads();
for(_jc = 0; _jc < _jn; _jc++){
if(_tx < _in){
_ja = (_jb+_jc)*tsize1 + (_ib+_tx);
arr0[_ja]=talshComplex8Add(arr0[_ja],talshComplex8Mul(lbuf[_tx],rbuf[_jc]));
}
}
__syncthreads();
}
}
return;
}
//---------------------------------------------------------------------------------------------------
// TENSOR TRANSPOSE-ADD (shared-memory version):
template <typename T>
__global__ void gpu_tensor_block_add_dlf__(int dmo, int drc, int dim_num, int const_args_pos,
const T * __restrict__ tens_in, T * __restrict__ tens_out)
/**
Shared-memory version of tensor transpose-add: tens_out+=TRN(tens_in):
INPUT:
# dmo - dimension extents order (0: normal, as it is in <const_args>; not 0: permuted dimension order will be imposed);
# drc - index permutation direction (0: normal, as it is in <const_args>; not 0: inversed permutation will be used);
# dim_num - tensor block rank;
# const_args_pos - entry in the __constant__ memory bank where tensor block dimension extents (const_args_dims)
and index permutation (const_args_prmn) are stored;
# tens_in[0:] - input tensor;
OUTPUT:
# tens_out[0:] - output (transposed) tensor in which accumulation is performed;
NOTES:
# Minimal CUDA execution configuration is <<<1,warpSize>>>
# Number of threads per block must be multiple of the warpSize!
**/
{
__shared__ T buf0[TENS_TRANSP_BUF_SIZE];
__shared__ float val;
__shared__ size_t base_in[MAX_TENSOR_RANK],base_out[MAX_TENSOR_RANK];
__shared__ size_t ftb[TENS_TRANSP_TAB_SIZE],gtb[TENS_TRANSP_TAB_SIZE];
__shared__ int htb[TENS_TRANSP_TAB_SIZE],stb[TENS_TRANSP_TAB_SIZE];
__shared__ int dim_in[MAX_TENSOR_RANK],dim_out[MAX_TENSOR_RANK],o2n[MAX_TENSOR_RANK],n2o[MAX_TENSOR_RANK];
__shared__ int pri[MAX_TENSOR_RANK],tmp0[MAX_TENSOR_RANK];
__shared__ int err_code,minor,minor_in,minor_out,s1_ind,s1_ond,s1_step,s1_dim,s2_ind,s2_ond,s2_step,s2_dim,ns1,ns2;
__shared__ size_t vol,vol_ext;
size_t _vol,_addr_in,_addr_out,_addr,_work_piece;
int i,j,k,l,m,n,_vol_minor,_vol_in,_vol_out,_s1,_s2;
/*
SHARED MEMORY USE (bytes) =
+ TENS_TRANSP_BUF_SIZE*sizeof(T)
+ MAX_TENSOR_RANK*(8+8+4+4+4+4+4+4)
+ TENS_TRANSP_TAB_SIZE*(8+8+4+4)
+ 4*15 + 8*2
MIN REGISTER USE (bytes) per thread =
+ 4*4 + 4*11 + 8*5 = 100
*/
//Determine the minor index set (only the master thread in each thread block):
if(threadIdx.x == 0){
err_code=0;
if(dim_num >= 0 && dim_num <= MAX_TENSOR_RANK && blockDim.x >= warpSize && blockDim.x%warpSize == 0){
s1_ind=dim_num+1; s2_ind=dim_num-1;
_vol=1; for(i=0;i<dim_num;i++){
_vol*=const_args_dims[const_args_pos][i]; if(const_args_prmn[const_args_pos][i] != i+1) s1_ind=0;
}; vol=_vol; //total volume (number of tensor elements)
if(s1_ind == 0){ //non-trivial permutation
// Set input/output permutations and dimension extents:
if(drc == 0){ //normal index permutation
for(i=0;i<dim_num;i++) o2n[i]=const_args_prmn[const_args_pos][i]-1; for(i=0;i<dim_num;i++) n2o[o2n[i]]=i;
}else{ //inversed index permutation
for(i=0;i<dim_num;i++) n2o[i]=const_args_prmn[const_args_pos][i]-1; for(i=0;i<dim_num;i++) o2n[n2o[i]]=i;
}
if(dmo == 0){ //normal dimension order
for(i=0;i<dim_num;i++) dim_in[i]=const_args_dims[const_args_pos][i];
for(i=0;i<dim_num;i++) dim_out[o2n[i]]=dim_in[i];
}else{ //inversed dimension order
for(i=0;i<dim_num;i++) dim_out[i]=const_args_dims[const_args_pos][i];
for(i=0;i<dim_num;i++) dim_in[n2o[i]]=dim_out[i];
}
s1_step=dim_in[s1_ind]; s2_step=dim_in[s2_ind];
if(_vol > TENS_TRANSP_BUF_SIZE){ //tensor block does not fit into the shared memory buffer
// Determine the input/output minor index sets and the combined minor index set:
l=(int)(sqrt((float)TENS_TRANSP_BUF_SIZE));
minor_in=0; _vol_in=1; for(i=0;i<dim_num;i++){j=_vol_in*dim_in[i]; if(j>l) break; minor_in++; _vol_in=j;}
minor_out=0; _vol_out=1; for(i=0;i<dim_num;i++){j=_vol_out*dim_out[i]; if(j>l) break; minor_out++; _vol_out=j;}
minor=minor_in; _vol_minor=_vol_in; for(i=0;i<minor_out;i++){if(n2o[i]>=minor_in){minor++; _vol_minor*=dim_out[i];}}
m=1; _s1=0; _s2=0;
while(_vol_minor < TENS_TRANSP_BUF_SIZE && m != 0){
m=0;
if(_s1 == 0){for(i=minor_in;i<dim_num;i++){if(o2n[i]<minor_out){minor_in++; _vol_in*=dim_in[i];}else{break;}}}
if(_s2 == 0){for(i=minor_out;i<dim_num;i++){if(n2o[i]<minor_in){minor_out++; _vol_out*=dim_out[i];}else{break;}}}
j=dim_in[minor_in]; l=dim_out[minor_out];
if(minor_in == n2o[minor_out] && _s1+_s2 == 0){ //same candidate index to both the input and output index sets
if(j > 1 && TENS_TRANSP_BUF_SIZE < _vol_minor*2) break;
if(_vol_minor*j > TENS_TRANSP_BUF_SIZE){s1_ind=minor_in; s1_step=TENS_TRANSP_BUF_SIZE/_vol_minor; _s1++; _s2++;}
minor_in++; _vol_in*=j; minor_out++; _vol_out*=j; minor++; _vol_minor*=j; m++;
}else{ //the input and output index sets consider two different candidates
if(_vol_minor*j*l <= TENS_TRANSP_BUF_SIZE && _s1+_s2 == 0){ //accept both, no splitting
minor_in++; _vol_in*=j; minor_out++; _vol_out*=l; minor+=2; _vol_minor*=(j*l); m++;
}else{ //try to accept either one of the two OR both with splitting
if(j == 1 || l == 1){
if(j == 1 && _s1 == 0){minor_in++; minor++; m++;}
if(l == 1 && _s2 == 0){minor_out++; minor++; m++;}
}else{
if(_vol_minor*j <= TENS_TRANSP_BUF_SIZE && _vol_minor*l > TENS_TRANSP_BUF_SIZE &&
_vol_out >= warpSize && _s1 == 0){ //accept the input index, no splitting
minor_in++; _vol_in*=j; minor++; _vol_minor*=j; m++;
}else if(_vol_minor*j > TENS_TRANSP_BUF_SIZE && _vol_minor*l <= TENS_TRANSP_BUF_SIZE &&
_vol_in >= warpSize && _s2 == 0){ //accept the output index, no splitting
minor_out++; _vol_out*=l; minor++; _vol_minor*=l; m++;
}else{ //splitting is unavoidable (both OR one OR none)
if(TENS_TRANSP_BUF_SIZE >= _vol_minor*2){
if(j >= 4 && l >= 4){ //dimension extents are large enough to be split
if(_vol_minor*4 > TENS_TRANSP_BUF_SIZE){ //impossible to split both indices
if(_vol_in <= _vol_out && _s1 == 0){ //split the input candidate index
s1_ind=minor_in; s1_step=TENS_TRANSP_BUF_SIZE/_vol_minor;
minor_in++; _vol_in*=j; minor++; _vol_minor*=j; _s1++; m++;
}else{ //split the output candidate index
if(_s2 == 0){
s1_ind=n2o[minor_out]; s1_step=TENS_TRANSP_BUF_SIZE/_vol_minor;
minor_out++; _vol_out*=l; minor++; _vol_minor*=l; _s2++; m++;
}
}
}else{ //possible to split both indices
i=(int)sqrt(((float)TENS_TRANSP_BUF_SIZE)/(float)_vol_minor); if(i < 2) i=2; //uniform splitting
s1_step=i; s2_step=i; val=(float)_vol_out/(float)_vol_in;
if(val < 1.0f){ //scale the initial uniform splitting to reflect the disbalance between _vol_in and _vol_out
if(val*(float)i < 1.0f) val=1.0f/(float)i; if(val*(float)l < (float)i) val=(float)i/(float)l;
}else{
if(val*(float)i > (float)j) val=(float)j/(float)i; if(val > float(i)) val=(float)i;
}
s1_step=(int)(((float)i)*val); s2_step=(int)(((float)i)/val);
if(s1_step >= 2 && _s1 == 0){ //&& s1_step <= dim_in[minor_in]
s1_ind=minor_in; minor_in++; _vol_in*=j; minor++; _vol_minor*=j; _s1++; m++;
}else{
s1_step=dim_in[s1_ind];
}
if(s2_step >= 2 && _s2 == 0){ //&& s2_step <= dim_out[minor_out]
s2_ind=n2o[minor_out]; minor_out++; _vol_out*=l; minor++; _vol_minor*=l; _s2++; m++;
}else{
s2_step=dim_in[s2_ind];
}
}
}else if(j >= 4 && l < 4 && _s1 == 0){ //split the input candidate index
s1_ind=minor_in; s1_step=TENS_TRANSP_BUF_SIZE/_vol_minor;
minor_in++; _vol_in*=j; minor++; _vol_minor*=j; _s1++; m++;
}else if(j < 4 && l >= 4 && _s2 == 0){ //split the output candidate index
s1_ind=n2o[minor_out]; s1_step=TENS_TRANSP_BUF_SIZE/_vol_minor;
minor_out++; _vol_out*=l; minor++; _vol_minor*=l; _s2++; m++;
}else{ //both candidate indices have too small extent to be split: try to add one of them fully
if(_vol_minor*j <= TENS_TRANSP_BUF_SIZE && _s1 == 0){
minor_in++; _vol_in*=j; minor++; _vol_minor*=j; m++;
}else if(_vol_minor*l <= TENS_TRANSP_BUF_SIZE && _s2 == 0){
minor_out++; _vol_out*=l; minor++; _vol_minor*=l; m++;
}
}
}else{ //unable to add more indices in the minor set
break;
}
}
}
}
}
}
if(s1_ind == dim_num-1 && s2_ind == dim_num-1){s2_ind=0; s2_step=dim_in[0];} //s1_ind was set while s2_ind was not
}else{ //tensor block fits into the shared memory buffer from the beginning
minor=dim_num; minor_in=dim_num; minor_out=dim_num; _vol_minor=_vol; _vol_in=_vol; _vol_out=_vol;
}
// Share the tensor transpose configuration with other threads in each block:
vol_ext=_vol/_vol_minor; s1_dim=dim_in[s1_ind]; s2_dim=dim_in[s2_ind];
// Set indexing bases (OUT:{out,in_c,ext_in}_new; IN:{in,out_c,ext_in}_old):
// OUTPUT indexing (dim_out[], base_out[]: prioritized new numeration):
for(i=0;i<dim_num;i++){tmp0[i]=dim_out[i];} //save output dimension extents (new numeration)
j=0; for(i=0;i<minor_out;i++){pri[j++]=i;} //output minor index set (new numeration))
for(i=0;i<dim_num;i++){if(o2n[i]>=minor_out) pri[j++]=o2n[i];} //{compl.input minor + external} index set (new numeration)
j=1; for(i=0;i<dim_num;i++){dim_out[i]=j; j*=tmp0[i];} //output bases (new numeration)
for(i=0;i<dim_num;i++){base_out[i]=dim_out[pri[i]];} //output bases (prioritized new numeration)
for(i=0;i<dim_num;i++){dim_out[i]=tmp0[pri[i]];} //output extents (prioritized new numeration)
for(i=0;i<dim_num;i++){if(n2o[pri[i]]==s1_ind){s1_ond=i;}else if(n2o[pri[i]]==s2_ind){s2_ond=i;}} //split indices (prioritized new numeration)
// INPUT indexing (dim_in[], base_in[]: prioritized old numeration):
for(i=0;i<dim_num;i++){tmp0[i]=dim_in[i];} //save input dimension extents (old numeration)
j=0; for(i=0;i<minor_in;i++){pri[j++]=i;} //input minor index set (old numeration)
for(i=0;i<minor_out;i++){if(n2o[i]>=minor_in) pri[j++]=n2o[i];} //compl.output minor idex set (old numeration)
for(i=j;i<dim_num;i++){pri[i]=n2o[pri[i]];} //external index set (just convert new numbers to old ones for consistency)
j=1; for(i=0;i<dim_num;i++){dim_in[i]=j; j*=tmp0[i];} //input bases (old numeration)
for(i=0;i<dim_num;i++){base_in[i]=dim_in[pri[i]];} //input bases (prioritized old numeration)
for(i=0;i<dim_num;i++){dim_in[i]=tmp0[pri[i]];} //input extents (prioritized old numeration)
for(i=0;i<dim_num;i++){if(pri[i]==s1_ind){_s1=i;}else if(pri[i]==s2_ind){_s2=i;}} //split indices (prioritized old numeration)
s1_ind=_s1; s2_ind=_s2;
ns1=1+(s1_dim-1)/s1_step; //number of segments from the 1st split minor index
ns2=1+(s2_dim-1)/s2_step; //number of segments from the 2nd split minor index
// Index position correspondence for the minor index set (pri-new --> pri-old):
j=0; for(i=0;i<minor_out;i++){if(n2o[i]<minor_in){pri[i]=n2o[i];}else{pri[i]=(minor_in+j); j++;}}
j=0; for(i=0;i<minor_in;i++){if(o2n[i]<minor_out){pri[o2n[i]]=i;}else{pri[minor_out+j]=i; j++;}}
// Check tensor transpose configuration parameters:
if(minor <= 0 || minor_in <= 0 || minor_out <= 0 || _vol <= 0 || _vol_minor <= 0) err_code+=5000; //trap
if(s1_ind >= dim_num || s2_ind >= dim_num || s1_ond >= dim_num || s2_ond >= dim_num ||
s1_ind == s2_ind || s1_ond == s2_ond || s1_step <= 0 || s2_step <= 0) err_code+=1000; //trap
if((s1_step != dim_in[s1_ind] && s1_ind != minor_in-1 && s1_ond != minor_out-1) ||
(s2_step != dim_in[s2_ind] && s2_ind != minor_in-1 && s2_ond != minor_out-1)) err_code+=500; //trap
if((_vol_minor*s1_step*s2_step)/(s1_dim*s2_dim) > TENS_TRANSP_BUF_SIZE) err_code+=100; //trap
} //endif: non-trivial permutation
}else{
err_code=1+2*blockDim.x%warpSize;
}
} //endif: Master thread.
#ifdef DEBUG_GPU
//DEBUG RECORD begin:
if(blockIdx.x == 0 && threadIdx.x == 0){
j=0; gpu_debug_dump[j++]=dim_num;
for(i=0;i<dim_num;i++) gpu_debug_dump[j++]=const_args_dims[const_args_pos][i];
for(i=0;i<dim_num;i++) gpu_debug_dump[j++]=const_args_prmn[const_args_pos][i];
for(i=0;i<dim_num;i++) gpu_debug_dump[j++]=base_in[i];
for(i=0;i<dim_num;i++) gpu_debug_dump[j++]=base_out[i];
gpu_debug_dump[j++]=vol; gpu_debug_dump[j++]=vol_ext; gpu_debug_dump[j++]=vol/vol_ext;
gpu_debug_dump[j++]=minor; gpu_debug_dump[j++]=minor_in; gpu_debug_dump[j++]=minor_out;
gpu_debug_dump[j++]=s1_ind; gpu_debug_dump[j++]=s1_ond; gpu_debug_dump[j++]=s1_step; gpu_debug_dump[j++]=s1_dim;
gpu_debug_dump[j++]=s2_ind; gpu_debug_dump[j++]=s2_ond; gpu_debug_dump[j++]=s2_step; gpu_debug_dump[j++]=s2_dim;
for(i=0;i<dim_num;i++) gpu_debug_dump[j++]=pri[i];
gpu_debug_dump[j++]=err_code; gpu_debug_dump[j++]=-1;
}
//DEBUG RECORD end.
#endif /*DEBUG_GPU*/
__syncthreads();
//Proceed:
if(err_code == 0){
if(s1_ind > dim_num){ //tag of a trivial permutation
// Direct copy:
_vol=vol; j=gridDim.x*blockDim.x; i=blockIdx.x*blockDim.x+threadIdx.x; _addr_in=_vol-_vol%j;
for(_addr=0;_addr<_addr_in;_addr+=j){
_addr_out=_addr+i;
tens_out[_addr_out]+=tens_in[_addr_out];