-
Notifications
You must be signed in to change notification settings - Fork 1
/
Copy pathLATfield2_PlanFFT.hpp
executable file
·1718 lines (1312 loc) · 57.3 KB
/
LATfield2_PlanFFT.hpp
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
#ifndef LATFIELD2_PLANFFT_HPP
#define LATFIELD2_PLANFFT_HPP
/*! \file LATfield2_PlanFFT.hpp
\brief LATfield2_PlanFFT.hpp contains the class PlanFFT definition.
\author David Daverio
*/
#include "fftw3.h"
#ifdef SINGLE
#define MPI_DATA_PREC MPI_FLOAT
#endif
#ifndef SINGLE
#define MPI_DATA_PREC MPI_DOUBLE
#endif
const int FFT_FORWARD = FFTW_FORWARD;
const int FFT_BACKWARD = FFTW_BACKWARD;
const int FFT_IN_PLACE = 16;
const int FFT_OUT_OF_PLACE = -16;
/*! \class temporaryMemFFT
\brief A class wich handle the additional memory needed by the class PlanFFT; Should never be used, internal temporary memroy handler, if real need, can be hacked, but could comflict with the PlanFFT!
*/
class temporaryMemFFT
{
public:
temporaryMemFFT();
~temporaryMemFFT();
temporaryMemFFT(long size);
int setTemp(long size);
#ifdef SINGLE
fftwf_complex* temp1(){return temp1_;}
fftwf_complex* temp2(){return temp2_;}
#endif
#ifndef SINGLE
fftw_complex * temp1(){return temp1_;}
fftw_complex * temp2(){return temp2_;}
#endif
private:
#ifdef SINGLE
fftwf_complex * temp1_;
fftwf_complex * temp2_;
#endif
#ifndef SINGLE
fftw_complex * temp1_;
fftw_complex * temp2_;
#endif
long allocated_; //number of variable stored (bit = allocated*sizeof(fftw(f)_complex))
};
temporaryMemFFT::temporaryMemFFT()
{
allocated_=0;
}
temporaryMemFFT::~temporaryMemFFT()
{
if(allocated_!=0){
#ifdef SINGLE
fftwf_free(temp1_);
fftwf_free(temp2_);
#endif
#ifndef SINGLE
fftw_free(temp1_);
fftw_free(temp2_);
#endif
}
}
temporaryMemFFT::temporaryMemFFT(long size)
{
#ifdef SINGLE
temp1_ = (fftwf_complex *)fftwf_malloc(size*sizeof(fftwf_complex));
temp2_ = (fftwf_complex *)fftwf_malloc(size*sizeof(fftwf_complex));
#endif
#ifndef SINGLE
temp1_ = (fftw_complex *)fftw_malloc(size*sizeof(fftw_complex));
temp2_ = (fftw_complex *)fftw_malloc(size*sizeof(fftw_complex));
#endif
allocated_=size;
}
int temporaryMemFFT::setTemp(long size)
{
if(size>allocated_)
{
#ifdef SINGLE
if(allocated_!=0){
fftwf_free(temp1_);
fftwf_free(temp2_);
}
temp1_ = (fftwf_complex *)fftwf_malloc(size*sizeof(fftwf_complex));
temp2_ = (fftwf_complex *)fftwf_malloc(size*sizeof(fftwf_complex));
//debug cout<<"("<< parallel.grid_rank()[0]<< ","<<parallel.grid_rank()[1] <<"): called temporary resize. Old size: " <<allocated_<<" , new size: "<< size<<endl;
#endif
#ifndef SINGLE
if(allocated_!=0){
fftw_free(temp1_);
fftw_free(temp2_);
}
temp1_ = (fftw_complex *)fftw_malloc(size*sizeof(fftw_complex));
temp2_ = (fftw_complex *)fftw_malloc(size*sizeof(fftw_complex));
#endif
allocated_ = size ;
}
return 1;
}
//////////////////////Temp memory///////////////////////////
temporaryMemFFT tempMemory;
/*! \class PlanFFT
\brief Class which handle Fourier transforms of fields (real/complex, single/double precision) on cubic lattices. (currently implemented only for 3d)
This class allow to perform Fourier transform of real and complex fields. See the FFTs example to have have a short intro of the usage. See the PoissonSolver example for more advanced usage (as linking several field to the same Fourier image)
One should understand that first a plan is created then executed (in the FFTW fashion). The plan links two fields, one on Fourier space, one on real space. Boths field will be allocated by the planer when the plan is initilized. But need to be initialized before passing them to the planer.
One need to be carefull to corretly define the lattice and field.
\sa void Lattice::initializeRealFFT(Lattice & lat_real, int halo);
\sa void Lattice::initializeComplexFFT(Lattice & lat_real, int halo);
\sa cKSite class
\sa rKSite class
*/
template<class compType>
class PlanFFT
{
public:
//! Constructor.
PlanFFT();
#ifndef SINGLE
/*!
Constructor with initialization for complex to complex transform.
\param rfield : real space field
\param kfield : Fourier space field
\param mem_type : memory type (FFT_OUT_OF_PLACE or FFT_IN_PLACE). In place mean that both Fourier and real space field point to the same data array.
\sa initialize(Field<compType>* rfield,Field<compType>* kfield,const int mem_type = FFT_OUT_OF_PLACE);
*/
PlanFFT(Field<compType>* rfield, Field<compType>* kfield,const int mem_type = FFT_OUT_OF_PLACE);
/*!
Initialization for complex to complex transform.
\param rfield : real space field
\param kfield : Fourier space field
\param mem_type : memory type (FFT_OUT_OF_PLACE or FFT_IN_PLACE). In place mean that both Fourier and real space field point to the same data array.
*/
void initialize(Field<compType>* rfield,Field<compType>* kfield,const int mem_type = FFT_OUT_OF_PLACE);
/*!
Constructor with initialization for real to complex transform.
\param rfield : real space field
\param kfield : Fourier space field
\param mem_type : memory type (FFT_OUT_OF_PLACE or FFT_IN_PLACE). In place mean that both Fourier and real space fields point to the same data array.
\sa initialize(Field<compType>* rfield,Field<compType>* kfield,const int mem_type = FFT_OUT_OF_PLACE);
*/
PlanFFT(Field<double>* rfield, Field<compType>* kfield,const int mem_type = FFT_OUT_OF_PLACE);
/*!
Initialization for real to complex transform.
\param rfield : real space field
\param kfield : Fourier space field
\param mem_type : memory type (FFT_OUT_OF_PLACE or FFT_IN_PLACE). In place mean that both Fourier and real space fields point to the same data array.
*/
void initialize(Field<double>* rfield,Field<compType>* kfield,const int mem_type = FFT_OUT_OF_PLACE);
#endif
#ifdef SINGLE
PlanFFT(Field<compType>* rfield, Field<compType>* kfield,const int mem_type = FFT_OUT_OF_PLACE);
void initialize(Field<compType> * rfield,Field<compType> * kfield,const int mem_type = FFT_OUT_OF_PLACE);
PlanFFT(Field<float>* rfield, Field<compType>* kfield,const int mem_type = FFT_OUT_OF_PLACE);
void initialize(Field<float>* rfield,Field<compType>* kfield,const int mem_type = FFT_OUT_OF_PLACE);
#endif
/*!
Execute the Fourier transform.
\param fft_type : dirrection of the transform. Can be FFT_BACKWARD or FFT_FORWARD.
*/
void execute(int fft_type);
private:
bool status_;
bool type_;
int mem_type_;
static bool R2C;
static bool C2C;
static bool initialized;
//data description variable, fftw plan, and temp
int components_;
int rSize_[3];
int kSize_[3];
int rJump_[3];
int kJump_[3];
int rSizeLocal_[3];
int kSizeLocal_[3];
int r2cSize_; //of 0 dimension
int r2cSizeLocal_; //of r2csize on proc_dim[1]
int r2cSizeLocal_as_;
int rHalo_;
int kHalo_;
#ifdef SINGLE
float * rData_; //pointer to start of data (halo skip)
fftwf_complex * cData_; //pointer to start of data (halo skip)
fftwf_complex * kData_; //pointer to start of data (halo skip)
fftwf_complex * temp_;
fftwf_complex * temp1_;//needed if field got more than 1 component
fftwf_plan fPlan_i_;
fftwf_plan fPlan_j_;
fftwf_plan fPlan_k_;
fftwf_plan fPlan_k_real_;
fftwf_plan bPlan_i_;
fftwf_plan bPlan_j_;
fftwf_plan bPlan_j_real_;
fftwf_plan bPlan_k_;
//transpostion fonction
// forward real to complex
// first transopsition
void transpose_0_2( fftwf_complex * in, fftwf_complex * out,int dim_i,int dim_j ,int dim_k);
void transpose_0_2_last_proc( fftwf_complex * in, fftwf_complex * out,int dim_i,int dim_j ,int dim_k);
void implement_local_0_last_proc( fftwf_complex * in, fftwf_complex * out,int proc_dim_i,int proc_dim_j,int proc_dim_k,int proc_size);
// second transposition
void transpose_1_2(fftwf_complex * in , fftwf_complex * out ,int dim_i,int dim_j ,int dim_k);
//third transposition
void transpose_back_0_3(fftwf_complex * in, fftwf_complex * out,int r2c,int local_r2c,int local_size_j,int local_size_k,int proc_size,int halo,int components,int comp);
void implement_0(fftwf_complex * in, fftwf_complex * out,int r2c_size,int local_size_j,int local_size_k,int halo,int components,int comp);
//backward real to complex
void b_arrange_data_0(fftwf_complex *in, fftwf_complex * out,int dim_i,int dim_j ,int dim_k, int khalo, int components, int comp);
void b_transpose_back_0_1(fftwf_complex * in, fftwf_complex * out,int r2c,int local_r2c,int local_size_j,int local_size_k,int proc_size);
void b_implement_0(fftwf_complex * in, fftwf_complex * out,int r2c_size,int local_size_j,int local_size_k);
#endif
#ifndef SINGLE
double * rData_; //pointer to start of data (halo skip)
fftw_complex * cData_; //pointer to start of data (halo skip)
fftw_complex * kData_; //pointer to start of data (halo skip)
fftw_complex * temp_;
fftw_complex * temp1_;//needed if field got more than 1 component
fftw_plan fPlan_i_;
fftw_plan fPlan_j_;
fftw_plan fPlan_k_;
fftw_plan fPlan_k_real_;
fftw_plan bPlan_i_;
fftw_plan bPlan_j_;
fftw_plan bPlan_j_real_;
fftw_plan bPlan_k_;
///transpostion fonction
/// forward real to complex
// first transopsition
void transpose_0_2( fftw_complex * in, fftw_complex * out,int dim_i,int dim_j ,int dim_k);
void transpose_0_2_last_proc( fftw_complex * in, fftw_complex * out,int dim_i,int dim_j ,int dim_k);
void implement_local_0_last_proc( fftw_complex * in, fftw_complex * out,int proc_dim_i,int proc_dim_j,int proc_dim_k,int proc_size);
// second transposition
void transpose_1_2(fftw_complex * in , fftw_complex * out ,int dim_i,int dim_j ,int dim_k);
//third transposition
void transpose_back_0_3(fftw_complex * in, fftw_complex * out,int r2c,int local_r2c,int local_size_j,int local_size_k,int proc_size,int halo,int components,int comp);
void implement_0(fftw_complex * in, fftw_complex * out,int r2c_size,int local_size_j,int local_size_k,int halo,int components,int comp);
////backward real to complex
void b_arrange_data_0(fftw_complex *in, fftw_complex * out,int dim_i,int dim_j ,int dim_k, int khalo, int components, int comp);
void b_transpose_back_0_1(fftw_complex * in, fftw_complex * out,int r2c,int local_r2c,int local_size_j,int local_size_k,int proc_size);
void b_implement_0(fftw_complex * in, fftw_complex * out,int r2c_size,int local_size_j,int local_size_k);
#endif
};
//constants
template<class compType>
bool PlanFFT<compType>::initialized = true;
template<class compType>
bool PlanFFT<compType>::R2C=false;
template<class compType>
bool PlanFFT<compType>::C2C=true;
template<class compType>
PlanFFT<compType>::PlanFFT()
{
status_ = false;
}
#ifdef SINGLE
template<class compType>
PlanFFT<compType>::PlanFFT(Field<compType>* rfield,Field<compType>* kfield,const int mem_type)
{
status_ = false;
initialize(rfield,kfield,mem_type);
}
template<class compType>
void PlanFFT<compType>::initialize(Field<compType>* rfield,Field<compType>* kfield,const int mem_type )
{
type_ = C2C;
mem_type_=mem_type;
//general variable
if(rfield->components() != kfield->components())
{
cerr<<"Latfield2d::PlanFFT::initialize : fft curently work only for fields with same number of components"<<endl;
cerr<<"Latfield2d::PlanFFT::initialize : coordinate and Fourier space fields have not the same number of components"<<endl;
cerr<<"Latfield2d : Abort Process Requested"<<endl;
}
else components_ = rfield->components();
for(int i = 0; i<3; i++)
{
rSize_[i]=rfield->lattice().size(i);
kSize_[i]=kfield->lattice().size(i);
rSizeLocal_[i]=rfield->lattice().sizeLocal(i);
kSizeLocal_[i]=kfield->lattice().sizeLocal(i);
rJump_[i]=rfield->lattice().jump(i);
kJump_[i]=kfield->lattice().jump(i);
}
rHalo_ = rfield->lattice().halo();
kHalo_ = kfield->lattice().halo();
tempMemory.setTemp((long)(rSize_[0]) * (long)(rSizeLocal_[1]) * (long)(rSizeLocal_[2]));
temp_ = tempMemory.temp1();
temp1_ = tempMemory.temp2();
if(rfield->lattice().dim()!=3)
{
if(parallel.isRoot())
{
cerr<<"Latfield2d::PlanFFT::initialize : fft curently work only for 3d cubic lattice"<<endl;
cerr<<"Latfield2d::PlanFFT::initialize : real lattice have not 3 dimensions"<<endl;
cerr<<"Latfield2d : Abort Process Requested"<<endl;
}
parallel.abortForce();
}
if(rSize_[0]!=rSize_[1] | rSize_[1]!=rSize_[2])
{
if(parallel.isRoot())
{
cerr<<"Latfield2d::PlanFFT::initialize : fft curently work only for 3d cubic lattice"<<endl;
cerr<<"Latfield2d::PlanFFT::initialize : real lattice is not cubic"<<endl;
cerr<<"Latfield2d : Abort Process Requested"<<endl;
}
parallel.abortForce();
}
//initialization of fftw plan
//Forward plan
fPlan_i_ = fftwf_plan_many_dft(1,&rSize_[0],rSizeLocal_[1] ,cData_,NULL,components_, rJump_[1]*components_,temp_,NULL,rSizeLocal_[1]*rSizeLocal_[2],1,FFTW_FORWARD,FFTW_ESTIMATE | FFTW_PRESERVE_INPUT);
fPlan_j_ = fftwf_plan_many_dft(1,&rSize_[0],rSizeLocal_[1]*rSizeLocal_[2],temp_,NULL,rSizeLocal_[1]*rSizeLocal_[2],1,temp_,NULL,rSizeLocal_[1]*rSizeLocal_[2],1,FFTW_FORWARD,FFTW_ESTIMATE);
fPlan_k_ = fftwf_plan_many_dft(1,&rSize_[0],rSizeLocal_[1] ,temp_,NULL,rSizeLocal_[1]*rSizeLocal_[2],1,kData_,NULL,components_, rJump_[1]*components_,FFTW_FORWARD,FFTW_ESTIMATE);
//Backward plan
bPlan_k_ = fftwf_plan_many_dft(1,&rSize_[0],rSizeLocal_[1] ,kData_,NULL,components_, rJump_[1]*components_,temp_,NULL,rSizeLocal_[1]*rSizeLocal_[2],1,FFTW_BACKWARD,FFTW_ESTIMATE | FFTW_PRESERVE_INPUT);
bPlan_j_ = fftwf_plan_many_dft(1,&rSize_[0],rSizeLocal_[1]*rSizeLocal_[2],temp_,NULL,rSizeLocal_[1]*rSizeLocal_[2],1,temp_,NULL,rSizeLocal_[1]*rSizeLocal_[2],1,FFTW_BACKWARD,FFTW_ESTIMATE);
bPlan_i_ = fftwf_plan_many_dft(1,&rSize_[0],rSizeLocal_[1] ,temp_,NULL,rSizeLocal_[1]*rSizeLocal_[2],1,kData_,NULL,components_, rJump_[1]*components_,FFTW_BACKWARD,FFTW_ESTIMATE);
//allocation of field
long rfield_size = rfield->lattice().sitesLocalGross();
long kfield_size = kfield->lattice().sitesLocalGross();
if(mem_type_ == FFT_IN_PLACE)
{
if(rfield_size>=kfield_size)
{
rfield->alloc();
kfield->data() = (Imag*)rfield->data();
}
else
{
kfield->alloc();
rfield->data() = (Imag*)kfield->data();
}
}
if(mem_type_ == FFT_OUT_OF_PLACE)
{
rfield->alloc();
kfield->alloc();
}
//Pointer to data
rData_ = (float*)rfield->data(); //to be sure that rData is instantiate !
cData_ = (fftwf_complex*)rfield->data();
cData_ += rfield->lattice().siteFirst()*components_;
kData_ = (fftwf_complex*)kfield->data();
kData_ += kfield->lattice().siteFirst()*components_;
}
template<class compType>
PlanFFT<compType>::PlanFFT(Field<float>* rfield, Field<compType>* kfield,const int mem_type )
{
status_ = false;
initialize(rfield,kfield,mem_type);
}
template<class compType>
void PlanFFT<compType>::initialize(Field<float>* rfield,Field<compType>* kfield,const int mem_type )
{
type_ = R2C;
mem_type_=mem_type;
//general variable
if(rfield->components() != kfield->components())
{
cerr<<"Latfield2d::PlanFFT::initialize : fft curently work only for fields with same number of components"<<endl;
cerr<<"Latfield2d::PlanFFT::initialize : coordinate and Fourier space fields have not the same number of components"<<endl;
cerr<<"Latfield2d : Abort Process Requested"<<endl;
}
else components_ = rfield->components();
for(int i = 0; i<3; i++)
{
rSize_[i]=rfield->lattice().size(i);
kSize_[i]=kfield->lattice().size(i);
rSizeLocal_[i]=rfield->lattice().sizeLocal(i);
kSizeLocal_[i]=kfield->lattice().sizeLocal(i);
rJump_[i]=rfield->lattice().jump(i);
kJump_[i]=kfield->lattice().jump(i);
}
r2cSize_ = rfield->lattice().size(0)/2 + 1;
r2cSizeLocal_as_ = rfield->lattice().sizeLocal(0)/(2*parallel.grid_size()[1]);
if(parallel.last_proc()[1]) r2cSizeLocal_ = r2cSizeLocal_as_ + 1;
else r2cSizeLocal_ = r2cSizeLocal_as_;
rHalo_ = rfield->lattice().halo();
kHalo_ = kfield->lattice().halo();
tempMemory.setTemp((r2cSize_+2) * (rSizeLocal_[1]+2) * (rSizeLocal_[2]+2));
temp_ = tempMemory.temp1();
temp1_ = tempMemory.temp2();
if(rfield->lattice().dim()!=3)
{
if(parallel.isRoot())
{
cerr<<"Latfield2d::PlanFFT::initialize : fft curently work only for 3d cubic lattice"<<endl;
cerr<<"Latfield2d::PlanFFT::initialize : real lattice have not 3 dimensions"<<endl;
cerr<<"Latfield2d : Abort Process Requested"<<endl;
}
parallel.abortForce();
}
//COUT << rSize_[0] << " "<< rSize_[1]<< " "<< rSize_[2] <<endl;
if(rSize_[0]!=rSize_[1] | rSize_[1]!=rSize_[2])
{
if(parallel.isRoot())
{
cerr<<"Latfield2d::PlanFFT::initialize : fft curently work only for 3d cubic lattice"<<endl;
cerr<<"Latfield2d::PlanFFT::initialize : real lattice is not cubic"<<endl;
cerr<<"Latfield2d : Abort Process Requested"<<endl;
}
parallel.abortForce();
}
//create the fftw_plan
fPlan_i_ = fftwf_plan_many_dft_r2c(1,&rSize_[0],rSizeLocal_[1] ,rData_,NULL,components_, rJump_[1]*components_,temp_,NULL,rSizeLocal_[1]*rSizeLocal_[2],1,FFTW_ESTIMATE | FFTW_PRESERVE_INPUT);
fPlan_j_ = fftwf_plan_many_dft(1,&rSize_[0],rSizeLocal_[2]*r2cSizeLocal_,temp_,NULL,rSizeLocal_[2]*r2cSizeLocal_,1,temp_,NULL,rSizeLocal_[2]*r2cSizeLocal_,1,FFTW_FORWARD,FFTW_ESTIMATE);
fPlan_k_ = fftwf_plan_many_dft(1,&rSize_[0],r2cSizeLocal_as_,temp_,NULL,rSizeLocal_[2]*r2cSizeLocal_,1,temp1_,NULL,rSizeLocal_[2]*r2cSizeLocal_as_,1,FFTW_FORWARD,FFTW_ESTIMATE);
fPlan_k_real_ = fftwf_plan_many_dft(1,&rSize_[0],rSizeLocal_[2],&temp_[r2cSizeLocal_as_],NULL,rSizeLocal_[2]*r2cSizeLocal_,r2cSizeLocal_,&temp1_[r2cSizeLocal_as_*rSizeLocal_[2]*rSize_[0]],NULL,rSizeLocal_[2],1,FFTW_FORWARD,FFTW_ESTIMATE);
bPlan_k_ = fftwf_plan_many_dft(1,&rSize_[0],kSizeLocal_[2]*r2cSizeLocal_,temp_,NULL,kSizeLocal_[2]*r2cSizeLocal_,1,temp_,NULL,kSizeLocal_[2]*r2cSizeLocal_,1,FFTW_BACKWARD,FFTW_ESTIMATE);
bPlan_j_ = fftwf_plan_many_dft(1,&rSize_[0],r2cSizeLocal_as_,temp_,NULL,rSizeLocal_[2]*r2cSizeLocal_,1,temp1_,NULL,rSizeLocal_[2]*r2cSizeLocal_as_,1,FFTW_BACKWARD,FFTW_ESTIMATE);
bPlan_j_real_ = fftwf_plan_many_dft(1,&rSize_[0],rSizeLocal_[2],&temp_[r2cSizeLocal_as_],NULL,rSizeLocal_[2]*r2cSizeLocal_,r2cSizeLocal_,&temp1_[r2cSizeLocal_as_*rSizeLocal_[2]*rSize_[0]],NULL,rSizeLocal_[2],1,FFTW_BACKWARD,FFTW_ESTIMATE);
bPlan_i_ = fftwf_plan_many_dft_c2r(1,&rSize_[0],rSizeLocal_[1] ,temp1_,NULL,1, r2cSize_,rData_,NULL,components_,rJump_[1]*components_,FFTW_ESTIMATE);
//allocation of field
long rfield_size = rfield->lattice().sitesLocalGross();
long kfield_size = kfield->lattice().sitesLocalGross()*2; //*2 for complex type
if(mem_type_==FFT_IN_PLACE)
{
if(rfield_size>kfield_size)
{
rfield->alloc();
kfield->data() = (Imag *)rfield->data();
}
else
{
kfield->alloc();
rfield->data() = (float *)kfield->data();
}
}
if(mem_type_ == FFT_OUT_OF_PLACE)
{
rfield->alloc();
kfield->alloc();
}
//Pointer to data
rData_ = rfield->data();
rData_ += rfield->lattice().siteFirst()*components_;//usefull point !!
cData_ = (fftwf_complex*)kfield->data(); //to be sure that cData is instantiate !
kData_ = (fftwf_complex*)kfield->data();
kData_ += kfield->lattice().siteFirst()*components_;
}
#endif
#ifndef SINGLE
template<class compType>
PlanFFT<compType>::PlanFFT(Field<compType>* rfield,Field<compType>* kfield,const int mem_type)
{
status_ = false;
initialize(rfield,kfield,mem_type);
}
template<class compType>
void PlanFFT<compType>::initialize(Field<compType>* rfield,Field<compType>* kfield,const int mem_type )
{
type_ = C2C;
mem_type_=mem_type;
//general variable
if(rfield->components() != kfield->components())
{
cerr<<"Latfield2d::PlanFFT::initialize : fft curently work only for fields with same number of components"<<endl;
cerr<<"Latfield2d::PlanFFT::initialize : coordinate and Fourier space fields have not the same number of components"<<endl;
cerr<<"Latfield2d : Abort Process Requested"<<endl;
}
else components_ = rfield->components();
for(int i = 0; i<3; i++)
{
rSize_[i]=rfield->lattice().size(i);
kSize_[i]=kfield->lattice().size(i);
rSizeLocal_[i]=rfield->lattice().sizeLocal(i);
kSizeLocal_[i]=kfield->lattice().sizeLocal(i);
rJump_[i]=rfield->lattice().jump(i);
kJump_[i]=kfield->lattice().jump(i);
}
rHalo_ = rfield->lattice().halo();
kHalo_ = kfield->lattice().halo();
tempMemory.setTemp((long)(rSize_[0]) * (long)(rSizeLocal_[1]) * (long)(rSizeLocal_[2]));
temp_ = tempMemory.temp1();
temp1_ = tempMemory.temp2();
if(rfield->lattice().dim()!=3)
{
if(parallel.isRoot())
{
cerr<<"Latfield2d::PlanFFT::initialize : fft curently work only for 3d cubic lattice"<<endl;
cerr<<"Latfield2d::PlanFFT::initialize : real lattice have not 3 dimensions"<<endl;
cerr<<"Latfield2d : Abort Process Requested"<<endl;
}
parallel.abortForce();
}
if(rSize_[0]!=rSize_[1] | rSize_[1]!=rSize_[2])
{
if(parallel.isRoot())
{
cerr<<"Latfield2d::PlanFFT::initialize : fft curently work only for 3d cubic lattice"<<endl;
cerr<<"Latfield2d::PlanFFT::initialize : real lattice is not cubic"<<endl;
cerr<<"Latfield2d : Abort Process Requested"<<endl;
}
parallel.abortForce();
}
//initialization of fftw plan
//Forward plan
fPlan_i_ = fftw_plan_many_dft(1,&rSize_[0],rSizeLocal_[1] ,cData_,NULL,components_, rJump_[1]*components_,temp_,NULL,rSizeLocal_[1]*rSizeLocal_[2],1,FFTW_FORWARD,FFTW_ESTIMATE | FFTW_PRESERVE_INPUT);
fPlan_j_ = fftw_plan_many_dft(1,&rSize_[0],rSizeLocal_[1]*rSizeLocal_[2],temp_,NULL,rSizeLocal_[1]*rSizeLocal_[2],1,temp_,NULL,rSizeLocal_[1]*rSizeLocal_[2],1,FFTW_FORWARD,FFTW_ESTIMATE);
fPlan_k_ = fftw_plan_many_dft(1,&rSize_[0],rSizeLocal_[1] ,temp_,NULL,rSizeLocal_[1]*rSizeLocal_[2],1,kData_,NULL,components_, rJump_[1]*components_,FFTW_FORWARD,FFTW_ESTIMATE);
//Backward plan
bPlan_k_ = fftw_plan_many_dft(1,&rSize_[0],rSizeLocal_[1] ,kData_,NULL,components_, rJump_[1]*components_,temp_,NULL,rSizeLocal_[1]*rSizeLocal_[2],1,FFTW_BACKWARD,FFTW_ESTIMATE | FFTW_PRESERVE_INPUT);
bPlan_j_ = fftw_plan_many_dft(1,&rSize_[0],rSizeLocal_[1]*rSizeLocal_[2],temp_,NULL,rSizeLocal_[1]*rSizeLocal_[2],1,temp_,NULL,rSizeLocal_[1]*rSizeLocal_[2],1,FFTW_BACKWARD,FFTW_ESTIMATE);
bPlan_i_ = fftw_plan_many_dft(1,&rSize_[0],rSizeLocal_[1] ,temp_,NULL,rSizeLocal_[1]*rSizeLocal_[2],1,kData_,NULL,components_, rJump_[1]*components_,FFTW_BACKWARD,FFTW_ESTIMATE);
//allocation of field
long rfield_size = rfield->lattice().sitesLocalGross();
long kfield_size = kfield->lattice().sitesLocalGross();
if(mem_type_ == FFT_IN_PLACE)
{
if(rfield_size>=kfield_size)
{
rfield->alloc();
kfield->data() = (Imag*)rfield->data();
}
else
{
kfield->alloc();
rfield->data() = (Imag*)kfield->data();
}
}
if(mem_type_ == FFT_OUT_OF_PLACE)
{
rfield->alloc();
kfield->alloc();
}
//Pointer to data
rData_ = (double*)rfield->data(); //to be sure that rData is instantiate !
cData_ = (fftw_complex*)rfield->data();
cData_ += rfield->lattice().siteFirst()*components_;
kData_ = (fftw_complex*)kfield->data();
kData_ += kfield->lattice().siteFirst()*components_;
}
template<class compType>
PlanFFT<compType>::PlanFFT(Field<double>* rfield, Field<compType>* kfield,const int mem_type )
{
status_ = false;
initialize(rfield,kfield,mem_type);
}
template<class compType>
void PlanFFT<compType>::initialize(Field<double>* rfield,Field<compType>* kfield,const int mem_type )
{
type_ = R2C;
mem_type_=mem_type;
//general variable
if(rfield->components() != kfield->components())
{
cerr<<"Latfield2d::PlanFFT::initialize : fft curently work only for fields with same number of components"<<endl;
cerr<<"Latfield2d::PlanFFT::initialize : coordinate and Fourier space fields have not the same number of components"<<endl;
cerr<<"Latfield2d : Abort Process Requested"<<endl;
}
else components_ = rfield->components();
for(int i = 0; i<3; i++)
{
rSize_[i]=rfield->lattice().size(i);
kSize_[i]=kfield->lattice().size(i);
rSizeLocal_[i]=rfield->lattice().sizeLocal(i);
kSizeLocal_[i]=kfield->lattice().sizeLocal(i);
rJump_[i]=rfield->lattice().jump(i);
kJump_[i]=kfield->lattice().jump(i);
}
r2cSize_ = rfield->lattice().size(0)/2 + 1;
r2cSizeLocal_as_ = rfield->lattice().sizeLocal(0)/(2*parallel.grid_size()[1]);
if(parallel.last_proc()[1]) r2cSizeLocal_ = r2cSizeLocal_as_ + 1;
else r2cSizeLocal_ = r2cSizeLocal_as_;
rHalo_ = rfield->lattice().halo();
kHalo_ = kfield->lattice().halo();
tempMemory.setTemp((r2cSize_+2) * (rSizeLocal_[1]+2) * (rSizeLocal_[2]+2));
temp_ = tempMemory.temp1();
temp1_ = tempMemory.temp2();
if(rfield->lattice().dim()!=3)
{
if(parallel.isRoot())
{
cerr<<"Latfield2d::PlanFFT::initialize : fft curently work only for 3d cubic lattice"<<endl;
cerr<<"Latfield2d::PlanFFT::initialize : real lattice have not 3 dimensions"<<endl;
cerr<<"Latfield2d : Abort Process Requested"<<endl;
}
parallel.abortForce();
}
//COUT << rSize_[0] << " "<< rSize_[1]<< " "<< rSize_[2] <<endl;
if(rSize_[0]!=rSize_[1] | rSize_[1]!=rSize_[2])
{
if(parallel.isRoot())
{
cerr<<"Latfield2d::PlanFFT::initialize : fft curently work only for 3d cubic lattice"<<endl;
cerr<<"Latfield2d::PlanFFT::initialize : real lattice is not cubic"<<endl;
cerr<<"Latfield2d : Abort Process Requested"<<endl;
}
parallel.abortForce();
}
//create the fftw_plan
fPlan_i_ = fftw_plan_many_dft_r2c(1,&rSize_[0],rSizeLocal_[1] ,rData_,NULL,components_, rJump_[1]*components_,temp_,NULL,rSizeLocal_[1]*rSizeLocal_[2],1,FFTW_ESTIMATE | FFTW_PRESERVE_INPUT);
fPlan_j_ = fftw_plan_many_dft(1,&rSize_[0],rSizeLocal_[2]*r2cSizeLocal_,temp_,NULL,rSizeLocal_[2]*r2cSizeLocal_,1,temp_,NULL,rSizeLocal_[2]*r2cSizeLocal_,1,FFTW_FORWARD,FFTW_ESTIMATE);
fPlan_k_ = fftw_plan_many_dft(1,&rSize_[0],r2cSizeLocal_as_,temp_,NULL,rSizeLocal_[2]*r2cSizeLocal_,1,temp1_,NULL,rSizeLocal_[2]*r2cSizeLocal_as_,1,FFTW_FORWARD,FFTW_ESTIMATE);
fPlan_k_real_ = fftw_plan_many_dft(1,&rSize_[0],rSizeLocal_[2],&temp_[r2cSizeLocal_as_],NULL,rSizeLocal_[2]*r2cSizeLocal_,r2cSizeLocal_,&temp1_[r2cSizeLocal_as_*rSizeLocal_[2]*rSize_[0]],NULL,rSizeLocal_[2],1,FFTW_FORWARD,FFTW_ESTIMATE);
bPlan_k_ = fftw_plan_many_dft(1,&rSize_[0],kSizeLocal_[2]*r2cSizeLocal_,temp_,NULL,kSizeLocal_[2]*r2cSizeLocal_,1,temp_,NULL,kSizeLocal_[2]*r2cSizeLocal_,1,FFTW_BACKWARD,FFTW_ESTIMATE);
bPlan_j_ = fftw_plan_many_dft(1,&rSize_[0],r2cSizeLocal_as_,temp_,NULL,rSizeLocal_[2]*r2cSizeLocal_,1,temp1_,NULL,rSizeLocal_[2]*r2cSizeLocal_as_,1,FFTW_BACKWARD,FFTW_ESTIMATE);
bPlan_j_real_ = fftw_plan_many_dft(1,&rSize_[0],rSizeLocal_[2],&temp_[r2cSizeLocal_as_],NULL,rSizeLocal_[2]*r2cSizeLocal_,r2cSizeLocal_,&temp1_[r2cSizeLocal_as_*rSizeLocal_[2]*rSize_[0]],NULL,rSizeLocal_[2],1,FFTW_BACKWARD,FFTW_ESTIMATE);
bPlan_i_ = fftw_plan_many_dft_c2r(1,&rSize_[0],rSizeLocal_[1] ,temp1_,NULL,1, r2cSize_,rData_,NULL,components_,rJump_[1]*components_,FFTW_ESTIMATE);
//allocation of field
long rfield_size = rfield->lattice().sitesLocalGross();
long kfield_size = kfield->lattice().sitesLocalGross()*2; //*2 for complex type
if(mem_type_==FFT_IN_PLACE)
{
if(rfield_size>kfield_size)
{
rfield->alloc();
kfield->data() = (Imag *)rfield->data();
}
else
{
kfield->alloc();
rfield->data() = (double *)kfield->data();
}
}
if(mem_type_ == FFT_OUT_OF_PLACE)
{
rfield->alloc();
kfield->alloc();
}
//Pointer to data
rData_ = rfield->data();
rData_ += rfield->lattice().siteFirst()*components_;//usefull point !!
cData_ = (fftw_complex*)kfield->data(); //to be sure that cData is instantiate !
kData_ = (fftw_complex*)kfield->data();
kData_ += kfield->lattice().siteFirst()*components_;
}
#endif
template<class compType>
void PlanFFT<compType>::execute(int fft_type)
{
temp_ = tempMemory.temp1();
temp1_ = tempMemory.temp2();
//#ifdef SINGLE
if(type_ == R2C)
{
if(fft_type == FFT_FORWARD)
{
int i,j,k;
int comp;
int comm_rank;
#ifdef SINGLE
float *p_in;
fftwf_complex *p_out;
#else
double *p_in;
fftw_complex *p_out;
#endif
for(comp=0;comp<components_;comp++)
{
//execute first dimension fft, prepar rData in temp_ to be send via AlltoAll + gather
#ifdef SINGLE
for(int l = 0;l< rSizeLocal_[2] ;l++)
{
p_in = &rData_[rJump_[2]*l*components_ + comp];
p_out = &temp_[l*rSizeLocal_[1]];
fftwf_execute_dft_r2c(fPlan_i_,p_in,p_out);
}
#else
for(int l = 0;l< rSizeLocal_[2] ;l++)
{
p_in = &rData_[rJump_[2]*l*components_ + comp];
p_out = &temp_[l*rSizeLocal_[1]];
fftw_execute_dft_r2c(fPlan_i_,p_in,p_out);
}
#endif
/*
//verif step 1
for(int proc=0; proc<parallel.size(); proc++)
{
if(proc==parallel.rank())
{
for(k=0 ; k < rSizeLocal_[2] ;k++)
{
for(j=0 ; j < rSizeLocal_[1] ;j++)
{
for(i=0 ; i < (r2cSize_) ;i++)
{
//cout << "("<<parallel.grid_rank()[0]<<","<< parallel.grid_rank()[1]<< "),("<< i <<","<< j + rSizeLocal_[1]*parallel.grid_rank()[1] <<","<< k+rSizeLocal_[2]*parallel.grid_rank()[0]<<")"<< temp_[j+rSizeLocal_[1]*(k+rSizeLocal_[2] *i)][0]<<" , "<<temp_[j+rSizeLocal_[1]*(k+rSizeLocal_[2] *i)][1]<<endl;
//temp_[j+rSizeLocal_[1]*(k+rSizeLocal_[2] *i)][0] = k + rSizeLocal_[2]*parallel.grid_rank()[0];
//temp_[j+rSizeLocal_[1]*(k+rSizeLocal_[2] *i)][1] = 0;//j + rSizeLocal_[1]*parallel.grid_rank()[1];
//cout <<rSizeLocal_[1]<<" "<< parallel.grid_rank()[1] <<"("<< i <<","<< j + rSizeLocal_[1]*parallel.grid_rank()[1] <<","<< k+rSizeLocal_[2]*parallel.grid_rank()[0]<<")"<< temp_[j+rSizeLocal_[1]*(k+rSizeLocal_[2] *i)][0]<<" , "<<temp_[j+rSizeLocal_[1]*(k+rSizeLocal_[2] *i)][1]<<endl;
}
}
}
}
MPI_Barrier(parallel.lat_world_comm());
} */
MPI_Alltoall(temp_, 2* rSizeLocal_[1]*rSizeLocal_[2]*r2cSizeLocal_as_, MPI_DATA_PREC, temp1_, 2* rSizeLocal_[1]*rSizeLocal_[2]*r2cSizeLocal_as_, MPI_DATA_PREC, parallel.dim1_comm()[parallel.grid_rank()[0]]);
MPI_Gather(&temp_[(r2cSize_-1)*rSizeLocal_[1]*rSizeLocal_[2]][0], 2*rSizeLocal_[1]*rSizeLocal_[2], MPI_DATA_PREC, &temp1_[rSize_[0]/2*rSizeLocal_[1]*rSizeLocal_[2]][0] , 2*rSizeLocal_[1]*rSizeLocal_[2], MPI_DATA_PREC ,parallel.grid_size()[1]-1, parallel.dim1_comm()[parallel.grid_rank()[0]]);
MPI_Barrier(parallel.dim1_comm()[parallel.grid_rank()[0]]);
if(parallel.last_proc()[1])
{
for(i=0;i<parallel.grid_size()[1];i++)transpose_0_2_last_proc(&temp1_[i*rSizeLocal_[1]*rSizeLocal_[2]*r2cSizeLocal_as_],&temp_[i*rSizeLocal_[1]*rSizeLocal_[2]*r2cSizeLocal_],rSizeLocal_[1],rSizeLocal_[2],r2cSizeLocal_as_);
implement_local_0_last_proc(&temp1_[rSize_[0]/2*rSizeLocal_[1]*rSizeLocal_[2]],temp_,rSizeLocal_[1],rSizeLocal_[2],r2cSizeLocal_as_,parallel.grid_size()[1]);
}
else for(i=0;i<parallel.grid_size()[1];i++)transpose_0_2(&temp1_[i*rSizeLocal_[1]*rSizeLocal_[2]*r2cSizeLocal_as_],&temp_[i*rSizeLocal_[1]*rSizeLocal_[2]*r2cSizeLocal_as_],rSizeLocal_[1],rSizeLocal_[2],r2cSizeLocal_as_);
#ifdef SINGLE
fftwf_execute(fPlan_j_);
#else
fftw_execute(fPlan_j_);
#endif
//verif step 2
/*
for(k=0 ; k < rSizeLocal_[2] ;k++)
{
for(j=0 ; j < rSize_[1] ;j++)
{
for(i=0 ; i < (r2cSizeLocal_) ;i++)
{
//cout << "("<< i + r2cSizeLocal_as_*parallel.grid_rank()[1] <<","<< j <<","<< k+rSizeLocal_[2]*parallel.grid_rank()[0]<<")"<< temp_[i+r2cSizeLocal_*(k+rSizeLocal_[2]*j)][0]<<" , "<<temp_[i+r2cSizeLocal_*(k+rSizeLocal_[2]*j)][1]<<endl;
temp_[i+r2cSizeLocal_*(k+rSizeLocal_[2]*j)][1] = k + rSizeLocal_[2]*parallel.grid_rank()[0] ;
temp_[i+r2cSizeLocal_*(k+rSizeLocal_[2]*j)][0] = i + r2cSizeLocal_as_*parallel.grid_rank()[1];
}
}
}*/
MPI_Barrier(parallel.lat_world_comm());
MPI_Alltoall(temp_, (2*rSizeLocal_[2]*rSizeLocal_[2]*r2cSizeLocal_), MPI_DATA_PREC, temp1_, (2*rSizeLocal_[2]*rSizeLocal_[2]*r2cSizeLocal_), MPI_DATA_PREC, parallel.dim0_comm()[parallel.grid_rank()[1]]);
MPI_Barrier(parallel.dim0_comm()[parallel.grid_rank()[1]]);
for(i=0;i<parallel.grid_size()[0];i++)transpose_1_2(&temp1_[i*rSizeLocal_[2]*rSizeLocal_[2]*r2cSizeLocal_],&temp_[i*rSizeLocal_[2]*rSizeLocal_[2]*r2cSizeLocal_], r2cSizeLocal_,rSizeLocal_[2],rSizeLocal_[2]);
/*
//verif transpos
for(int proc=0; proc<parallel.size(); proc++)
{
if(proc==parallel.rank())
{
for(k=0 ; k < rSize_[2] ;k++)
{
for(j=0 ; j < rSizeLocal_[2] ;j++)
{
for(i=0 ; i < r2cSizeLocal_ ;i++)
{
cout << "("<<parallel.grid_rank()[0]<<","<< parallel.grid_rank()[1]<< "),("<< i + r2cSizeLocal_as_*parallel.grid_rank()[1] <<","<< j+rSizeLocal_[2]*parallel.grid_rank()[0] <<","<< k<<")"<< temp_[i+r2cSizeLocal_*(j+rSizeLocal_[2]*k)][0]<<" , "<<temp_[i+r2cSizeLocal_*(j+rSizeLocal_[2]*k)][1]<<endl;
temp_[i+r2cSizeLocal_*(j+rSizeLocal_[2]*k)][0] = 1;
temp_[i+r2cSizeLocal_*(j+rSizeLocal_[2]*k)][1] = 0;
}
}
}
}
MPI_Barrier(parallel.lat_world_comm());
}*/
#ifdef SINGLE
for(int l=0;l<rSizeLocal_[2];l++)
{
fftwf_execute_dft(fPlan_k_,&temp_[l*r2cSizeLocal_],&temp1_[l*r2cSizeLocal_as_]);
}
if(parallel.last_proc()[1])fftwf_execute(fPlan_k_real_);
#else
for(int l=0;l<rSizeLocal_[2];l++)
{
fftw_execute_dft(fPlan_k_,&temp_[l*r2cSizeLocal_],&temp1_[l*r2cSizeLocal_as_]);
}
if(parallel.last_proc()[1])fftw_execute(fPlan_k_real_);
#endif