forked from dealias/fftwpp
-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathconvolution.h
2017 lines (1730 loc) · 57.8 KB
/
convolution.h
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
/* Implicitly dealiased convolution routines.
Copyright (C) 2010-2015 John C. Bowman and Malcolm Roberts, Univ. of Alberta
This program is free software; you can redistribute it and/or modify
it under the terms of the GNU Lesser General Public License as published by
the Free Software Foundation; either version 3 of the License, or
(at your option) any later version.
This program is distributed in the hope that it will be useful,
but WITHOUT ANY WARRANTY; without even the implied warranty of
MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
GNU Lesser General Public License for more details.
You should have received a copy of the GNU Lesser General Public License
along with this program; if not, write to the Free Software
Foundation, Inc., 675 Mass Ave, Cambridge, MA 02139, USA. */
#include "Complex.h"
#include "fftw++.h"
#include "cmult-sse2.h"
#include "transposeoptions.h"
namespace fftwpp {
#ifndef __convolution_h__
#define __convolution_h__ 1
extern const double sqrt3;
extern const double hsqrt3;
extern const Complex hSqrt3;
extern const Complex mhsqrt3;
extern const Complex mhalf;
extern const Complex zeta3;
extern const double twopi;
inline unsigned int min(unsigned int a, unsigned int b)
{
return (a < b) ? a : b;
}
inline unsigned int max(unsigned int a, unsigned int b)
{
return (a > b) ? a : b;
}
// Build the factored zeta tables.
unsigned int BuildZeta(double arg, unsigned int m,
Complex *&ZetaH, Complex *&ZetaL,
unsigned int threads=1);
unsigned int BuildZeta(unsigned int n, unsigned int m,
Complex *&ZetaH, Complex *&ZetaL,
unsigned int threads=1);
struct convolveOptions {
unsigned int nx,ny,nz; // |
unsigned int stride2,stride3; // | Used internally by the MPI interface.
utils::mpiOptions mpi; // |
bool toplevel;
convolveOptions(unsigned int nx, unsigned int ny, unsigned int nz,
unsigned int stride2, unsigned int stride3) :
nx(nx), ny(ny), nz(nz), stride2(stride2), stride3(stride3),
toplevel(true) {}
convolveOptions(unsigned int nx, unsigned int ny, unsigned int stride2,
utils::mpiOptions mpi, bool toplevel=true) :
nx(nx), ny(ny), stride2(stride2), mpi(mpi), toplevel(toplevel) {}
convolveOptions(unsigned int ny, unsigned int nz,
unsigned int stride2, unsigned int stride3,
utils::mpiOptions mpi, bool toplevel=true) :
ny(ny), nz(nz), stride2(stride2), stride3(stride3), mpi(mpi),
toplevel(toplevel) {}
convolveOptions(bool toplevel=true) : nx(0), ny(0), nz(0),
toplevel(toplevel) {}
};
static const convolveOptions defaultconvolveOptions;
typedef void multiplier(Complex **, unsigned int m,
const unsigned int indexsize,
const unsigned int *index,
unsigned int r, unsigned int threads);
typedef void realmultiplier(double **, unsigned int m,
const unsigned int indexsize,
const unsigned int *index,
unsigned int r, unsigned int threads);
// Multipliers for binary convolutions.
multiplier multautoconvolution;
multiplier multautocorrelation;
multiplier multbinary;
multiplier multcorrelation;
multiplier multbinary2;
multiplier multbinary3;
multiplier multbinary4;
multiplier multbinary8;
realmultiplier multbinary;
realmultiplier multbinary2;
realmultiplier multadvection2;
struct general {};
struct pretransform1 {};
struct pretransform2 {};
struct pretransform3 {};
struct pretransform4 {};
// In-place implicitly dealiased 1D complex convolution using
// function pointers for multiplication
class ImplicitConvolution : public ThreadBase {
private:
unsigned int m;
Complex **U;
unsigned int A;
unsigned int B;
Complex *u;
unsigned int s;
Complex *ZetaH, *ZetaL;
fft1d *BackwardsO,*ForwardsO;
fft1d *Backwards,*Forwards;
bool pointers;
bool allocated;
unsigned int indexsize;
public:
unsigned int *index;
void initpointers(Complex **&U, Complex *u) {
unsigned int C=max(A,B);
U=new Complex *[C];
for(unsigned int a=0; a < C; ++a)
U[a]=u+a*m;
pointers=true;
}
void deletepointers(Complex **&U) {
delete [] U;
}
void allocateindex(unsigned int n, unsigned int *i) {
indexsize=n;
index=i;
}
void init() {
indexsize=0;
Complex* U0=U[0];
Complex* U1=A == 1 ? utils::ComplexAlign(m) : U[1];
BackwardsO=new fft1d(m,1,U0,U1);
ForwardsO=new fft1d(m,-1,U0,U1);
threads=std::min(threads,max(BackwardsO->Threads(),ForwardsO->Threads()));
if(A == B) {
Backwards=new fft1d(m,1,U0);
threads=std::min(threads,Backwards->Threads());
}
if(A <= B) {
Forwards=new fft1d(m,-1,U0);
threads=std::min(threads,Forwards->Threads());
}
if(A == 1) utils::deleteAlign(U1);
s=BuildZeta(2*m,m,ZetaH,ZetaL,threads);
}
// m is the number of Complex data values.
// U is an array of C distinct work arrays each of size m, where C=max(A,B)
// A is the number of inputs.
// B is the number of outputs.
ImplicitConvolution(unsigned int m, Complex **U, unsigned int A=2,
unsigned int B=1,
unsigned int threads=fftw::maxthreads)
: ThreadBase(threads), m(m), U(U), A(A), B(B), pointers(false),
allocated(false) {
init();
}
// m is the number of Complex data values.
// u is a work array of C*m Complex values.
// A is the number of inputs.
// B is the number of outputs.
ImplicitConvolution(unsigned int m, Complex *u,
unsigned int A=2, unsigned int B=1,
unsigned int threads=fftw::maxthreads)
: ThreadBase(threads), m(m), A(A), B(B), u(u), allocated(false) {
initpointers(U,u);
init();
}
// m is the number of Complex data values.
// A is the number of inputs.
// B is the number of outputs.
ImplicitConvolution(unsigned int m,
unsigned int A=2, unsigned int B=1,
unsigned int threads=fftw::maxthreads)
: ThreadBase(threads), m(m), A(A), B(B), allocated(true) {
u=utils::ComplexAlign(max(A,B)*m);
initpointers(U,u);
init();
}
~ImplicitConvolution() {
utils::deleteAlign(ZetaH);
utils::deleteAlign(ZetaL);
if(pointers) deletepointers(U);
if(allocated) utils::deleteAlign(u);
if(A == B)
delete Backwards;
if(A <= B)
delete Forwards;
delete ForwardsO;
delete BackwardsO;
}
// F is an array of A pointers to distinct data blocks each of size m,
// shifted by offset (contents not preserved).
void convolve(Complex **F, multiplier *pmult, unsigned int i=0,
unsigned int offset=0);
void autoconvolve(Complex *f) {
Complex *F[]={f};
convolve(F,multautoconvolution);
}
void autocorrelate(Complex *f) {
Complex *F[]={f};
convolve(F,multautocorrelation);
}
// Binary convolution:
void convolve(Complex *f, Complex *g) {
Complex *F[]={f,g};
convolve(F,multbinary);
}
// Binary correlation:
void correlate(Complex *f, Complex *g) {
Complex *F[]={f,g};
convolve(F,multcorrelation);
}
template<class T>
inline void pretransform(Complex **F, unsigned int k, Vec& Zetak);
template<class T>
void pretransform(Complex **F);
void posttransform(Complex *f, Complex *u);
};
// In-place implicitly dealiased 1D Hermitian convolution.
class ImplicitHConvolution : public ThreadBase {
protected:
unsigned int m;
unsigned int c;
bool compact;
Complex **U;
unsigned int A;
unsigned int B;
Complex *u;
unsigned int s;
Complex *ZetaH,*ZetaL;
rcfft1d *rc,*rco,*rcO;
crfft1d *cr,*cro,*crO;
Complex *w; // Work array of size max(A,B) to hold f[c] in even case.
bool pointers;
bool allocated;
bool even;
unsigned int indexsize;
public:
unsigned int *index;
void initpointers(Complex **&U, Complex *u) {
unsigned int C=max(A,B);
U=new Complex *[C];
unsigned stride=c+1;
for(unsigned int a=0; a < C; ++a)
U[a]=u+a*stride;
pointers=true;
}
void deletepointers(Complex **&U) {
delete [] U;
}
void allocateindex(unsigned int n, unsigned int *i) {
indexsize=n;
index=i;
}
void init() {
even=m == 2*c;
indexsize=0;
Complex* U0=U[0];
rc=new rcfft1d(m,U0);
cr=new crfft1d(m,U0);
Complex* U1=A == 1 ? utils::ComplexAlign(m) : U[1];
rco=new rcfft1d(m,(double *) U0,U1);
cro=new crfft1d(m,U1,(double *) U0);
if(A == 1) utils::deleteAlign(U1);
if(A != B) {
rcO=rco;
crO=cro;
} else {
rcO=rc;
crO=cr;
}
threads=std::min(threads,std::max(rco->Threads(),cro->Threads()));
s=BuildZeta(3*m,c+2,ZetaH,ZetaL,threads);
w=even ? utils::ComplexAlign(max(A,B)) : u;
}
// m is the number of independent data values
// U is an array of max(A,B) distinct work arrays of size c+1, where c=m/2
// A is the number of inputs.
// B is the number of outputs.
ImplicitHConvolution(unsigned int m, Complex **U, unsigned int A=2,
unsigned int B=1, unsigned int threads=fftw::maxthreads)
: ThreadBase(threads), m(m), c(m/2), compact(true), U(U), A(A), B(B),
pointers(false), allocated(false) {
init();
}
ImplicitHConvolution(unsigned int m, bool compact, Complex **U,
unsigned int A=2, unsigned int B=1,
unsigned int threads=fftw::maxthreads)
: ThreadBase(threads), m(m), c(m/2), compact(compact), U(U), A(A), B(B),
pointers(false), allocated(false) {
init();
}
// m is the number of independent data values
// u is a work array of max(A,B)*(c+1) Complex values, where c=m/2
// A is the number of inputs.
// B is the number of outputs.
ImplicitHConvolution(unsigned int m, Complex *u,
unsigned int A=2, unsigned int B=1,
unsigned int threads=fftw::maxthreads)
: ThreadBase(threads), m(m), c(m/2), compact(true), A(A), B(B), u(u),
allocated(false) {
initpointers(U,u);
init();
}
ImplicitHConvolution(unsigned int m, bool compact, Complex *u,
unsigned int A=2, unsigned int B=1,
unsigned int threads=fftw::maxthreads)
: ThreadBase(threads), m(m), c(m/2), compact(compact), A(A), B(B), u(u),
allocated(false) {
initpointers(U,u);
init();
}
// m is the number of independent data values
// u is a work array of max(A,B)*(c+1) Complex values, where c=m/2
// A is the number of inputs.
// B is the number of outputs.
ImplicitHConvolution(unsigned int m, bool compact=true, unsigned int A=2,
unsigned int B=1, unsigned int threads=fftw::maxthreads)
: ThreadBase(threads), m(m), c(m/2), compact(compact), A(A), B(B),
u(utils::ComplexAlign(max(A,B)*(c+1))), allocated(true) {
initpointers(U,u);
init();
}
virtual ~ImplicitHConvolution() {
if(even) utils::deleteAlign(w);
utils::deleteAlign(ZetaH);
utils::deleteAlign(ZetaL);
if(pointers) deletepointers(U);
if(allocated) utils::deleteAlign(u);
if(A != B) {
delete cro;
delete rco;
}
delete cr;
delete rc;
}
// F is an array of A pointers to distinct data blocks each of size m,
// shifted by offset (contents not preserved).
void convolve(Complex **F, realmultiplier *pmult, unsigned int i=0,
unsigned int offset=0);
void pretransform(Complex *F, Complex *f1c, Complex *U);
void posttransform(Complex *F, const Complex& f1c, Complex *U);
// Binary convolution:
void convolve(Complex *f, Complex *g) {
Complex *F[]={f,g};
convolve(F,multbinary);
}
};
// Compute the scrambled implicitly m-padded complex Fourier transform of M
// complex vectors, each of length m.
// The arrays in and out (which may coincide), along with the array u, must
// be allocated as Complex[M*m].
//
// fftpad fft(m,M,stride);
// fft.backwards(in,u);
// fft.forwards(in,u);
//
// Notes:
// stride is the spacing between the elements of each Complex vector.
//
class fftpad {
unsigned int m;
unsigned int M;
unsigned int stride;
unsigned int dist;
unsigned int s;
Complex *ZetaH, *ZetaL;
unsigned int threads;
public:
mfft1d *Backwards;
mfft1d *Forwards;
fftpad(unsigned int m, unsigned int M,
unsigned int stride, Complex *u=NULL,
unsigned int Threads=fftw::maxthreads)
: m(m), M(M), stride(stride), threads(Threads) {
Backwards=new mfft1d(m,1,M,stride,1,u,NULL,threads);
Forwards=new mfft1d(m,-1,M,stride,1,u,NULL,threads);
threads=std::max(Backwards->Threads(),Forwards->Threads());
s=BuildZeta(2*m,m,ZetaH,ZetaL,threads);
}
~fftpad() {
utils::deleteAlign(ZetaL);
utils::deleteAlign(ZetaH);
delete Forwards;
delete Backwards;
}
void expand(Complex *f, Complex *u);
void reduce(Complex *f, Complex *u);
void backwards(Complex *f, Complex *u);
void forwards(Complex *f, Complex *u);
};
// Compute the scrambled implicitly m-padded complex Fourier transform of M
// complex vectors, each of length 2m-1 with the origin at index m-1,
// containing physical data for wavenumbers -m+1 to m-1.
// The arrays in and out (which may coincide) must be allocated as
// Complex[M*(2m-1)]. The array u must be allocated as Complex[M*(m+1)].
//
// fft0pad fft(m,M,stride,u);
// fft.backwards(in,u);
// fft.forwards(in,u);
//
// Notes:
// stride is the spacing between the elements of each Complex vector.
//
class fft0pad {
protected:
unsigned int m;
unsigned int M;
unsigned int s;
unsigned int stride;
Complex *ZetaH, *ZetaL;
unsigned int threads;
public:
mfft1d *Forwards;
mfft1d *Backwards;
fft0pad(unsigned int m, unsigned int M, unsigned int stride, Complex *u=NULL,
unsigned int Threads=fftw::maxthreads)
: m(m), M(M), stride(stride), threads(Threads) {
Backwards=new mfft1d(m,1,M,stride,1,u,NULL,threads);
Forwards=new mfft1d(m,-1,M,stride,1,u,NULL,threads);
s=BuildZeta(3*m,m,ZetaH,ZetaL);
}
virtual ~fft0pad() {
utils::deleteAlign(ZetaL);
utils::deleteAlign(ZetaH);
delete Forwards;
delete Backwards;
}
// Unscramble indices, returning spatial index stored at position i
inline static unsigned findex(unsigned i, unsigned int m) {
return i < m-1 ? 3*i : 3*i+4-3*m; // for i >= m-1: j=3*(i-(m-1))+1
}
inline static unsigned uindex(unsigned i, unsigned int m) {
return i > 0 ? (i < m ? 3*i-1 : 3*m-3) : 3*m-1;
}
virtual void expand(Complex *f, Complex *u);
virtual void reduce(Complex *f, Complex *u);
void backwards(Complex *f, Complex *u);
virtual void forwards(Complex *f, Complex *u);
virtual void Backwards1(Complex *f, Complex *u);
virtual void Forwards0(Complex *f);
virtual void Forwards1(Complex *f, Complex *u);
};
// Compute the scrambled implicitly m-padded complex Fourier transform of M
// complex vectors, each of length 2m with the origin at index m,
// corresponding to wavenumbers -m to m-1.
// The arrays in and out (which may coincide) must be allocated as
// Complex[M*2m]. The array u must be allocated as Complex[M*m].
//
// fft1pad fft(m,M,stride,u);
// fft.backwards(in,u);
// fft.forwards(in,u);
//
// Notes:
// stride is the spacing between the elements of each Complex vector.
//
class fft1pad : public fft0pad {
public:
fft1pad(unsigned int m, unsigned int M, unsigned int stride,
Complex *u=NULL, unsigned int threads=fftw::maxthreads) :
fft0pad(m,M,stride,u,threads) {}
// Unscramble indices, returning spatial index stored at position i
inline static unsigned findex(unsigned i, unsigned int m) {
return i < m ? 3*i : 3*(i-m)+1;
}
inline static unsigned uindex(unsigned i, unsigned int m) {
return i > 0 ? 3*i-1 : 3*m-1;
}
void expand(Complex *f, Complex *u);
void reduce(Complex *f, Complex *u);
void forwards(Complex *f, Complex *u);
void Backwards1(Complex *f, Complex *u);
void Forwards0(Complex *f);
void Forwards1(Complex *f, Complex *u);
};
// In-place implicitly dealiased 2D complex convolution.
class ImplicitConvolution2 : public ThreadBase {
protected:
unsigned int mx,my;
Complex *u1;
Complex *u2;
unsigned int A,B;
fftpad *xfftpad;
ImplicitConvolution **yconvolve;
Complex **U2;
bool allocated;
unsigned int indexsize;
bool toplevel;
public:
unsigned int *index;
void initpointers2(Complex **&U2, Complex *u2, unsigned int stride) {
U2=new Complex *[A];
for(unsigned int a=0; a < A; ++a)
U2[a]=u2+a*stride;
if(toplevel) allocateindex(1,new unsigned int[1]);
}
void deletepointers2(Complex **&U2) {
if(toplevel) {
delete [] index;
for(unsigned int t=1; t < threads; ++t)
delete [] yconvolve[t]->index;
}
delete [] U2;
}
void allocateindex(unsigned int n, unsigned int *i) {
indexsize=n;
index=i;
yconvolve[0]->allocateindex(n,i);
for(unsigned int t=1; t < threads; ++t)
yconvolve[t]->allocateindex(n,new unsigned int[n]);
}
void init(const convolveOptions& options) {
toplevel=options.toplevel;
xfftpad=new fftpad(mx,options.ny,options.ny,u2,threads);
unsigned int C=max(A,B);
yconvolve=new ImplicitConvolution*[threads];
for(unsigned int t=0; t < threads; ++t)
yconvolve[t]=new ImplicitConvolution(my,u1+t*my*C,A,B,innerthreads);
initpointers2(U2,u2,options.stride2);
}
void set(convolveOptions& options) {
if(options.nx == 0) options.nx=mx;
if(options.ny == 0) {
options.ny=my;
options.stride2=mx*my;
}
}
// u1 is a temporary array of size my*C*threads.
// u2 is a temporary array of size mx*my*C.
// A is the number of inputs.
// B is the number of outputs.
// Here C=max(A,B).
ImplicitConvolution2(unsigned int mx, unsigned int my,
Complex *u1, Complex *u2,
unsigned int A=2, unsigned int B=1,
unsigned int threads=fftw::maxthreads,
convolveOptions options=defaultconvolveOptions) :
ThreadBase(threads), mx(mx), my(my), u1(u1), u2(u2), A(A), B(B),
allocated(false) {
set(options);
multithread(options.nx);
init(options);
}
ImplicitConvolution2(unsigned int mx, unsigned int my,
unsigned int A=2, unsigned int B=1,
unsigned int threads=fftw::maxthreads,
convolveOptions options=defaultconvolveOptions) :
ThreadBase(threads), mx(mx), my(my), A(A), B(B), allocated(true) {
set(options);
multithread(options.nx);
unsigned int C=max(A,B);
u1=utils::ComplexAlign(my*C*threads);
u2=utils::ComplexAlign(options.stride2*C);
init(options);
}
virtual ~ImplicitConvolution2() {
deletepointers2(U2);
for(unsigned int t=0; t < threads; ++t)
delete yconvolve[t];
delete [] yconvolve;
delete xfftpad;
if(allocated) {
utils::deleteAlign(u2);
utils::deleteAlign(u1);
}
}
void backwards(Complex **F, Complex **U2, unsigned int offset) {
for(unsigned int a=0; a < A; ++a)
xfftpad->backwards(F[a]+offset,U2[a]);
}
void subconvolution(Complex **F, multiplier *pmult,
unsigned int r, unsigned int M, unsigned int stride,
unsigned int offset=0) {
if(threads > 1) {
#ifndef FFTWPP_SINGLE_THREAD
#pragma omp parallel for num_threads(threads)
#endif
for(unsigned int i=0; i < M; ++i)
yconvolve[get_thread_num()]->convolve(F,pmult,2*i+r,offset+i*stride);
} else {
ImplicitConvolution *yconvolve0=yconvolve[0];
for(unsigned int i=0; i < M; ++i)
yconvolve0->convolve(F,pmult,2*i+r,offset+i*stride);
}
}
void forwards(Complex **F, Complex **U2, unsigned int offset) {
for(unsigned int b=0; b < B; ++b)
xfftpad->forwards(F[b]+offset,U2[b]);
}
// F is a pointer to A distinct data blocks each of size mx*my,
// shifted by offset (contents not preserved).
virtual void convolve(Complex **F, multiplier *pmult, unsigned int i=0,
unsigned int offset=0) {
if(!toplevel) {
index[indexsize-2]=i;
if(threads > 1) {
for(unsigned int t=1; t < threads; ++t) {
unsigned int *Index=yconvolve[t]->index;
for(unsigned int i=0; i < indexsize; ++i)
Index[i]=index[i];
}
}
}
backwards(F,U2,offset);
subconvolution(F,pmult,0,mx,my,offset);
subconvolution(U2,pmult,1,mx,my);
forwards(F,U2,offset);
}
// Binary convolution:
void convolve(Complex *f, Complex *g) {
Complex *F[]={f,g};
convolve(F,multbinary);
}
// Binary correlation:
void correlate(Complex *f, Complex *g) {
Complex *F[]={f, g};
convolve(F,multcorrelation);
}
void autoconvolve(Complex *f) {
Complex *F[]={f};
convolve(F,multautoconvolution);
}
void autocorrelate(Complex *f) {
Complex *F[]={f};
convolve(F,multautocorrelation);
}
};
inline void HermitianSymmetrizeX(unsigned int mx, unsigned int my,
unsigned int xorigin, Complex *f)
{
unsigned int offset=xorigin*my;
unsigned int stop=mx*my;
f[offset].im=0.0;
for(unsigned int i=my; i < stop; i += my)
f[offset-i]=conj(f[offset+i]);
}
// Enforce 3D Hermiticity using specified (x,y > 0,z=0) and (x >= 0,y=0,z=0)
// data.
inline void HermitianSymmetrizeXY(unsigned int mx, unsigned int my,
unsigned int mz, unsigned int xorigin,
unsigned int yorigin, Complex *f,
unsigned int threads=fftw::maxthreads)
{
int stride=(yorigin+my)*mz;
int mxstride=mx*stride;
unsigned int myz=my*mz;
unsigned int origin=xorigin*stride+yorigin*mz;
f[origin].im=0.0;
for(int i=stride; i < mxstride; i += stride)
f[origin-i]=conj(f[origin+i]);
PARALLEL(
for(int i=stride-mxstride; i < mxstride; i += stride) {
int stop=i+myz;
for(int j=i+mz; j < stop; j += mz) {
f[origin-j]=conj(f[origin+j]);
}
}
);
}
typedef unsigned int IndexFunction(unsigned int, unsigned int m);
class ImplicitHConvolution2 : public ThreadBase {
protected:
unsigned int mx,my;
bool xcompact,ycompact;
Complex *u1;
Complex *u2;
unsigned int A,B;
fft0pad *xfftpad;
ImplicitHConvolution **yconvolve;
Complex **U2;
bool allocated;
unsigned int indexsize;
bool toplevel;
public:
unsigned int *index;
void initpointers2(Complex **&U2, Complex *u2, unsigned int stride)
{
unsigned int C=max(A,B);
U2=new Complex *[C];
for(unsigned int a=0; a < C; ++a)
U2[a]=u2+a*stride;
if(toplevel) allocateindex(1,new unsigned int[1]);
}
void deletepointers2(Complex **&U2) {
if(toplevel) {
delete [] index;
for(unsigned int t=1; t < threads; ++t)
delete [] yconvolve[t]->index;
}
delete [] U2;
}
void allocateindex(unsigned int n, unsigned int *i) {
indexsize=n;
index=i;
yconvolve[0]->allocateindex(n,i);
for(unsigned int t=1; t < threads; ++t)
yconvolve[t]->allocateindex(n,new unsigned int[n]);
}
void init(const convolveOptions& options) {
unsigned int C=max(A,B);
toplevel=options.toplevel;
xfftpad=xcompact ? new fft0pad(mx,options.ny,options.ny,u2) :
new fft1pad(mx,options.ny,options.ny,u2);
yconvolve=new ImplicitHConvolution*[threads];
for(unsigned int t=0; t < threads; ++t)
yconvolve[t]=new ImplicitHConvolution(my,ycompact,u1+t*(my/2+1)*C,A,B,
innerthreads);
initpointers2(U2,u2,options.stride2);
}
void set(convolveOptions& options) {
if(options.nx == 0) options.nx=mx;
if(options.ny == 0) {
options.ny=my+!ycompact;
options.stride2=(mx+xcompact)*options.ny;
}
}
// u1 is a temporary array of size (my/2+1)*C*threads.
// u2 is a temporary array of size (mx+xcompact)*(my+!ycompact)*C;
// A is the number of inputs.
// B is the number of outputs.
// Here C=max(A,B).
ImplicitHConvolution2(unsigned int mx, unsigned int my,
Complex *u1, Complex *u2,
unsigned int A=2, unsigned int B=1,
unsigned int threads=fftw::maxthreads,
convolveOptions options=defaultconvolveOptions) :
ThreadBase(threads), mx(mx), my(my), xcompact(true), ycompact(true),
u1(u1), u2(u2), A(A), B(B), allocated(false) {
set(options);
multithread(options.nx);
init(options);
}
ImplicitHConvolution2(unsigned int mx, unsigned int my,
bool xcompact, bool ycompact,
Complex *u1, Complex *u2,
unsigned int A=2, unsigned int B=1,
unsigned int threads=fftw::maxthreads,
convolveOptions options=defaultconvolveOptions) :
ThreadBase(threads), mx(mx), my(my),
xcompact(xcompact), ycompact(ycompact), u1(u1), u2(u2), A(A), B(B),
allocated(false) {
set(options);
multithread(options.nx);
init(options);
}
ImplicitHConvolution2(unsigned int mx, unsigned int my,
bool xcompact=true, bool ycompact=true,
unsigned int A=2, unsigned int B=1,
unsigned int threads=fftw::maxthreads,
convolveOptions options=defaultconvolveOptions) :
ThreadBase(threads), mx(mx), my(my),
xcompact(xcompact), ycompact(ycompact), A(A), B(B), allocated(true) {
set(options);
multithread(options.nx);
unsigned int C=max(A,B);
u1=utils::ComplexAlign((my/2+1)*C*threads);
u2=utils::ComplexAlign(options.stride2*C);
init(options);
}
virtual ~ImplicitHConvolution2() {
deletepointers2(U2);
for(unsigned int t=0; t < threads; ++t)
delete yconvolve[t];
delete [] yconvolve;
delete xfftpad;
if(allocated) {
utils::deleteAlign(u2);
utils::deleteAlign(u1);
}
}
void backwards(Complex **F, Complex **U2, unsigned int ny,
bool symmetrize, unsigned int offset) {
for(unsigned int a=0; a < A; ++a) {
Complex *f=F[a]+offset;
if(symmetrize)
HermitianSymmetrizeX(mx,ny,mx-xcompact,f);
xfftpad->backwards(f,U2[a]);
}
}
void subconvolution(Complex **F, realmultiplier *pmult,
IndexFunction indexfunction,
unsigned int M, unsigned int stride,
unsigned int offset=0) {
if(threads > 1) {
#ifndef FFTWPP_SINGLE_THREAD
#pragma omp parallel for num_threads(threads)
#endif
for(unsigned int i=0; i < M; ++i)
yconvolve[get_thread_num()]->convolve(F,pmult,indexfunction(i,mx),
offset+i*stride);
} else {
ImplicitHConvolution *yconvolve0=yconvolve[0];
for(unsigned int i=0; i < M; ++i)
yconvolve0->convolve(F,pmult,indexfunction(i,mx),offset+i*stride);
}
}
void forwards(Complex **F, Complex **U2, unsigned int offset) {
for(unsigned int b=0; b < B; ++b)
xfftpad->forwards(F[b]+offset,U2[b]);
}
// F is a pointer to A distinct data blocks each of size
// (2mx-compact)*(my+!ycompact), shifted by offset (contents not preserved).
virtual void convolve(Complex **F, realmultiplier *pmult,
bool symmetrize=true, unsigned int i=0,
unsigned int offset=0) {
if(!toplevel) {
index[indexsize-2]=i;
if(threads > 1) {
for(unsigned int t=1; t < threads; ++t) {
unsigned int *Index=yconvolve[t]->index;
for(unsigned int i=0; i < indexsize; ++i)
Index[i]=index[i];
}
}
}
unsigned stride=my+!ycompact;
backwards(F,U2,stride,symmetrize,offset);
subconvolution(F,pmult,xfftpad->findex,2*mx-xcompact,stride,offset);
subconvolution(U2,pmult,xfftpad->uindex,mx+xcompact,stride);
forwards(F,U2,offset);
}
// Binary convolution:
void convolve(Complex *f, Complex *g, bool symmetrize=true) {
Complex *F[]={f,g};
convolve(F,multbinary,symmetrize);
}
};
// In-place implicitly dealiased 3D complex convolution.
class ImplicitConvolution3 : public ThreadBase {
protected:
unsigned int mx,my,mz;
Complex *u1;
Complex *u2;
Complex *u3;
unsigned int A,B;
fftpad *xfftpad;
ImplicitConvolution2 **yzconvolve;
Complex **U3;
bool allocated;
unsigned int indexsize;
bool toplevel;
public:
unsigned int *index;
void initpointers3(Complex **&U3, Complex *u3, unsigned int stride) {
unsigned int C=max(A,B);
U3=new Complex *[C];
for(unsigned int a=0; a < C; ++a)
U3[a]=u3+a*stride;
if(toplevel) allocateindex(2,new unsigned int[2]);
}
void deletepointers3(Complex **&U3) {
if(toplevel) {
delete [] index;
for(unsigned int t=1; t < threads; ++t)
delete [] yzconvolve[t]->index;
}
delete [] U3;
}