forked from rajeshmondal18/FoF-Halo-finder
-
Notifications
You must be signed in to change notification settings - Fork 0
/
nbody_funcs.c
1454 lines (1143 loc) · 39.4 KB
/
nbody_funcs.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
#include<stdio.h>
#include<math.h>
#include<stdlib.h>
#include<string.h>
#include<fftw3.h>
#include"nbody.h"
#include<omp.h>
/*---------------------GLOBAL VARIABLES declaration---------------------*/
// Cosmological parameters read from input file "input.zel"
extern float vhh, // Hubble parameter
vomegam, // Omega_matter; total matter density (baryons+CDM) parameter
vomegalam, // Cosmological Constant
vomegab, //Omega_baryon
tcmb, // CMB temperature
sigma_8_present ,// Last updated value of sigma_8 (Presently WMAP)
vnn; // Spectral index of primordial Power spectrum
extern long N1,N2,N3;// box dimension (grid)
extern int NF, // Fill every NF grid point
Nbin; // Number of bins to calculate final P(k) (output)
extern long MM; // Number of particles
extern float LL; // grid spacing in Mpc
// global variables (calculated )
extern int zel_flag, // memory allocation for zel is 3 times that for nbody
fourier_flag; //for fourier transfrom
extern float DM_m, // Darm matter mass of simulation particle (calculated)
norm, // normalize Pk
pi;
extern io_header header1; // structure for header
extern float ***ro; // for density or potential
extern fftwf_plan p_ro; // for FFT
extern fftwf_plan q_ro; // for FFT
// arrays for storing data
static float ***va; // for grad phi
static float norml, CC; // normalization constants
static float rho_b_inv, // (Mean number of particles/grid cell)^{-1}
Cx,Cy,Cz, // conversion from i,j,k to kx,ky,kz
Lcube,vol, // L^3, Volume=L^3*N1*N2*N3
tpibyL;// 2 pi /LL
static fftwf_plan q_va; // for FFT
static omp_lock_t *ro_lck;
/*********************************************************************************************************/
void Setting_Up_Memory_For_Ro(float av)
{
// for multiple threads
fftwf_init_threads();
fftwf_plan_with_nthreads(omp_get_max_threads());
// fftwf_plan_with_nthreads(1);
omp_set_num_threads(omp_get_max_threads());
// omp_set_num_threads(1);
printf("No of threads = %d\n",omp_get_max_threads());
// done multi thread
ro = allocate_fftwf_3d(N1,N2,N3+2);
va = allocate_fftwf_3d(N1,N2,N3+2); /* The last dimension gets padded because of using REAL FFT */
ro_lck=(omp_lock_t*)&(va[0][0][0]);
/****************************************************/
/* Creating the plans for forward and reverse FFT's */
p_ro = fftwf_plan_dft_r2c_3d(N1, N2, N3, &(ro[0][0][0]), (fftwf_complex*)&(ro[0][0][0]), FFTW_ESTIMATE);
q_ro = fftwf_plan_dft_c2r_3d(N1, N2, N3, (fftwf_complex*)&(ro[0][0][0]), &(ro[0][0][0]), FFTW_ESTIMATE);
q_va = fftwf_plan_dft_c2r_3d(N1, N2, N3, (fftwf_complex*)&(va[0][0][0]), &(va[0][0][0]), FFTW_ESTIMATE);
/****************************************************/
/* Fwd FFT followed immediately by Reverse FFT gives results reduced by a product N1*N2*N3 */
/* The FFT is unnormalized.. So we introduce a normalization factor */
norml=1./((1.0*N1)*(1.0*N2)*(1.0*N3));
Cx=2.*pi/(N1*LL); Cy=2.*pi/(N2*LL); Cz=2.*pi/(N3*LL);
Lcube=powf(LL,3.);
vol=Lcube/norml;
CC=pi*pi*vol*powf(Df(av),2.); // (2*pi^2 Volume DD^2)/2 DD-> initial growing mode
rho_b_inv=1.0/(norml*MM); // inverse of density in units of particles per unit cell
tpibyL=2.0*pi/LL;// 2 pi /LL
printf("norml=%e\tCC=%e\trho_b_inv=%e\n",norml,CC,rho_b_inv);
}
/*******************************************************************************************/
void delta_fill(long* seed)
{
long i,j,k,m,jj1,jj2,kk;
long index1,index2;
float amp,phasev,DD;
fftwf_complex *A;
float ran1(long*);
float gasdev(long*);
float p(long,long,long); // returns amplitude of delta(k)
A=(fftwf_complex*)&(ro[0][0][0]);
/********************FILLING POINTS***********************/
for(i=0;i<2;++i)
for(j=0;j<2;++j)
for(k=0;k<2;++k)
{
if (i==0 && j==0 && k==0)
{
index1=0;
A[index1][0]=0.0;
A[index1][1]=0.0;
}
else
{
index1 = i*(N1/2)*N2*(N3/2+1) + j*(N2/2)*(N3/2+1) + k*(N3/2);
amp=p(i*N1/2,j*N2/2,k*N3/2);
A[index1][0]=amp*gasdev(seed);
A[index1][1]=0.0;
}
}
/********************FILLING LINES***********************/
for(i=0;i<2;++i)
for(j=0;j<2;++j)
{
for(k=1;k<N1/2;k++)
{
amp=p(k,i*N2/2,j*N3/2);
index1 =k*N2*(N3/2+1) +i*(N2/2)*(N3/2+1) + j*(N3/2);
A[index1][0]=amp*gasdev(seed);
A[index1][1]=amp*gasdev(seed);
index2 =(N1-k)*N2*(N3/2+1) +i*(N2/2)*(N3/2+1) + j*(N3/2);
A[index2][1]=A[index1][0];
A[index2][1]=-A[index1][1];
}
for(k=1;k<N2/2;k++)
{
amp=p(i*N1/2,k,j*N3/2);
index1 =i*(N1/2)*N2*(N3/2+1)+ k*(N3/2+1) + j*(N3/2);
index2 =i*(N1/2)*N2*(N3/2+1)+(N2-k)*(N3/2+1) + j*(N3/2);
A[index1][0]=amp*gasdev(seed);
A[index1][1]=amp*gasdev(seed);
A[index2][0]=A[index1][0];
A[index2][1]=-A[index1][1];
}
for(k=1;k<N3/2;k++)
{
amp=p(i*N1/2,j*N2/2,k);
index1 =i*(N1/2)*N2*(N3/2+1)+ j*(N2/2)*(N3/2+1) + k;
A[index1][0]=amp*gasdev(seed);
A[index1][1]=amp*gasdev(seed);
}
}
/**********************FILLING PLANES***********************/
for(i=0;i<2;i++)
for(m=0;m<2;m++)
{
for(j=1;j<N2/2;j++)
for(k=1;k<N3/2;k++)
{
amp=p(i*N1/2,j,k);
jj1=(1-m)*j + m * (N2-j);
index1 =i*(N1/2)*N2*(N3/2+1)+ jj1*(N3/2+1) + k;
A[index1][0]=amp*gasdev(seed);
A[index1][1]=amp*gasdev(seed);
}
for(j=1;j<N1/2;j++)
for(k=1;k<N3/2;k++)
{
amp=p(j,i*N2/2,k);
jj1=(1-m)*j + m * (N1-j);
index1 =jj1*N2*(N3/2+1)+ i*(N2/2)*(N3/2+1) + k;
A[index1][0]=amp*gasdev(seed);
A[index1][1]=amp*gasdev(seed);
}
for(j=1;j<N1/2;j++)
for(k=1;k<N2/2;k++)
{
amp=p(j,k,i*N3/2);
jj1=(1-m)*j + m * (N1-j);
kk=k;
index1 =jj1*N2*(N3/2+1)+ kk*(N3/2+1)+i*N3/2;
jj2=m*j + (1-m) * (N1-j);
kk=N2-k;
index2 =jj2*N2*(N3/2+1)+ kk*(N3/2+1)+i*N3/2;
A[index1][0]=amp*gasdev(seed);
A[index1][1]=amp*gasdev(seed);
A[index2][0]=A[index1][0];
A[index2][1]=-A[index1][1];
}
}
/***************************FILLING CUBES***********************/
for(i=1;i<N1/2;i++)
for(j=1;j<N2/2;j++)
for(k=1;k<N3/2;k++)
{
amp=p(i,j,k);
index1 =i*N2*(N3/2+1) + j*(N3/2+1) + k;
A[index1][0]=amp*gasdev(seed);
A[index1][1]=amp*gasdev(seed);
index1 =(N1-i)*N2*(N3/2+1) + j*(N3/2+1) + k;
A[index1][0]=amp*gasdev(seed);
A[index1][1]=amp*gasdev(seed);
index1 =i*N2*(N3/2+1) + (N2-j)*(N3/2+1) + k;
A[index1][0]=amp*gasdev(seed);
A[index1][1]=amp*gasdev(seed);
index1 = (N1-i)*N2*(N3/2+1) + (N2-j)*(N3/2+1) + k;
A[index1][0]=amp*gasdev(seed);
A[index1][1]=amp*gasdev(seed);
}
}
/**************************************************************************************************/
float p(long i,long j,long k)
// returns amplitude of delta(k)
{
float kk,val,Pk(float),kx,ky,kz;
kx=Cx*i;
ky=Cy*j;
kz=Cz*k;
kk = (float)sqrt((double)(kx*kx+ ky*ky + kz*kz));
val = (float)sqrt(CC*Pk(kk));
return(val);
}
/*******************************************************************************************************/
/*******************************************************************************************************/
/******************** TO FIND POWER SPECTRUM **************************/
void calpow(int f_flag, int Nbin, double* power, double* powerk, double* kmode, long *no)
{
long i, j, k, a, b, c;
int d;
fftwf_complex *comp_ro;
float fac1,fac2,fac3, fac, m, scale;
long index,index1,index2,c1=0;
// f-flag==0 -> ro contains delta(k)
// f-flag==1 -> ro contains delta(x)
/*************** TAKING FOURIER TRANSFORM OF RO. **************/
// do forward transform to get delta(k)
if(f_flag==1)
{
for(i=0;i<N1;i++)
for(j=0;j<N2;j++)
for(k=0;k<N3;k++)
ro[i][j][k]=ro[i][j][k]*Lcube;
fftwf_execute(p_ro); // ro mow contains delta(k)
}
comp_ro = (fftwf_complex *)&(ro[0][0][0]);
/*********** TO FIND POWER SPECTRUM OF RO. **************/
fac1=1./(1.*N1*N1);
fac2=1./(1.*N2*N2);
fac3=1./(1.*N3*N3);
for(i=0;i<Nbin;++i)
{
power[i]=0.0;
powerk[i]=0.0;
kmode[i]=0.0;
no[i]=0;
}
// dor logarithmic bins
scale=log10(0.5*N1)/Nbin;
// ----------------------
/**************** BINNING POWER SPECTRA **********************/
/*-------------------------- half lines ------------------------- */
for(i=1;i<=N1/2;i++)
for(j=0;j<=N2/2;j=j+N2/2)
for(k=0;k<=N3/2;k=k+N3/2)
{
a=i;
b=j;
c=k;
index = i*N2*(N3/2+1) + j*(N3/2+1) + k;
m = sqrt(fac1*a*a + fac2*b*b + fac3*c*c);
d=(int)floorf(log10(m*N1)/scale);//logarithmic bins
if(d>=0 && d<Nbin)
{
power[d]+= ((comp_ro[index][0]* comp_ro[index][0]) + (comp_ro[index][1]* comp_ro[index][1]))*1.;
powerk[d] += Pk((float)(tpibyL*m))*1.0;
kmode[d] += m*1.;
no[d]=no[d]+1;
}
}
/*----------------------- half planes -----------------------*/
for(i=0;i<N1;i++)
{
a=(i>N1/2)? N1-i: i;
index1 = i*N2*(N3/2+1) ;
for(j=1;j<N2/2;j++)
{
b=j;
index2 = index1 + j*(N3/2+1) ;
for(k=0;k<=N3/2;k=k+N3/2)
{
c=k;
index = index2 + k;
m = sqrt(fac1*a*a + fac2*b*b + fac3*c*c);
d=(int)floorf(log10(m*N1)/scale);//logarithmic bins
if(d>=0 && d<Nbin)
{
power[d]+= ((comp_ro[index][0]* comp_ro[index][0]) + (comp_ro[index][1]* comp_ro[index][1]))*1.;
powerk[d] += Pk((float)(tpibyL*m))*1.0;
kmode[d] += m*1.;
no[d]=no[d]+1;
}
}
}
}
/**************** half cube **********************/
for(i=0;i<N1;i++)
{
a=(i>N1/2)? N1-i: i;
index1 = i*N2*(N3/2+1) ;
for(j=0;j<N2;j++)
{
b=(j>N2/2)? N2-j: j;
index2 = index1 + j*(N3/2+1) ;
for(k=1;k<N3/2;k++)
{
c=k;
index = index2 + k;
m = sqrt(fac1*a*a + fac2*b*b + fac3*c*c);
/* m*(2 * pi/LL) is |k| */
/* m=1/2 corresponds to kmode[Nbin-1] i.e. Nyquits */
d=(int)floorf(log10(m*N1)/scale);//logarithmic bins
if (d>=0 && d<Nbin)
{
power[d] += ((comp_ro[index][0]* comp_ro[index][0]) + (comp_ro[index][1]* comp_ro[index][1]))*1.;
powerk[d] += Pk((float)(tpibyL*m))*1.0;
kmode[d] += m*1.;
no[d]=no[d]+1;
}
} /* end k for loop */
}
}
for(i=0;i<Nbin;i++)
{
if (no[i]>0)
{
power[i] = power[i]/(1.*no[i]*vol);
powerk[i] = powerk[i]/(1.*no[i]);
kmode[i]=tpibyL*kmode[i]/(1.*no[i]);
}
}
// do back transform to get delta(x)
if(f_flag==1)
{
for(i=0;i<N1;i++)
{
index1 = i*N2*(N3/2+1) ;
for(j=0;j<N2;j++)
{
index2=index1 + j*(N3/2+1) ;
for(k=0;k<N3/2;k++)
{
index=index2 + k;
comp_ro[index][0]=comp_ro[index][0]/vol;
comp_ro[index][1]=comp_ro[index][1]/vol;
}
}
}
/* now convert the array back to real space */
fftwf_execute(q_ro);
}
}
/**********************************************************************/
/******************** TO FIND POWER SPECTRUM **************************/
void calpow_k(int f_flag, float kmin, float kmax, int Nbin,double* power, double* powerk, double* kmode, long *no)
{
long i, j, k, a, b, c;
int d;
fftwf_complex *comp_ro;
float fac1,fac2,fac3, fac, m, scale;
long index,index1,index2;
// f-flag==0 -> ro contains delta(k)
// f-flag==1 -> ro contains delta(x)
/*************** TAKING FOURIER TRANSFORM OF RO. **************/
// do forward transform to get delta(k)
if(f_flag==1)
{
for(i=0;i<N1;i++)
for(j=0;j<N2;j++)
for(k=0;k<N3;k++)
ro[i][j][k]=ro[i][j][k]*Lcube;
fftwf_execute(p_ro); // ro mow contains delta(k)
}
comp_ro = (fftwf_complex *)&(ro[0][0][0]);
/*********** TO FIND POWER SPECTRUM OF RO. **************/
fac1=1./(1.*N1*N1);
fac2=1./(1.*N2*N2);
fac3=1./(1.*N3*N3);
for(i=0;i<Nbin;++i)
{
power[i]=0.0;
powerk[i]=0.0;
kmode[i]=0.0;
no[i]=0;
}
// dor logarithmic bins
scale=(log10(kmax)-log10(kmin))/Nbin;
// ----------------------
/**************** BINNING POWER SPECTRA **********************/
/*-------------------------- half lines ------------------------- */
for(i=1;i<=N1/2;i++)
for(j=0;j<=N2/2;j=j+N2/2)
for(k=0;k<=N3/2;k=k+N3/2)
{
a=i;
b=j;
c=k;
index = i*N2*(N3/2+1) + j*(N3/2+1) + k;
m = tpibyL*sqrt(fac1*a*a + fac2*b*b + fac3*c*c); /* m is |k| */
d=(int)floorf((log10(m)-log10(kmin))/scale); //logarithmic bins
if(d>=0 && d<Nbin)
{
power[d]+= ((comp_ro[index][0]* comp_ro[index][0]) + (comp_ro[index][1]* comp_ro[index][1]))*1.;
powerk[d] += Pk((float)(m))*1.0;
kmode[d] += m*1.;
no[d]=no[d]+1;
}
}
/*----------------------- half planes -----------------------*/
for(i=0;i<N1;i++)
{
a=(i>N1/2)? N1-i: i;
index1 = i*N2*(N3/2+1) ;
for(j=1;j<N2/2;j++)
{
b=j;
index2 = index1 + j*(N3/2+1) ;
for(k=0;k<=N3/2;k=k+N3/2)
{
c=k;
index = index2 + k;
m = tpibyL*sqrt(fac1*a*a + fac2*b*b + fac3*c*c); /* m is |k| */
d=(int)floorf((log10(m)-log10(kmin))/scale); //logarithmic bins
if(d>=0 && d<Nbin)
{
power[d]+= ((comp_ro[index][0]* comp_ro[index][0]) + (comp_ro[index][1]* comp_ro[index][1]))*1.;
powerk[d] += Pk((float)(m))*1.0;
kmode[d] += m*1.;
no[d]=no[d]+1;
}
}
}
}
/**************** half cube **********************/
for(i=0;i<N1;i++)
{
a=(i>N1/2)? N1-i: i;
index1 = i*N2*(N3/2+1) ;
for(j=0;j<N2;j++)
{
b=(j>N2/2)? N2-j: j;
index2 = index1 + j*(N3/2+1) ;
for(k=1;k<N3/2;k++)
{
c=k;
index = index2 + k;
m = tpibyL*sqrt(fac1*a*a + fac2*b*b + fac3*c*c); /* m is |k| */
d=(int)floorf((log10(m)-log10(kmin))/scale); //logarithmic bins
if (d>=0 && d<Nbin)
{
power[d] += ((comp_ro[index][0]* comp_ro[index][0]) + (comp_ro[index][1]* comp_ro[index][1]))*1.;
powerk[d] += Pk((float)(m))*1.0;
kmode[d] += m*1.;
no[d]=no[d]+1;
}
} /* end k for loop */
}
}
for(i=0;i<Nbin;i++)
{
if (no[i]>0)
{
power[i] = power[i]/(1.*no[i]*vol);
powerk[i] = powerk[i]/(1.*no[i]);
kmode[i] = kmode[i]/(1.*no[i]);
}
}
// do back transform to get delta(x)
if(f_flag==1)
{
for(i=0;i<N1;i++)
{
index1 = i*N2*(N3/2+1) ;
for(j=0;j<N2;j++)
{
index2=index1 + j*(N3/2+1) ;
for(k=0;k<N3/2;k++)
{
index=index2 + k;
comp_ro[index][0]=comp_ro[index][0]/vol;
comp_ro[index][1]=comp_ro[index][1]/vol;
}
}
}
/* now convert the array back to real space */
fftwf_execute(q_ro);
}
}
/*******************************************************************************************************/
/*******************************************************************************************************/
void grad_phi(int ix)
{
// at start of this function: ro contains delta(k)
// at end of function va contains the ix component of Grad[Laplacian^-1[delta(x)]]
// this is used in ZA
long ii,jj,kk;
long index;
float AA,a0,a1,a2;
fftwf_complex *vc,*roc;
vc=(fftwf_complex*)&(va[0][0][0]);
roc=(fftwf_complex*)&(ro[0][0][0]);
#pragma omp parallel for private(jj,kk,a0,a1,a2,AA,index)
for(ii=0;ii<N1;ii++)
{
a0 = (ii>N1/2)? Cx*(ii-N1) : Cx*ii; // kx
for(jj=0;jj<N2;jj++)
{
a1 = (jj>N2/2)? Cy*(jj-N2) : Cy*jj; // ky
kk=0;
a2=0.;
index = (ii*N2+jj)*(N3/2+1) + kk;
AA = pow(a0,2.)+pow(a1,2.)+pow(a2,2.); // AA=K^2
AA = (fabs(AA)>1.e-10)? (((ix-1)*(ix-2)/2)*a0 - ix*(ix-2)*a1 + (ix*(ix-1)/2)*a2)/AA : 0.0; // AA=k_ix/K^2 or 0
vc[index][0]=-1.0*AA*roc[index][1]/vol;
vc[index][1]=AA*roc[index][0]/vol;
for(kk=1;kk<N3/2+1;kk++)
{
a2=kk*Cz; // kz
index = (ii*N2+jj)*(N3/2+1) + kk;
AA= pow(a0,2.)+pow(a1,2.)+pow(a2,2.); // AA=K^2
AA=(((ix-1)*(ix-2)/2)*a0 - ix*(ix-2)*a1 + (ix*(ix-1)/2)*a2)/AA; // AA=k_ix/K^2
vc[index][0] = -1.0*AA*roc[index][1]/vol;
vc[index][1] = AA*roc[index][0]/vol;
}
}
}
fftwf_execute(q_va); /* Take fourier transform */
}
//############################################################################################
//############################################################################################
void Zel_move_gradphi(float av,float **rra,float **vva)
{
float DD, // growing mode
vfac; // converts to peculiar velocity
int ii;
long jj, kk, ll, N, pos;
long pin;
vfac = av*av*Hf(av)*ff(av)/LL;
for(ii=0;ii<3;++ii)
{
grad_phi(ii); // Calculate ix component of -Grad[Laplacian^{_1}[delta]] in real space
#pragma omp parallel for private(jj,kk,ll,pos,pin,N)
for(jj=0;jj<N1/NF;jj++)
for(kk=0;kk<N2/NF;kk++)
for(ll=0;ll<N3/NF;ll++)
{
pin = jj*(N2/NF)*(N3/NF) + kk*(N3/NF) + ll;
pos= (ii-1)*(ii-2)*jj/2 - ii*(ii-2)*kk + ii*(ii-1)*ll/2;
rra[pin][ii]=(float)(NF*pos)+ va[NF*jj][NF*kk][NF*ll]/LL;
vva[pin][ii]=vfac*va[NF*jj][NF*kk][NF*ll];
// periodic boundary condition
N = (ii-1)*(ii-2)*N1/2 - ii*(ii-2)*N2 + ii*(ii-1)*N3/2;
/* to ensure the values are not negative */
rra[pin][ii]=rra[pin][ii]+N;
rra[pin][ii]= rra[pin][ii]-1.0*N*(int)(floor(rra[pin][ii])/(1.*N));
}
/*------------------- done ---------------------*/
}
} /* end of Zel_move_gradphi */
//######################################################################################
void cic(float **rra)
// //For a given particle distribution
// This uses Cloud in Cell(CIC) to calculate (1 + delta) on the Grid
{
long i, j, k, ix, jy, kz;
int ii, jj, kk;
long pin,index;
float wx,wy,wz;
/* Clear out the array ro. ******/
#pragma omp parallel for private(i,j,k,index)
for(i=0;i<N1;i++)
for(j=0;j<N2;j++)
for(k=0;k<N3;k++)
{
index = (i*N2+j)*N3 + k;
ro[i][j][k] = 0.0;
omp_init_lock(&ro_lck[index]);
}
/********************************/
#pragma omp parallel for private(i, j, k, ii, jj, kk, wx, ix, wy, jy, wz, kz, index)
for(pin=0;pin<MM;pin++)
{ /* begin particle index loop */
/* (a/b/c)[0] or (a/b/c)[1] can never be greater than (N1/N2/N3) */
// left bottom corner of cell containing the particle
i = (long)floor(rra[pin][0]);
j = (long)floor(rra[pin][1]);
k = (long)floor(rra[pin][2]);
/* for each of the 8 corner points */
for(ii=0;ii<=1;ii++)
{
wx=fabs(1.-rra[pin][0]+i-ii)*rho_b_inv; // divide by mean density
ix=(i+ii)%N1;
for(jj=0;jj<=1;jj++)
{
wy=fabs(1.-rra[pin][1]+j-jj);
jy=(j+jj)%N2;
for(kk=0;kk<=1;kk++)
{
wz=fabs(1.-rra[pin][2]+k-kk);
kz=(k+kk)%N3;
index=(ix*N2+jy)*N3 + kz;
omp_set_lock(&ro_lck[index]);
ro[ix][jy][kz]+=wx*wy*wz;
omp_unset_lock(&ro_lck[index]);
}
}
} /* end of 8 grid corners loop> */
} /* end of each particle loop */
} /* end of function cic */
//############################################################################################
//############################################################################################
void cic_sampled(float **rra, int *s_indx)
// //For a given particle distribution
// This uses Cloud in Cell(CIC) to calculate (1 + delta) on the Grid
{
long i, j, k, ix, jy, kz;
int ii, jj, kk;
long pin,index;
float wx,wy,wz;
/* Clear out the array ro. ******/
#pragma omp parallel for private(i,j,k,index)
for(i=0;i<N1;i++)
for(j=0;j<N2;j++)
for(k=0;k<N3;k++)
{
index = (i*N2+j)*N3 + k;
ro[i][j][k] = 0.0;
omp_init_lock(&ro_lck[index]);
}
/********************************/
#pragma omp parallel for private(i, j, k, ii, jj, kk, wx, ix, wy, jy, wz, kz, index)
for(pin=0;pin<MM;pin++)
if(s_indx[pin]==-1)
{ /* begin particle index loop */
/* (a/b/c)[0] or (a/b/c)[1] can never be greater than (N1/N2/N3) */
// left bottom corner of cell containing the particle
i = (long)floor(rra[pin][0]);
j = (long)floor(rra[pin][1]);
k = (long)floor(rra[pin][2]);
/* for each of the 8 corner points */
for(ii=0;ii<=1;ii++)
{
wx=fabs(1.-rra[pin][0]+i-ii)*rho_b_inv; // divide by mean density
ix=(i+ii)%N1;
for(jj=0;jj<=1;jj++)
{
wy=fabs(1.-rra[pin][1]+j-jj);
jy=(j+jj)%N2;
for(kk=0;kk<=1;kk++)
{
wz=fabs(1.-rra[pin][2]+k-kk);
kz=(k+kk)%N3;
index=(ix*N2+jy)*N3 + kz;
omp_set_lock(&ro_lck[index]);
ro[ix][jy][kz]+=wx*wy*wz;
omp_unset_lock(&ro_lck[index]);
}
}
} /* end of 8 grid corners loop> */
} /* end of each particle loop */
} /* end of function cic_sampled */
//##########################################################################################
void Get_phi(int f_flag) // calculates Laplacian^{-1}[ ro]]
{
long ii,jj,kk;
long index;
float AA,a0,a1,a2;
fftwf_complex *roc;
/*************** TAKING FOURIER TRANSFORM OF RO. **************/
// do forward transform to get ro(k)
if(f_flag==1)
{
for(ii=0;ii<N1;ii++)
for(jj=0;jj<N2;jj++)
for(kk=0;kk<N3;kk++)
ro[ii][jj][kk]*=Lcube;
fftwf_execute(p_ro); // ro mow contains delta(k)
}
//------------------------------------------
roc=(fftwf_complex*)&(ro[0][0][0]);
#pragma omp parallel for private(jj,kk,a0,a1,a2,AA,index)
for(ii=0;ii<N1;ii++)
{
a0=ii*pi/N1; // kx *LL/2
a0=pow(2.*sin(a0)/LL,2.);
for(jj=0;jj<N2;jj++)
{
a1=jj*pi/N2; // ky *LL/2
a1=pow(2.*sin(a1)/LL,2.);
kk=0;
index = (ii*N2+jj)*(N3/2+1) + kk;
AA=a0+a1;
AA=(fabs(AA)>1.e-10)? 1/AA:0.; // AA=k_ix/K^2 or 0
roc[index][0] = -1.*AA*roc[index][0]/vol;
roc[index][1] = -1.*AA*roc[index][1]/vol;
for(kk=1;kk<N3/2+1;kk++)
{
a2=kk*pi/N3; // kz *LL/2
a2=pow(2.*sin(a2)/LL,2.);
index = (ii*N2+jj)*(N3/2+1) + kk;
AA=1./(a0+a1+a2);
roc[index][0] = -1.*AA*roc[index][0]/vol;
roc[index][1] = -1.*AA*roc[index][1]/vol;
}
}
}
fftwf_execute(q_ro); /* Take fourier transform */
}
//#############################################################################################
void Update_v(float av,float delta_aa,float **rra,float **vva) // udates v if (delta a >0) else intializes v
{
long pin;
int ii,jj,kk;
long xp,yp,zp,xn,yn,zn; /* (x/y/z)(p/n) = (x/y/z)(previous/next) */
long a,b,c,ix,jy,kz;
float wx,wy,wz,g0,g1,g2;
float xx,yy,zz;
float coeff,Hf(float aa),vflag;
cic(rra); // calculate delta(x) from particle distribution
Get_phi(1);// Calculate [Laplacian^{_1}[delta]] in real space
if(delta_aa>0.)
{
coeff=delta_aa*1.5* vomegam/(LL*av*av*Hf(av)); // delta_aa scaled
vflag=1.0;
}
else
{
coeff = av*av*Hf(av)*ff(av)/LL;
vflag=0.0;
}
coeff/=(2.*LL); // divide by 2 LL for Grad on grid
#pragma omp parallel for private(a,b,c,xx,yy,zz,g0,g1,g2,ii,jj,kk,ix,wx,xp,xn,jy,wy,yn,yp,kz,wz,zp,zn)
for(pin=0;pin<MM;pin++)
{
/* left most corner of the cube enclosing the particle */
a = (long)floor(rra[pin][0]);
b = (long)floor(rra[pin][1]);
c = (long)floor(rra[pin][2]);
/* particle co-ordinates itself */
xx = rra[pin][0];
yy = rra[pin][1];
zz = rra[pin][2];
g0=0;
g1=0;
g2=0;
/* begin 8 corners loop */
for(ii=0;ii<=1;ii++)
{
ix = (a+ii)%N1;
wx=fabs(1.- xx +a-ii);
xp = ((ix-1)==-1) ? N1-1: ix-1;
xn = ((ix+1)==N1) ? 0 : ix+1;
for(jj=0;jj<=1;jj++)
{
jy = (b+jj)%N2;
wy=fabs(1.- yy +b-jj);
yp = ((jy-1)==-1) ? N2-1: jy-1;
yn = ((jy+1)==N2) ? 0 : jy+1;
for(kk=0;kk<=1;kk++)
{
kz = (c+kk)%N3;
wz=fabs(1.- zz +c-kk);
zp = ((kz-1)==-1) ? N3-1: kz-1;
zn = ((kz+1)==N3) ? 0 : kz+1;
/* ix,jy,kz are the current co-ordinates of the cube vertex point */
/* calculating the difference from the respective corner */