-
Notifications
You must be signed in to change notification settings - Fork 58
/
Copy pathfft.c
1554 lines (1460 loc) · 63.3 KB
/
fft.c
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
/* Initialization of all FFT for matrix-vector products; and FFT procedures themselves; not used in sparse mode
* TODO: A lot of indirect indexing used - way to optimize.
*
* Copyright (C) ADDA contributors
* This file is part of ADDA.
*
* ADDA is free software: you can redistribute it and/or modify it under the terms of the GNU General Public License as
* published by the Free Software Foundation, either version 3 of the License, or (at your option) any later version.
*
* ADDA 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 General Public License for more details.
*
* You should have received a copy of the GNU General Public License along with ADDA. If not, see
* <http://www.gnu.org/licenses/>.
*/
/* The following tests for compilation inconsistencies, but also helps proper syntax checking in IDE, such as Eclipse.
* Otherwise, a lot of unresolved-symbol errors are produced, when another build configuration is selected.
*/
#ifdef SPARSE
# error "This file is incompatible with SPARSE"
# undef SPARSE
#endif
#include "const.h" // keep this first
#include "fft.h" // corresponding header
// project headers
#include "cmplx.h"
#include "comm.h"
#include "debug.h"
#include "function.h"
#include "io.h"
#include "interaction.h"
#include "memory.h"
#include "oclcore.h"
#include "prec_time.h"
#include "vars.h"
// system headers
#include <math.h>
#include <stdlib.h>
#include <string.h>
#ifdef CLFFT
IGNORE_WARNING(-Wstrict-prototypes) // no way to change the library header (tested with v.2.12.2)
# include <clFFT.h> //external library
STOP_IGNORE
// Defines precision of clFFT transforms.
// CLFFT_DOUBLE_FAST should be tested when becomes operational - in v.2.12.2 it is identical to CLFFT_DOUBLE
# define PRECISION_CLFFT CLFFT_DOUBLE
#elif defined(CLFFT_APPLE)
# include "cpp/clFFT.h" //nearly unmodified APPLE FFT header file
# ifdef NO_CPP
# error "Apple clFFT relies on C++ sources, hence is incompatible with NO_CPP option"
# endif
#endif
/* standard FFT routines (FFTW3 of FFT_TEMPERTON) are required even when OpenCL is used, since they are used for Fourier
* transform of the D-matrix
*/
#ifdef FFTW3
# include <fftw3.h> // types.h or cmplx.h should be defined before (to match C99 complex type)
/* define level of planning for usual and Dmatrix (DM) FFT: FFTW_ESTIMATE (heuristics), FFTW_MEASURE (default),
* FFTW_PATIENT, or FFTW_EXHAUSTIVE
*/
# define PLAN_FFTW FFTW_MEASURE
# define PLAN_FFTW_DM FFTW_ESTIMATE
# define ONLY_FOR_FFTW3 // this is used in function argument declarations
#else
# define ONLY_FOR_FFTW3 ATT_UNUSED
#endif
#ifdef FFT_TEMPERTON
# define ONLY_FOR_TEMPERTON // this is used in function argument declarations
#else
# define ONLY_FOR_TEMPERTON ATT_UNUSED
#endif
// SEMI-GLOBAL VARIABLES
// defined and initialized in interaction.c
extern const int local_Nz_Rm;
// defined and initialized in timing.c
extern TIME_TYPE Timing_FFT_Init,Timing_Dm_Init;
// used in comm.c
double * restrict BT_buffer, * restrict BT_rbuffer; // buffers for BlockTranspose
// used in matvec.c; in OpenCL mode some of those are not used at all, others - only locally
doublecomplex * restrict Dmatrix; // holds FFT of the interaction matrix
doublecomplex * restrict Rmatrix; // holds FFT of the reflection matrix
#ifndef OPENCL
// holds input vector (on expanded grid) to matvec, also used as storage space in iterative.c
doublecomplex * restrict Xmatrix;
doublecomplex * restrict slices; // used in inner cycle of matvec - holds 3 components (for fixed x)
doublecomplex * restrict slices_tr; // additional storage space for slices to accelerate transpose
doublecomplex * restrict slicesR,* restrict slicesR_tr; // same as above, but for reflected interaction
#endif
size_t DsizeY,DsizeZ,DsizeYZ; // size of the 'matrix' D
size_t RsizeY; // size of the 'matrix' R; in OpenCL mode it is used in oclmatvec.c
// used in oclmatvec.c
#ifdef OPENCL
/* clxslices is the number of slices required in MatVec,
* so that more GPU memory is used and kernels can run over larger arrays
*/
size_t clxslices;
size_t local_gridX; // 'thickness' in x direction of a slice*
size_t slicesize; // total number of doublecomplex values per slice
#endif
// LOCAL VARIABLES
// D2 matrix and its two slices; used only temporary for InitDmatrix
static doublecomplex * restrict slice,* restrict slice_tr,* restrict D2matrix;
static doublecomplex * restrict R2matrix; // same for surface (slice and slice_tr are reused from Dmatrix)
static size_t D2sizeY; // size of the 'matrix' D2 (x-size is gridX), Z size is not used
static size_t R2sizeY; // size of the 'matrix' R2 (x- and z-sizes are corresponding grids)
static size_t lz_Dm,lz_Rm; // local sizes along z for D(2) and R(2) matrices
// the following two lines are defined in InitDmatrix but used in InitRmatrix, they are analogous to Dm values
static size_t Rsize,R2sizeTot; // sizes of R and R2 matrices
static int jstartR; // starting index for y
static bool weird_nprocs; // whether weird number of processors is used
#ifdef OPENCL
// clFFT plans
# ifdef CLFFT
static clfftPlanHandle clplanX,clplanY,clplanZ;
static size_t clfftBufSize=0;
# elif defined(CLFFT_APPLE)
static clFFT_Plan clplanX,clplanY,clplanZ;
# endif
#endif
#ifdef FFTW3
// FFTW3 plans: f - FFT_FORWARD; b - FFT_BACKWARD
static fftw_plan planXf_Dm,planYf_slice,planZf_slice,planXf_Rm;
# ifndef OPENCL // these plans are used only if OpenCL is not used
static fftw_plan planXf,planXb,planYf,planYb,planZf,planZb,planYRf,planZRf; // last two for reflected interaction
# endif
#elif defined(FFT_TEMPERTON)
# ifdef NO_FORTRAN
# error "Tempertron FFT is implemented in Fortran, hence incompatible with NO_FORTRAN option"
# endif
# define IFAX_SIZE 20
// arrays for Temperton FFT
static double * restrict trigsX,* restrict trigsY,* restrict trigsZ,* restrict work;
static int ifaxX[IFAX_SIZE],ifaxY[IFAX_SIZE],ifaxZ[IFAX_SIZE];
// Fortran routines from cfft99D.f
void cftfax_(const int *nn,int * restrict ifax,double * restrict trigs);
void cfft99_(double * restrict data,double * restrict _work,const double * restrict trigs,const int * restrict ifax,
const int *inc,const int *jump,const int *nn,const int *lot,const int *isign);
#endif
#ifdef CLFFT
// Test clFFT version
# define CLFFT_VER_REQ 2
# define CLFFT_SUBVER_REQ 12
# if !GREATER_EQ2(clfftVersionMajor,clfftVersionMinor,CLFFT_VER_REQ,CLFFT_SUBVER_REQ)
# error "clFFT version is too old"
# endif
// Error-checking functionality for clFFT
# define CLFFT_CH_ERR(a) Check_clFFT_Err(a,ALL_POS)
//======================================================================================================================
static const char *Print_clFFT_Errstring(clfftStatus err)
// produces meaningful error message from the clFFT-specific error code, based on clFFT.h v.2.12.2 (NULL if not found)
{
switch (err) {
case CLFFT_BUGCHECK: return "Bugcheck";
case CLFFT_NOTIMPLEMENTED: return "Functionality not implemented";
case CLFFT_TRANSPOSED_NOTIMPLEMENTED: return "Transposed functionality not implemented";
case CLFFT_FILE_NOT_FOUND: return "File not found";
case CLFFT_FILE_CREATE_FAILURE: return "File creation failure";
case CLFFT_VERSION_MISMATCH: return "Version conflict between client and library";
case CLFFT_INVALID_PLAN: return "Invalid plan";
case CLFFT_DEVICE_NO_DOUBLE: return "Device does not support double precision";
case CLFFT_DEVICE_MISMATCH: return "Mismatch of plan and device";
default: return NULL;
}
}
//======================================================================================================================
static void Check_clFFT_Err(const clfftStatus err,ERR_LOC_DECL)
/* Checks error code for clFFT calls and prints error if necessary. First searches among clFFT specific errors. If not
* found, uses general error processing for CL calls (since clFFT error codes can take standard cl values as well).
*/
{
if (err != CLFFT_SUCCESS) {
const char *str=Print_clFFT_Errstring(err);
if (str!=NULL) LogError(ERR_LOC_CALL,"clFFT error code %d: %s\n",err,str);
else PrintCLErr((cl_int)err,ERR_LOC_CALL,NULL);
}
}
#endif
//======================================================================================================================
static inline size_t IndexDmatrix(const size_t x,size_t y,size_t z)
// index D matrix to store final result (symmetric with respect to center for y and z)
{
if (y>=DsizeY) y=gridY-y;
if (z>=DsizeZ) z=gridZ-z;
return(NDCOMP*((x*DsizeZ+z)*DsizeY+y));
}
//======================================================================================================================
static inline size_t IndexGarbledD(const size_t x,int y,int z)
// index D2 matrix after BlockTranspose (periodic over y and z)
{
if (y<0) y+=gridY;
if (z<0) z+=gridZ;
#ifdef PARALLEL
return(((z%lz_Dm)*D2sizeY+y)*gridX+(z/lz_Dm)*local_Nx+x%local_Nx);
#else
return((z*D2sizeY+y)*gridX+x);
#endif
}
//======================================================================================================================
static inline size_t IndexSliceD2matrix(int y,int z)
// index slice of D2 matrix (periodic over y and z)
{
if (y<0) y+=gridY;
if (z<0) z+=gridZ;
return(y*gridZ+z);
}
//======================================================================================================================
static inline size_t Index2matrix(int x,int y,const int z,const int sizeY)
// index D2 or R2 matrix to store calculated elements (periodic over x and y), z should already be shifted
{
if (x<0) x+=gridX;
if (y<0) y+=gridY;
return((z*sizeY+y)*gridX+x);
}
//======================================================================================================================
static inline size_t IndexSlice_zy(const size_t y,const size_t z)
// index transposed slice of D2 (or R2) matrix
{
return (z*gridY+y);
}
//======================================================================================================================
static inline size_t IndexRmatrix(const size_t x,size_t y,const size_t z)
// index R matrix to store final result (symmetric with respect to center for y)
{
if (y>=RsizeY) y=gridY-y;
return(NDCOMP*((x*gridZ+z)*RsizeY+y));
}
//======================================================================================================================
static inline size_t IndexGarbledR(const size_t x,int y,const int z)
// index R2 matrix after BlockTranspose (periodic over y)
{
if (y<0) y+=gridY;
#ifdef PARALLEL
return(((z%lz_Rm)*R2sizeY+y)*gridX+(z/lz_Rm)*local_Nx+x%local_Nx);
#else
return((z*R2sizeY+y)*gridX+x);
#endif
}
//======================================================================================================================
static inline size_t IndexSliceR2matrix(int y,const int z)
// index slice of R2 matrix (periodic over y)
{
if (y<0) y+=gridY;
return(y*gridZ+z);
}
//======================================================================================================================
static void transpose(const doublecomplex * restrict data,doublecomplex * restrict trans,const size_t Y,const size_t Z)
// optimized routine to transpose complex matrix with dimensions YxZ: data -> trans
{
size_t y,z,y1,y2,z1,z2,i,j,y0,z0;
doublecomplex *t1,*t2,*t3,*t4;
const doublecomplex *w1,*w2,*w3;
const size_t blockTr=64; // block size
/* Intel compiler 11.1 seems to produce broken code for this function whenever blockTr>1, at least when the function
* is called in row of three from TransposeYZ(). So maybe the bug is due to incorrect inlining. This bug appears
* only for -O3 compilation, but not for -O2. Moreover, the bug is not present for icc 13.0.
* We have no idea what can be done on our side.
*/
y1=Y/blockTr;
y2=Y%blockTr;
z1=Z/blockTr;
z2=Z%blockTr;
w1=data;
t1=trans-Y;
for(i=0;i<=y1;i++) {
if (i==y1) y0=y2;
else y0=blockTr;
w2=w1;
t2=t1;
for(j=0;j<=z1;j++) {
if (j==z1) z0=z2;
else z0=blockTr;
w3=w2;
t3=t2;
for (y=0;y<y0;y++) {
t4=t3+y;
for (z=0;z<z0;z++) *(t4+=Y)=w3[z];
w3+=Z;
}
w2+=blockTr;
t2+=blockTr*Y;
}
w1+=blockTr*Z;
t1+=blockTr;
}
}
//======================================================================================================================
void TransposeYZ(const int direction)
/* optimized routine to transpose y and z; forward: slices->slices_tr; backward: slices_tr->slices; direction can be
* made boolean but this contradicts with existing definitions of FFT_FORWARD and FFT_BACKWARD, which themselves are
* determined by FFT routines invocation format
*/
{
#ifdef OPENCL
const size_t blocksize=16; //this corresponds to BLOCK_DIM in oclkernels.cl
const size_t tblock[3]={blocksize,blocksize,1};
size_t enqtglobalzy[3]={gridZ,gridY,3*local_gridX};
size_t enqtglobalyz[3]={gridY,gridZ,3*local_gridX};
// if the grid is not divisible by blocksize, extend it. Kernel takes care of borders
size_t tgridZ = (gridZ%blocksize==0) ? gridZ : (gridZ/blocksize+1)*blocksize;
size_t tgridY = (gridY%blocksize==0) ? gridY : (gridY/blocksize+1)*blocksize;
enqtglobalzy[0]=tgridZ;
enqtglobalzy[1]=tgridY;
enqtglobalyz[0]=tgridY;
enqtglobalyz[1]=tgridZ;
if (direction==FFT_FORWARD) {
CL_CH_ERR(clEnqueueNDRangeKernel(command_queue,cltransposeof,3,NULL,enqtglobalzy,tblock,0,NULL,NULL));
if (surface)
CL_CH_ERR(clEnqueueNDRangeKernel(command_queue,cltransposeofR,3,NULL,enqtglobalzy,tblock,0,NULL,NULL));
}
else CL_CH_ERR(clEnqueueNDRangeKernel(command_queue,cltransposeob,3,NULL,enqtglobalyz,tblock,0,NULL,NULL));
#else
size_t Xcomp,ind;
if (direction==FFT_FORWARD) for (Xcomp=0;Xcomp<3;Xcomp++) {
ind=Xcomp*gridYZ;
transpose(slices+ind,slices_tr+ind,gridY,gridZ);
if (surface) transpose(slicesR+ind,slicesR_tr+ind,gridY,gridZ);
}
else for (Xcomp=0;Xcomp<3;Xcomp++) { // direction==FFT_BACKWARD
ind=Xcomp*gridYZ;
transpose(slices_tr+ind,slices+ind,gridZ,gridY);
}
#endif
}
//======================================================================================================================
void fftX(const int isign)
// FFT three components of (buf)Xmatrix(x) for all y,z; called from matvec
{
#ifdef OPENCL
# ifdef CLFFT
CLFFT_CH_ERR(clfftEnqueueTransform(clplanX,(clfftDirection)isign,1,&command_queue,0,NULL,NULL,&bufXmatrix,NULL,
NULL));
# elif defined(CLFFT_APPLE)
CL_CH_ERR(clFFT_ExecuteInterleaved(command_queue,clplanX,(int)3*local_Nz*smallY,(clFFT_Direction)isign,bufXmatrix,
bufXmatrix,0,NULL,NULL));
# endif
#elif defined(FFTW3)
if (isign==FFT_FORWARD) fftw_execute(planXf);
else fftw_execute(planXb);
#elif defined(FFT_TEMPERTON)
int nn=gridX,inc=1,jump=nn,lot=boxY;
size_t z;
/* Calls to Temperton FFT cause warnings for translation from doublecomplex to double pointers. However, such a cast
* is perfectly valid in C99. So we set pragmas to remove these warnings.
*
* !!! TODO: Another (ultimate) solution is to remove this routine altogether, since FFTW is perfect in all
* respects. This is also reasonable considering future switch to tgmath.h
*/
IGNORE_WARNING(-Wstrict-aliasing);
for (z=0;z<3*local_Nz;z++) cfft99_((double *)(Xmatrix+z*gridX*smallY),work,trigsX,ifaxX,&inc,&jump,&nn,&lot,&isign);
STOP_IGNORE;
#endif
}
//======================================================================================================================
void fftY(const int isign)
// FFT three components of slices_tr(y) for all z; called from matvec
{
#ifdef OPENCL
# ifdef CLFFT
CLFFT_CH_ERR(clfftEnqueueTransform(clplanY,(clfftDirection)isign,1,&command_queue,0,NULL,NULL,&bufslices_tr,NULL,
NULL));
if (surface && isign==FFT_FORWARD)
CLFFT_CH_ERR(clfftEnqueueTransform(clplanY,(clfftDirection)isign,1,&command_queue,0,NULL,NULL,&bufslicesR_tr,
NULL,NULL));
# elif defined(CLFFT_APPLE)
CL_CH_ERR(clFFT_ExecuteInterleaved(command_queue,clplanY,(int)3*gridZ*local_gridX,(clFFT_Direction)isign,
bufslices_tr,bufslices_tr,0,NULL,NULL));
if (surface && isign==FFT_FORWARD)
CL_CH_ERR(clFFT_ExecuteInterleaved(command_queue,clplanY,(int)3*gridZ*local_gridX,(clFFT_Direction)isign,
bufslicesR_tr,bufslicesR_tr,0,NULL,NULL));
# endif
#elif defined(FFTW3)
if (isign==FFT_FORWARD) {
fftw_execute(planYf);
if (surface) fftw_execute(planYRf);
}
else fftw_execute(planYb);
#elif defined(FFT_TEMPERTON)
int nn=gridY,inc=1,jump=nn,lot=3*gridZ;
IGNORE_WARNING(-Wstrict-aliasing);
cfft99_((double *)(slices_tr),work,trigsY,ifaxY,&inc,&jump,&nn,&lot,&isign);
// the same operation is applied to sliceR_tr, when required
if (surface && isign==FFT_FORWARD) cfft99_((double *)(slicesR_tr),work,trigsY,ifaxY,&inc,&jump,&nn,&lot,&isign);
STOP_IGNORE;
#endif
}
//======================================================================================================================
void fftZ(const int isign)
// FFT three components of slices(z) for all y; called from matvec
{
#ifdef OPENCL
# ifdef CLFFT
CLFFT_CH_ERR(clfftEnqueueTransform(clplanZ,(clfftDirection)isign,1,&command_queue,0,NULL,NULL,&bufslices,NULL,
NULL));
if (surface && isign==FFT_FORWARD) // the same operation is applied to bufslicesR, but with inverse transform
CLFFT_CH_ERR(clfftEnqueueTransform(clplanZ,(clfftDirection)FFT_BACKWARD,1,&command_queue,0,NULL,NULL,
&bufslicesR,NULL,NULL));
# elif defined(CLFFT_APPLE)
CL_CH_ERR(clFFT_ExecuteInterleaved(command_queue,clplanZ,(int)3*gridY*local_gridX,(clFFT_Direction)isign,bufslices,
bufslices,0,NULL,NULL));
if (surface && isign==FFT_FORWARD) // the same operation is applied to bufslicesR, but with inverse transform
CL_CH_ERR(clFFT_ExecuteInterleaved(command_queue,clplanZ,(int)3*gridY*local_gridX,(clFFT_Direction)FFT_BACKWARD,
bufslicesR,bufslicesR,0,NULL,NULL));
# endif
#elif defined(FFTW3)
if (isign==FFT_FORWARD) {
fftw_execute(planZf);
if (surface) fftw_execute(planZRf);
}
else fftw_execute(planZb);
#elif defined(FFT_TEMPERTON)
int nn=gridZ,inc=1,jump=nn,lot=boxY,Xcomp;
IGNORE_WARNING(-Wstrict-aliasing);
for (Xcomp=0;Xcomp<3;Xcomp++) cfft99_((double *)(slices+gridYZ*Xcomp),work,trigsZ,ifaxZ,&inc,&jump,&nn,&lot,&isign);
if (surface && isign==FFT_FORWARD) { // the same operation is applied to slicesR, but with inverse transform
const int invSign=FFT_BACKWARD;
for (Xcomp=0;Xcomp<3;Xcomp++)
cfft99_((double *)(slicesR+gridYZ*Xcomp),work,trigsZ,ifaxZ,&inc,&jump,&nn,&lot,&invSign);
}
STOP_IGNORE;
#endif
}
//======================================================================================================================
static void fftX_Dm(void)
// FFT(forward) D2matrix(x) for all y,z; used for Dmatrix calculation
{
#ifdef FFTW3
fftw_execute(planXf_Dm);
#elif defined(FFT_TEMPERTON)
int nn=gridX,inc=1,jump=nn,lot=D2sizeY,isign=FFT_FORWARD;
size_t z;
IGNORE_WARNING(-Wstrict-aliasing);
for (z=0;z<lz_Dm;z++) cfft99_((double *)(D2matrix+z*gridX*D2sizeY),work,trigsX,ifaxX,&inc,&jump,&nn,&lot,&isign);
STOP_IGNORE;
#endif
}
//======================================================================================================================
static void fftX_Rm(void)
// FFT(forward) D2matrix(x) for all y,z; used for Rmatrix calculation
{
#ifdef FFTW3
fftw_execute(planXf_Rm);
#elif defined(FFT_TEMPERTON)
int nn=gridX,inc=1,jump=nn,lot=R2sizeY,isign=FFT_FORWARD;
size_t z;
const size_t zlim=local_Nz_Rm; // can be smaller by 1 than lz_Rm
IGNORE_WARNING(-Wstrict-aliasing);
for (z=0;z<zlim;z++) cfft99_((double *)(R2matrix+z*gridX*R2sizeY),work,trigsX,ifaxX,&inc,&jump,&nn,&lot,&isign);
STOP_IGNORE;
#endif
}
//======================================================================================================================
static void fftY_slice(void)
// FFT(forward) slice_tr(y) for all z; used for Dmatrix and Rmatrix calculation
{
#ifdef FFTW3
fftw_execute(planYf_slice);
#elif defined(FFT_TEMPERTON)
int nn=gridY,inc=1,jump=nn,lot=gridZ,isign=FFT_FORWARD;
IGNORE_WARNING(-Wstrict-aliasing);
cfft99_((double *)slice_tr,work,trigsY,ifaxY,&inc,&jump,&nn,&lot,&isign);
STOP_IGNORE;
#endif
}
//======================================================================================================================
static void fftZ_slice(void)
// FFT(forward) slice(z) for all y; used for Dmatrix and Rmatrix calculation
{
#ifdef FFTW3
fftw_execute(planZf_slice);
#elif defined(FFT_TEMPERTON)
int nn=gridZ,inc=1,jump=nn,lot=gridY,isign=FFT_FORWARD;
IGNORE_WARNING(-Wstrict-aliasing);
cfft99_((double *)slice,work,trigsZ,ifaxZ,&inc,&jump,&nn,&lot,&isign);
STOP_IGNORE;
#endif
}
//======================================================================================================================
void CheckNprocs(void)
// checks for consistency the specified number of processors; called in the beginning from InitComm
{
int y=nprocs;
// initialize weird_nprocs
weird_nprocs=false;
// remove simple prime divisors of y
while (y%2==0) y/=2;
#ifdef CLFFT_APPLE
// this is redundant (and CLFTT below), since OpenCL is not currently intended to run in parallel
if (y!=1) PrintError("Specified number of processors (%d) is incompatible with Apple clFFT, since the latter "
"supports only FFTs with size 2^n. Revise the number of processors or recompile with clFFT support.",nprocs);
#else
while (y%3==0) y/=3;
while (y%5==0) y/=5;
# ifdef FFT_TEMPERTON
if (y!=1) PrintError("Specified number of processors (%d) is weird (has prime divisors larger than 5). That is "
"incompatible with Temperton FFT. Revise the number of processors (recommended) or recompile with FFTW 3 "
"support.",nprocs);
# else
while (y%7==0) y/=7;
// one multiplier of either 11 or 13 is allowed
if (y%11==0) y/=11;
else if (y%13==0) y/=13;
if (y!=1) {
LogWarning(EC_WARN,ONE_POS,"Specified number of processors (%d) is weird (has prime divisors larger than 13 or "
"more than one divisor of either 11 or 13). FFTW3 will work less efficiently. It is strongly recommended "
"to revise the number of processors.",nprocs);
weird_nprocs=true;
}
# ifdef CLFFT // it allows more factors of 11 and 13, but fails completely otherwise
while (y%11==0) y/=11;
while (y%13==0) y/=13;
if (y!=1) PrintError("Specified number of processors (%d) is weird (has prime divisors larger than 13). That is "
"incompatible with clFFT. Revise the number of processors",nprocs);
# endif
# endif
#endif
}
//======================================================================================================================
int fftFit(int x,int divis)
/* find the first number >=x divisible by 2 only (Apple clFFT) or 2,3,5 only (Temperton FFT) or also
* allowing 7 and one of 11 or 13 (FFTW3) or any number of 11 and 13 factors (clFFT), and also divisible by 2 and divis.
* If weird_nprocs is used, only the latter condition is required.
* In OpenCL mode both CPU and GPU FFT routines are tested (since both are used), otherwise - only CPU ones
*/
{
int y;
if (weird_nprocs) {
if (IS_ODD(divis)) divis*=2;
return (divis*((x+divis-1)/divis));
}
else while (true) { // not very efficient but robust way
if (IS_EVEN(x) && x%divis==0) {
y=x;
while (y%2==0) y/=2; // here Apple clFFT ends
#ifndef CLFFT_APPLE // we assume that it is defined only in OpenCL mode
while (y%3==0) y/=3;
while (y%5==0) y/=5; // here Temperton FFT ends
# ifndef FFT_TEMPERTON
// the following is for FFTW3 and clFFT
while (y%7==0) y/=7;
// one multiplier of either 11 or 13 is allowed - limited by FFTW3
if (y%11==0) y/=11;
else if (y%13==0) y/=13;
# endif
#endif
if (y==1) return(x);
}
x++;
}
}
//======================================================================================================================
static void fftInitBeforeD(void)
// initialize fft before initialization of Dmatrix
{
#ifdef FFTW3
int grXint=gridX,grYint=gridY,grZint=gridZ; // this is needed to provide 'int *' to grids
/* For some reason, the following FFTW strings cannot be found by Visual Studio linker. They are not present in .def
* file for the corresponding library (in v. 3.3.5), while the default way to produce .lib file for linking with this
* linker is through the .def file (http://www.fftw.org/install/windows.html). Potentially, this can be solved by using
* 'dumpbin /exports ...' on the provided DLL, but it will probably require some manual editing.
*/
#ifndef _MSC_VER
D("FFTW library version: %s\n compiler: %s\n codelet optimizations: %s",fftw_version,fftw_cc,
fftw_codelet_optim);
#endif
planYf_slice=fftw_plan_many_dft(1,&grYint,gridZ,slice_tr,NULL,1,gridY,slice_tr,NULL,1,gridY,FFT_FORWARD,
PLAN_FFTW_DM);
planZf_slice=fftw_plan_many_dft(1,&grZint,gridY,slice,NULL,1,gridZ,slice,NULL,1,gridZ,FFT_FORWARD,PLAN_FFTW_DM);
planXf_Dm=fftw_plan_many_dft(1,&grXint,lz_Dm*D2sizeY,D2matrix,NULL,1,gridX,D2matrix,NULL,1,gridX,FFT_FORWARD,
PLAN_FFTW_DM);
// very similar to Dm, but local_Nz_Rm can be smaller by 1 than lz_Rm
if (surface) planXf_Rm=fftw_plan_many_dft(1,&grXint,local_Nz_Rm*R2sizeY,R2matrix,NULL,1,gridX,R2matrix,NULL,1,gridX,
FFT_FORWARD,PLAN_FFTW_DM);
#elif defined(FFT_TEMPERTON)
int nn;
size_t size;
// allocate memory
MALLOC_VECTOR(trigsX,double,2*gridX,ALL);
MALLOC_VECTOR(trigsY,double,2*gridY,ALL);
MALLOC_VECTOR(trigsZ,double,2*gridZ,ALL);
size=MAX(gridX*D2sizeY,3*gridYZ);
if (surface) size=MAX(size,gridX*R2sizeY);
MALLOC_VECTOR(work,double,2*size,ALL);
// initialize ifax and trigs
nn=gridX;
cftfax_(&nn,ifaxX,trigsX);
nn=gridY;
cftfax_(&nn,ifaxY,trigsY);
nn=gridZ;
cftfax_(&nn,ifaxZ,trigsZ);
#endif
}
//======================================================================================================================
static void fftInitAfterD(void)
/* second part of fft initialization
* completely separate code is used for OpenCL and FFTW3, because even precise-timing output is significantly different.
* In particular, FFTW3 uses separate plans for forward and backward, while clFFT uses one plan for both directions.
*
* clFft access the OpenCL buffers directly, so they are not anyhow affected by the definition of complex numbers in the
* main part of the code (although, it is consistent with it)
*/
{
#ifdef OPENCL
# ifdef CLFFT_APPLE
cl_int err; // error code
# else
size_t bufsize;
# endif
# ifdef PRECISE_TIMING
SYSTEM_TIME tvp[4];
# endif
if (IFROOT) PRINTFB("Initializing clFFT\n");
# ifdef PRECISE_TIMING
GET_SYSTEM_TIME(tvp);
# endif
# ifdef CLFFT
// test library version at runtime
cl_uint major,minor,patch;
CLFFT_CH_ERR(clfftGetVersion(&major,&minor,&patch));
if (!GREATER_EQ2(major,minor,CLFFT_VER_REQ,CLFFT_SUBVER_REQ)) LogError(ONE_POS,
"clFFT library version (%u.%u) is too old. Version %d.%d or newer is required",
major,minor,CLFFT_VER_REQ,CLFFT_SUBVER_REQ);
D("clFFT library version - %u.%u.%u",major,minor,patch);
// initialize library
CLFFT_CH_ERR(clfftSetup(NULL));
/* Here and further we explicitly set all plan parameters for clFFT, even those that are equal to the default
* values (as recommended in clFFT manual)
*/
/* Unfortunately, clFFT (and Apple clFFT as well) currently supports only simple regular batches of transforms
* (similar to fftw_plan_many_dft) but not fully flexible configurations, like offered by fftw_plan_guru_dft.
* Basically the problem is due to lack of multi-dimensional (non-tightly packed) batches. So to make X transform as
* a single plan we have have to cycle over the whole smallY instead of (possibly smaller) boxY. This incurs a small
* performance hit for "non-standard" values of boxY, but should be overall faster than making an explicit loop over
* smaller kernels (like is now done with Temperton FFT).
*/
CLFFT_CH_ERR(clfftCreateDefaultPlan(&clplanX,context,CLFFT_1D,&gridX));
CLFFT_CH_ERR(clfftSetPlanBatchSize(clplanX,3*local_Nz*smallY));
CLFFT_CH_ERR(clfftSetPlanPrecision(clplanX,PRECISION_CLFFT));
CLFFT_CH_ERR(clfftSetResultLocation(clplanX,CLFFT_INPLACE));
CLFFT_CH_ERR(clfftSetLayout(clplanX,CLFFT_COMPLEX_INTERLEAVED,CLFFT_COMPLEX_INTERLEAVED));
// more robust to use clfftDirection(FFT_FORWARD) instead of CLFFT_FORWARD, but here the scale is symmetric
CLFFT_CH_ERR(clfftSetPlanScale(clplanX,CLFFT_FORWARD,1));
CLFFT_CH_ERR(clfftSetPlanScale(clplanX,CLFFT_BACKWARD,1)); // override the default 1/N scale for backward direction
CLFFT_CH_ERR(clfftBakePlan(clplanX,1,&command_queue,NULL,NULL));
CLFFT_CH_ERR(clfftGetTmpBufSize(clplanX,&bufsize));
clfftBufSize+=bufsize;
# elif defined(CLFFT_APPLE)
clFFT_Dim3 xdimen;
xdimen.x=(unsigned int)gridX;
xdimen.y=1;
xdimen.z=1;
clplanX=clFFT_CreatePlan(context,xdimen,clFFT_2D,clFFT_InterleavedComplexFormat,&err);
CL_CH_ERR(err);
# endif
# ifdef PRECISE_TIMING
GET_SYSTEM_TIME(tvp+1);
# endif
# ifdef CLFFT
CLFFT_CH_ERR(clfftCreateDefaultPlan(&clplanY,context,CLFFT_1D,&gridY));
CLFFT_CH_ERR(clfftSetPlanBatchSize(clplanY,3*gridZ*local_gridX));
CLFFT_CH_ERR(clfftSetPlanPrecision(clplanY,PRECISION_CLFFT));
CLFFT_CH_ERR(clfftSetResultLocation(clplanY,CLFFT_INPLACE));
CLFFT_CH_ERR(clfftSetLayout(clplanY,CLFFT_COMPLEX_INTERLEAVED,CLFFT_COMPLEX_INTERLEAVED));
CLFFT_CH_ERR(clfftSetPlanScale(clplanY,CLFFT_FORWARD,1));
CLFFT_CH_ERR(clfftSetPlanScale(clplanY,CLFFT_BACKWARD,1)); // override the default 1/N scale for backward direction
CLFFT_CH_ERR(clfftBakePlan(clplanY,1,&command_queue,NULL,NULL));
CLFFT_CH_ERR(clfftGetTmpBufSize(clplanY,&bufsize));
clfftBufSize+=bufsize;
# elif defined(CLFFT_APPLE)
clFFT_Dim3 ydimen;
ydimen.x=(unsigned int)gridY;
ydimen.y=1;
ydimen.z=1;
clplanY=clFFT_CreatePlan(context,ydimen,clFFT_1D,clFFT_InterleavedComplexFormat,&err);
CL_CH_ERR(err);
# endif
# ifdef PRECISE_TIMING
GET_SYSTEM_TIME(tvp+2);
# endif
# ifdef CLFFT
/* Here the issue is similar to clplanX described above. However, we are using full gridY instead of boxY, which
* incurs at least a-factor-of-two performance hit. To solve this problem one need to execute separate plans for 3
* components of vectors. Unfortunately, this cannot be simply done using 3-element loop (like in Temperton FFT),
* because a part of cl_mem object can't be addressed independently (as can be done with simple C arrays). The only
* way to address this issue is to either create three separate cl_mem objects or to change the indexing of levels
* inside the array, so that 3 components are stored together.
*/
CLFFT_CH_ERR(clfftCreateDefaultPlan(&clplanZ,context,CLFFT_1D,&gridZ));
/* TODO: last slices can be very slightly thinner than the previous ones, but since the batchsize is part of the
* plan, another plan would be needed to address this. However, we ignore this currently and assume that every
* slice has a thickness of local_gridX.
* This issue also applies to clplanY.
*/
CLFFT_CH_ERR(clfftSetPlanBatchSize(clplanZ,3*gridY*local_gridX));
CLFFT_CH_ERR(clfftSetPlanPrecision(clplanZ,PRECISION_CLFFT));
CLFFT_CH_ERR(clfftSetResultLocation(clplanZ,CLFFT_INPLACE));
CLFFT_CH_ERR(clfftSetLayout(clplanZ,CLFFT_COMPLEX_INTERLEAVED,CLFFT_COMPLEX_INTERLEAVED));
CLFFT_CH_ERR(clfftSetPlanScale(clplanZ,CLFFT_FORWARD,1));
CLFFT_CH_ERR(clfftSetPlanScale(clplanZ,CLFFT_BACKWARD,1)); // override the default 1/N scale for backward direction
CLFFT_CH_ERR(clfftBakePlan(clplanZ,1,&command_queue,NULL,NULL));
CLFFT_CH_ERR(clfftGetTmpBufSize(clplanZ,&bufsize));
clfftBufSize+=bufsize;
/* In most cases clfftBufSize is zero, except some weird grid sizes like 2x2x60000. Still, we rigorously account
* for this memory. However, we do not update oclMemMaxObj, since even single plan is not guaranteed to allocate a
* single object. So we assume that clFFT will either handle maximum object size itself or produce a meaningful
* error.
*/
oclMem+=clfftBufSize;
MAXIMIZE(oclMemPeak,oclMem);
# elif defined(CLFFT_APPLE)
clFFT_Dim3 zdimen;
zdimen.x=(unsigned int)gridZ;
zdimen.y=1;
zdimen.z=1;
clplanZ=clFFT_CreatePlan(context,zdimen,clFFT_1D,clFFT_InterleavedComplexFormat,&err);
CL_CH_ERR(err);
# endif
# ifdef PRECISE_TIMING
GET_SYSTEM_TIME(tvp+3);
// print precise timing of FFT planning
if (IFROOT) PrintBoth(logfile,
"~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~\n"
" clFFT planning \n"
"~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~\n"
"X = "FFORMPT" Total = "FFORMPT"\n"
"Y = "FFORMPT"\n"
"Z = "FFORMPT"\n\n",
DiffSystemTime(tvp,tvp+1),DiffSystemTime(tvp,tvp+3),DiffSystemTime(tvp+1,tvp+2),DiffSystemTime(tvp+2,tvp+3));
# endif
#elif defined(FFTW3) // this is not needed when OpenCL is used
int lot;
fftw_iodim dims,howmany_dims[2];
int grYint=gridY; // this is needed to provide 'int *' to gridY
# ifdef PRECISE_TIMING
SYSTEM_TIME tvp[7];
# endif
if (IFROOT) PRINTFB("Initializing FFTW3\n");
# ifdef PRECISE_TIMING
GET_SYSTEM_TIME(tvp);
# endif
lot=3*gridZ;
planYf=fftw_plan_many_dft(1,&grYint,lot,slices_tr,NULL,1,gridY,slices_tr,NULL,1,gridY,FFT_FORWARD,PLAN_FFTW);
if (surface) // same operation, but applied to slicesR_tr
planYRf=fftw_plan_many_dft(1,&grYint,lot,slicesR_tr,NULL,1,gridY,slicesR_tr,NULL,1,gridY,FFT_FORWARD,PLAN_FFTW);
# ifdef PRECISE_TIMING
GET_SYSTEM_TIME(tvp+1);
# endif
planYb=fftw_plan_many_dft(1,&grYint,lot,slices_tr,NULL,1,gridY,slices_tr,NULL,1,gridY,FFT_BACKWARD,PLAN_FFTW);
# ifdef PRECISE_TIMING
GET_SYSTEM_TIME(tvp+2);
# endif
dims.n=gridZ;
dims.is=dims.os=1;
howmany_dims[0].n=3;
howmany_dims[0].is=howmany_dims[0].os=gridZ*gridY;
howmany_dims[1].n=boxY;
howmany_dims[1].is=howmany_dims[1].os=gridZ;
planZf=fftw_plan_guru_dft(1,&dims,2,howmany_dims,slices,slices,FFT_FORWARD,PLAN_FFTW);
// same operation but for slicesR and inverse transform (since correlation is computed instead of convolution)
if (surface) planZRf=fftw_plan_guru_dft(1,&dims,2,howmany_dims,slicesR,slicesR,FFT_BACKWARD,PLAN_FFTW);
# ifdef PRECISE_TIMING
GET_SYSTEM_TIME(tvp+3);
# endif
planZb=fftw_plan_guru_dft(1,&dims,2,howmany_dims,slices,slices,FFT_BACKWARD,PLAN_FFTW);
# ifdef PRECISE_TIMING
GET_SYSTEM_TIME(tvp+4);
# endif
dims.n=gridX;
dims.is=dims.os=1;
howmany_dims[0].n=3*local_Nz;
howmany_dims[0].is=howmany_dims[0].os=smallY*gridX;
howmany_dims[1].n=boxY;
howmany_dims[1].is=howmany_dims[1].os=gridX;
planXf=fftw_plan_guru_dft(1,&dims,2,howmany_dims,Xmatrix,Xmatrix,FFT_FORWARD,PLAN_FFTW);
# ifdef PRECISE_TIMING
GET_SYSTEM_TIME(tvp+5);
# endif
planXb=fftw_plan_guru_dft(1,&dims,2,howmany_dims,Xmatrix,Xmatrix,FFT_BACKWARD,PLAN_FFTW);
# ifdef PRECISE_TIMING
GET_SYSTEM_TIME(tvp+6);
// print precise timing of FFT planning
if (IFROOT) PrintBoth(logfile,
"~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~\n"
" FFTW3 planning \n"
"~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~\n"
"Yf = "FFORMPT" Total = "FFORMPT"\n"
"Yb = "FFORMPT"\n"
"Zf = "FFORMPT"\n"
"Zb = "FFORMPT"\n"
"Xf = "FFORMPT"\n"
"Xb = "FFORMPT"\n\n",
DiffSystemTime(tvp,tvp+1),DiffSystemTime(tvp,tvp+6),DiffSystemTime(tvp+1,tvp+2),DiffSystemTime(tvp+2,tvp+3),
DiffSystemTime(tvp+3,tvp+4),DiffSystemTime(tvp+4,tvp+5),DiffSystemTime(tvp+5,tvp+6));
# endif
#endif
#ifdef FFTW3
// destroy old (D,R-matrix) plans; also in OpenCL mode
fftw_destroy_plan(planXf_Dm);
fftw_destroy_plan(planYf_slice);
fftw_destroy_plan(planZf_slice);
if (surface) fftw_destroy_plan(planXf_Rm);
# ifdef OPENCL // in this case, FFTW ends here
fftw_cleanup();
# endif
#endif
}
//======================================================================================================================
static void InitRmatrix(const double invNgrid)
/* Initializes the matrix R. R[i][j][k]=GR[i1-i2][j1-j2][k1+k2]. Actually R=-FFT(GR)/Ngrid. Then -GR.x=invFFT(R*FFT(x))
* for practical implementation of FFT such that invFFT(FFT(x))=Ngrid*x. GR is exactly reflected Green's tensor. The
* routine is very similar to the corresponding part of InitDmatrix. Moreover, some initialization is delegated to the
* latter function.
*/
{
int i,j,k,Rcomp;
size_t x,y,z,indexfrom,indexto,ind,index;
// allocate memory for Rmatrix (R2matrix is allocated earlier in InitDmatrix)
MALLOC_VECTOR(Rmatrix,complex,Rsize,ALL);
#ifdef PARALLEL
// allocate buffer for BlockTranspose_DRm
size_t bufsize = 2*lz_Rm*R2sizeY*local_Nx;
MALLOC_VECTOR(BT_buffer,double,bufsize,ALL);
MALLOC_VECTOR(BT_rbuffer,double,bufsize,ALL);
#endif
if (IFROOT) PRINTFB("Calculating reflected Green's function (Rmatrix)\n");
/* Interaction matrix values are calculated all at once for performance reasons. They are stored in Rmatrix with
* indexing corresponding to R2matrix (to facilitate copying) but NDCOMP elements instead of one. Afterwards they
* are replaced by Fourier transforms (with different indexing) component-wise (in cycle over NDCOMP)
*/
/* fill Rmatrix with 0, this if to fill the possible gap between e.g. boxY and gridY/2; (and for R=0) probably
* faster than using a lot of conditionals
*/
for (ind=0;ind<Rsize;ind++) Rmatrix[ind]=0;
// fill Rmatrix with values of reflected Green's tensor
for(k=0;k<local_Nz_Rm;k++) for (j=jstartR;j<boxY;j++) for (i=1-boxX;i<boxX;i++) {
index=NDCOMP*Index2matrix(i,j,k,R2sizeY);
(*ReflTerm_int)(i,j,k,Rmatrix+index);
} // end of i,j,k loop
if (IFROOT) PRINTFB("Fourier transform of Rmatrix\n");
for(Rcomp=0;Rcomp<NDCOMP;Rcomp++) { // main cycle over components of Rmatrix
// fill R2matrix with precomputed values from Rmatrix
for (ind=0;ind<R2sizeTot;ind++) R2matrix[ind]=Rmatrix[NDCOMP*ind+Rcomp];
fftX_Rm(); // fftX R2matrix
BlockTranspose_DRm(R2matrix,R2sizeY,lz_Rm);
for(x=local_x0;x<local_x1;x++) {
for (ind=0;ind<gridYZ;ind++) slice[ind]=0.0; // fill slice with 0.0
for(j=jstartR;j<boxY;j++) for(k=0;k<2*boxZ-1;k++) {
indexfrom=IndexGarbledR(x,j,k);
indexto=IndexSliceR2matrix(j,k);
slice[indexto]=R2matrix[indexfrom];
}
/* here a specific symmetry is used, that elements of R depend on direction y/|rho| either as even order
* (0 or 2) or as odd (1) - the latter are elements 1 and 4 (see GetSomIntegral in interaction.c)
*/
if (reduced_FFT) for(j=1;j<boxY;j++) for(k=0;k<2*boxZ-1;k++) {
// mirror along y
indexfrom=IndexSliceR2matrix(j,k);
indexto=IndexSliceR2matrix(-j,k);
if (Rcomp==1 || Rcomp==4) slice[indexto]=-slice[indexfrom];
else slice[indexto]=slice[indexfrom];
}
fftZ_slice(); // fftZ slice
transpose(slice,slice_tr,gridY,gridZ);
fftY_slice(); // fftY slice_tr
for(z=0;z<gridZ;z++) for(y=0;y<RsizeY;y++) {
indexto=IndexRmatrix(x-local_x0,y,z)+Rcomp;
indexfrom=IndexSlice_zy(y,z);
Rmatrix[indexto]=-invNgrid*slice_tr[indexfrom];
}
} // end slice X
} // end of Rcomp
#ifdef PARALLEL
// deallocate buffers for BlockTranspose_DRm
Free_general(BT_buffer);
Free_general(BT_rbuffer);
#endif
#ifdef OPENCL
// Setting kernel arguments which are always the same
// for arith3_surface
CL_CH_ERR(clSetKernelArg(clarith3_surface,0,sizeof(cl_mem),&bufslices_tr));
CL_CH_ERR(clSetKernelArg(clarith3_surface,1,sizeof(cl_mem),&bufDmatrix));
CL_CH_ERR(clSetKernelArg(clarith3_surface,2,sizeof(size_t),&smallY));
CL_CH_ERR(clSetKernelArg(clarith3_surface,3,sizeof(size_t),&smallZ));
CL_CH_ERR(clSetKernelArg(clarith3_surface,4,sizeof(size_t),&gridX));
CL_CH_ERR(clSetKernelArg(clarith3_surface,5,sizeof(size_t),&DsizeY));
CL_CH_ERR(clSetKernelArg(clarith3_surface,6,sizeof(size_t),&DsizeZ));
CL_CH_ERR(clSetKernelArg(clarith3_surface,11,sizeof(cl_mem),&bufslicesR_tr));
CL_CH_ERR(clSetKernelArg(clarith3_surface,12,sizeof(cl_mem),&bufRmatrix));
// for transpose forward (backward are not needed for surface)
CL_CH_ERR(clSetKernelArg(cltransposeofR,0,sizeof(cl_mem),&bufslicesR));
CL_CH_ERR(clSetKernelArg(cltransposeofR,1,sizeof(cl_mem),&bufslicesR_tr));
CL_CH_ERR(clSetKernelArg(cltransposeofR,2,sizeof(size_t),&gridZ));
CL_CH_ERR(clSetKernelArg(cltransposeofR,3,sizeof(size_t),&gridY));
CL_CH_ERR(clSetKernelArg(cltransposeofR,4,17*16*sizeof(doublecomplex),NULL));
// copy Rmatrix to OpenCL buffer, blocking to ensure completion before function end
CL_CH_ERR(clEnqueueWriteBuffer(command_queue,bufRmatrix,CL_TRUE,0,Rsize*sizeof(*Rmatrix),Rmatrix,0,NULL,NULL));
Free_cVector(Rmatrix);
#endif
}
//======================================================================================================================
void InitDmatrix(void)
/* Initializes the matrix D. D[i][j][k]=A[i1-i2][j1-j2][k1-k2]. Actually D=-FFT(G)/Ngrid. Then -G.x=invFFT(D*FFT(x)) for
* practical implementation of FFT such that invFFT(FFT(x))=Ngrid*x. G is exactly Green's tensor. The routine is called
* only once, so does not need to be very fast, however we tried to optimize it.
*/
{
int i,j,k,kcor,Dcomp;
size_t x,y,z,indexfrom,indexto,ind,index,Dsize,D2sizeTot;
double invNgrid;
int nnn; // multiplier used for reduced_FFT or not reduced; 1 or 2
int jstart,kstart;
TIME_TYPE start,time1;
#ifdef PRECISE_TIMING
// precise timing of the Dmatrix computation
SYSTEM_TIME tvp[15];
SYSTEM_TIME Timing_fftX,Timing_fftY,Timing_fftZ,Timing_Gcalc,Timing_ar1,Timing_ar2,Timing_ar3,Timing_BT,Timing_TYZ,
Timing_beg,Timing_InitMV;
double t_fftX,t_fftY,t_fftZ,t_ar1,t_ar2,t_ar3,t_TYZ,t_beg,t_Gcalc,t_Arithm,t_FFT,t_BT,t_InitMV,t_Rm,t_Tot;
// This should be the first occurrence of PRECISE_TIMING in the program
SetTimerFreq();