ngenic.cc 40 KB
Newer Older
Volker Springel's avatar
Volker Springel committed
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
/*******************************************************************************
 * \copyright   This file is part of the GADGET4 N-body/SPH code developed
 * \copyright   by Volker Springel. Copyright (C) 2014-2020 by Volker Springel
 * \copyright   (vspringel@mpa-garching.mpg.de) and all contributing authors.
 *******************************************************************************/

/*! \file  ngenic.cc
 *
 *  \brief sets up cosmological initial conditions
 */

#include "gadgetconfig.h"

#ifdef NGENIC

#include <gsl/gsl_rng.h>
#include <math.h>
#include <mpi.h>
#include <stdlib.h>
#include <algorithm>

#include "../data/allvars.h"
#include "../data/dtypes.h"
#include "../data/intposconvert.h"
#include "../data/mymalloc.h"
#include "../logs/logs.h"
#include "../main/simulation.h"
#include "../mpi_utils/mpi_utils.h"
#include "../mpi_utils/shared_mem_handler.h"
#include "../ngenic/ngenic.h"
#include "../pm/pm_mpi_fft.h"
#include "../system/system.h"

#ifdef GRIDX
#undef GRIDX
#undef GRIDY
#undef GRIDZ
#undef INTCELL
#endif

#define GRIDX (NGENIC)
#define GRIDY (NGENIC)
#define GRIDZ (NGENIC)

#define INTCELL ((~((MyIntPosType)0)) / GRIDX + 1)

#define GRIDz (GRIDZ / 2 + 1)
#define GRID2 (2 * GRIDz)

#define FI(x, y, z) (((large_array_offset)GRID2) * (GRIDY * (x) + (y)) + (z))
#define FC(c, z) (((large_array_offset)GRID2) * ((c)-myplan.firstcol_XY) + (z))

#if(GRIDZ > 1024)
typedef long long large_array_offset; /* use a larger data type in this case so that we can always address all cells of the 3D grid
                                         with a single index */
#else
typedef unsigned int large_array_offset;
#endif

#ifdef NUMPART_PER_TASK_LARGE
typedef long long large_numpart_type; /* if there is a risk that the local particle number times 8 overflows a 32-bit integer, this
                                         data type should be used */
#else
typedef int large_numpart_type;
#endif

void ngenic::ngenic_displace_particles(void)
{
  TIMER_START(CPU_NGENIC);

  mpi_printf("NGENIC: computing displacement fields...\n");

  All.set_cosmo_factors_for_current_time();

  double vel_prefac1 = All.cf_atime * All.cf_hubble_a * ngenic_f1_omega(All.cf_atime);
  double vel_prefac2 = All.cf_atime * All.cf_hubble_a * ngenic_f2_omega(All.cf_atime);

  vel_prefac1 /= sqrt(All.cf_atime); /* converts to Gadget velocity */
  vel_prefac2 /= sqrt(All.cf_atime); /* converts to Gadget velocity */

  mpi_printf("NGENIC: vel_prefac1= %g  hubble_a=%g   fom1=%g\n", vel_prefac1, All.cf_hubble_a, ngenic_f1_omega(All.cf_atime));
  mpi_printf("NGENIC: vel_prefac2= %g  hubble_a=%g   fom2=%g\n", vel_prefac2, All.cf_hubble_a, ngenic_f2_omega(All.cf_atime));

  rnd_generator_conjugate = gsl_rng_alloc(gsl_rng_ranlxd1);
  rnd_generator           = gsl_rng_alloc(gsl_rng_ranlxd1);
  gsl_rng_set(rnd_generator, All.NgenicSeed);

  ngenic_initialize_powerspectrum();

  ngenic_initialize_ffts();

  if(!(seedtable = (unsigned int *)Mem.mymalloc("seedtable", NGENIC * NGENIC * sizeof(unsigned int))))
    Terminate("could not allocate seed table");

  for(int i = 0; i < NGENIC / 2; i++)
    {
      for(int j = 0; j < i; j++)
        seedtable[i * NGENIC + j] = 0x7fffffff * gsl_rng_uniform(rnd_generator);

      for(int j = 0; j < i + 1; j++)
        seedtable[j * NGENIC + i] = 0x7fffffff * gsl_rng_uniform(rnd_generator);

      for(int j = 0; j < i; j++)
        seedtable[(NGENIC - 1 - i) * NGENIC + j] = 0x7fffffff * gsl_rng_uniform(rnd_generator);

      for(int j = 0; j < i + 1; j++)
        seedtable[(NGENIC - 1 - j) * NGENIC + i] = 0x7fffffff * gsl_rng_uniform(rnd_generator);

      for(int j = 0; j < i; j++)
        seedtable[i * NGENIC + (NGENIC - 1 - j)] = 0x7fffffff * gsl_rng_uniform(rnd_generator);

      for(int j = 0; j < i + 1; j++)
        seedtable[j * NGENIC + (NGENIC - 1 - i)] = 0x7fffffff * gsl_rng_uniform(rnd_generator);

      for(int j = 0; j < i; j++)
        seedtable[(NGENIC - 1 - i) * NGENIC + (NGENIC - 1 - j)] = 0x7fffffff * gsl_rng_uniform(rnd_generator);

      for(int j = 0; j < i + 1; j++)
        seedtable[(NGENIC - 1 - j) * NGENIC + (NGENIC - 1 - i)] = 0x7fffffff * gsl_rng_uniform(rnd_generator);
    }

  if(Shmem.Island_NTask != Shmem.World_NTask)
    {
      // We actually have multiple shared memory nodes in which we set aside one MPI rank for shared memory communictiona.
      // In this casem move the seedtable to the communication rank in order to consume this memory only once on the node

      if(Shmem.Island_ThisTask == 0)
        {
          size_t tab_len = NGENIC * NGENIC * sizeof(unsigned int);

          MPI_Send(&tab_len, sizeof(tab_len), MPI_BYTE, Shmem.MyShmRankInGlobal, TAG_TABLE_ALLOC, MPI_COMM_WORLD);
          MPI_Send(seedtable, tab_len, MPI_BYTE, Shmem.MyShmRankInGlobal, TAG_DMOM, MPI_COMM_WORLD);
        }

      Mem.myfree(seedtable);

      ptrdiff_t off;
      MPI_Bcast(&off, sizeof(ptrdiff_t), MPI_BYTE, Shmem.Island_NTask - 1, Shmem.SharedMemComm);

      seedtable = (unsigned int *)((char *)Shmem.SharedMemBaseAddr[Shmem.Island_NTask - 1] + off);
    }

  ngenic_distribute_particles();

  /* allocate displacement vectors */
  Pdisp = (disp_data *)Mem.mymalloc_clear("disp_data", Sp->NumPart * sizeof(disp_data));

#if defined(MULTICOMPONENTGLASSFILE) && defined(DIFFERENT_TRANSFER_FUNC)
  for(Type = MinType; Type <= MaxType; Type++)
#endif
    {
#ifdef NGENIC_2LPT

      /* allocate temporary buffers for second derivatives */
      fft_real *d2phi1[3];
      for(int axes = 0; axes < 3; axes++)
        d2phi1[axes] = (fft_real *)Mem.mymalloc_clear("d2Phi1", maxfftsize * sizeof(fft_real));

      for(int axes = 0; axes < 3; axes++)
        {
          mpi_printf("NGENIC_2LPT: Computing secondary source term, derivatices %d %d\n", axes, axes);

          fft_real *disp = (fft_real *)Mem.mymalloc("disp", maxfftsize * sizeof(fft_real));

          ngenic_setup_modes_in_kspace((fft_complex *)disp);
          ngenic_get_derivate_from_fourier_field(axes, axes, (fft_complex *)disp);

          memcpy(d2phi1[axes], disp, maxfftsize * sizeof(fft_real));

          Mem.myfree(disp);
        }

      /* allocate second source potential */
      fft_real *Phi2 = (fft_real *)Mem.mymalloc_movable(&Phi2, "Phi2", maxfftsize * sizeof(fft_real));

      for(size_t n = 0; n < maxfftsize; n++)
        Phi2[n] = d2phi1[0][n] * d2phi1[1][n] + d2phi1[0][n] * d2phi1[2][n] + d2phi1[1][n] * d2phi1[2][n];

      for(int axes = 2; axes >= 0; axes--)
        Mem.myfree_movable(d2phi1[axes]);

      for(int i = 0; i < 3; i++)
        for(int j = i + 1; j < 3; j++)
          {
            mpi_printf("NGENIC_2LPT: Computing secondary source term, derivatices %d %d\n", i, j);

            fft_real *disp = (fft_real *)Mem.mymalloc("disp", maxfftsize * sizeof(fft_real));

            ngenic_setup_modes_in_kspace((fft_complex *)disp);
            ngenic_get_derivate_from_fourier_field(i, j, (fft_complex *)disp);

            for(size_t n = 0; n < maxfftsize; n++)
              Phi2[n] -= disp[n] * disp[n];

            Mem.myfree(disp);
          }

      mpi_printf("NGENIC_2LPT: Secondary source term computed in real space\n");

      /* Do a forward inplace-FFT to get Phi2 in Fourier space */
      ngenic_compute_transform_of_source_potential(Phi2);

      mpi_printf("NGENIC_2LPT: Done transforming it to k-space\n");

      for(int axes = 0; axes < 3; axes++)
        {
          mpi_printf("NGENIC_2LPT: Obtaining second order displacements for axes=%d\n", axes);

          fft_real *disp = (fft_real *)Mem.mymalloc("disp", maxfftsize * sizeof(fft_real));

          memcpy(disp, Phi2, maxfftsize * sizeof(fft_real));

          ngenic_get_derivate_from_fourier_field(axes, -1, (fft_complex *)disp);

          ngenic_readout_disp(disp, axes, 3.0 / 7, 3.0 / 7 * vel_prefac2);

          Mem.myfree(disp);
        }

      Mem.myfree(Phi2);
#endif

      /* now carry out Zeldovich approximation, yielding first order displacements */
      for(int axes = 0; axes < 3; axes++)
        {
          mpi_printf("NGENIC_2LPT: Obtaining Zeldovich displacements for axes=%d\n", axes);

          fft_real *disp = (fft_real *)Mem.mymalloc("disp", maxfftsize * sizeof(fft_real));

          ngenic_setup_modes_in_kspace((fft_complex *)disp);

          ngenic_get_derivate_from_fourier_field(axes, -1, (fft_complex *)disp);

          ngenic_readout_disp(disp, axes, 1.0, vel_prefac1);

          Mem.myfree(disp);
        }
    }

  /* now add displacement to Lagrangian coordinates  */
  double maxdisp = 0;
  double maxvel  = 0;
  for(int n = 0; n < Sp->NumPart; n++)
    {
      double posdiff[3] = {Pdisp[n].deltapos[0], Pdisp[n].deltapos[1], Pdisp[n].deltapos[2]};

      MyIntPosType delta[3];
      Sp->pos_to_signedintpos(posdiff, (MySignedIntPosType *)delta);

      for(int axes = 0; axes < 3; axes++)
        {
          Sp->P[n].IntPos[axes] += delta[axes];

          if(Pdisp[n].deltapos[axes] > maxdisp)
            maxdisp = Pdisp[n].deltapos[axes];

          if(Sp->P[n].Vel[axes] > maxvel)
            maxvel = Sp->P[n].Vel[axes];
        }
    }

  double max_disp_global, maxvel_global;
  MPI_Reduce(&maxdisp, &max_disp_global, 1, MPI_DOUBLE, MPI_MAX, 0, Communicator);
  MPI_Reduce(&maxvel, &maxvel_global, 1, MPI_DOUBLE, MPI_MAX, 0, Communicator);

  mpi_printf("\nNGENIC: Maximum displacement: %g, in units of the part-spacing= %g\n\n", max_disp_global,
             max_disp_global / (All.BoxSize / NGENIC));
  mpi_printf("\nNGENIC: Maximum velocity component: %g\n\n", maxvel_global);

  Mem.myfree(Pdisp);

  Mem.myfree(partin);
  Mem.myfree(Rcvpm_offset);
  Mem.myfree(Rcvpm_count);
  Mem.myfree(Sndpm_offset);
  Mem.myfree(Sndpm_count);

  if(Shmem.Island_NTask != Shmem.World_NTask)
    {
      if(Shmem.Island_ThisTask == 0)
        {
          // need to send this flag to the correct processor rank (our shared memory handler) so that the table is freed there
          size_t tab_len = NGENIC * NGENIC * sizeof(unsigned int);
          MPI_Send(&tab_len, sizeof(tab_len), MPI_BYTE, Shmem.MyShmRankInGlobal, TAG_TABLE_FREE, MPI_COMM_WORLD);
        }
    }
  else
    {
      Mem.myfree(seedtable);
    }

#ifndef FFT_COLUMN_BASED
  my_slab_based_fft_free(&myplan);
#else
  my_column_based_fft_free(&myplan);
#endif

  FFTW(destroy_plan)(myplan.forward_plan_zdir);
  FFTW(destroy_plan)(myplan.forward_plan_ydir);
  FFTW(destroy_plan)(myplan.forward_plan_xdir);

  FFTW(destroy_plan)(myplan.backward_plan_zdir);
  FFTW(destroy_plan)(myplan.backward_plan_ydir);
  FFTW(destroy_plan)(myplan.backward_plan_xdir);

306
307
  print_spec();

Volker Springel's avatar
Volker Springel committed
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
1001
1002
1003
1004
1005
1006
1007
1008
1009
1010
1011
1012
1013
1014
1015
1016
1017
1018
1019
1020
1021
1022
1023
1024
1025
1026
1027
1028
1029
1030
1031
1032
1033
1034
1035
1036
1037
1038
1039
1040
1041
1042
1043
1044
1045
1046
1047
1048
1049
1050
1051
1052
1053
1054
1055
1056
1057
1058
1059
1060
1061
1062
1063
1064
1065
1066
1067
1068
1069
1070
1071
1072
1073
1074
1075
1076
1077
1078
1079
1080
1081
1082
1083
1084
1085
1086
1087
1088
1089
1090
1091
1092
1093
1094
1095
1096
1097
1098
1099
1100
1101
1102
1103
1104
1105
1106
1107
1108
1109
1110
1111
1112
1113
1114
1115
1116
1117
1118
1119
1120
1121
1122
1123
1124
1125
1126
1127
1128
1129
1130
1131
1132
1133
1134
1135
1136
1137
1138
1139
1140
1141
1142
1143
1144
1145
1146
1147
1148
1149
1150
1151
1152
1153
1154
1155
1156
1157
1158
1159
1160
1161
1162
1163
1164
1165
1166
1167
1168
1169
1170
1171
1172
1173
1174
1175
1176
1177
1178
1179
1180
1181
1182
1183
1184
1185
1186
1187
1188
  if(All.PowerSpectrumType == 2)
    free_power_table();

  gsl_rng_free(rnd_generator);
  gsl_rng_free(rnd_generator_conjugate);

  TIMER_STOP(CPU_NGENIC);
}

void ngenic::ngenic_distribute_particles(void)
{
  Sndpm_count  = (size_t *)Mem.mymalloc("Sndpm_count", NTask * sizeof(size_t));
  Sndpm_offset = (size_t *)Mem.mymalloc("Sndpm_offset", NTask * sizeof(size_t));
  Rcvpm_count  = (size_t *)Mem.mymalloc("Rcvpm_count", NTask * sizeof(size_t));
  Rcvpm_offset = (size_t *)Mem.mymalloc("Rcvpm_offset", NTask * sizeof(size_t));

#ifdef FFT_COLUMN_BASED
  int columns         = GRIDX * GRIDY;
  int avg             = (columns - 1) / NTask + 1;
  int exc             = NTask * avg - columns;
  int tasklastsection = NTask - exc;
  int pivotcol        = tasklastsection * avg;
#endif

  /* determine the slabs/columns each particles accesses */
  {
    size_t *send_count = Sndpm_count;

    for(int j = 0; j < NTask; j++)
      send_count[j] = 0;

    for(int i = 0; i < Sp->NumPart; i++)
      {
        int slab_x  = Sp->P[i].IntPos[0] / INTCELL;
        int slab_xx = slab_x + 1;

        if(slab_xx >= GRIDX)
          slab_xx = 0;

#ifndef FFT_COLUMN_BASED
        int task0 = myplan.slab_to_task[slab_x];
        int task1 = myplan.slab_to_task[slab_xx];

        send_count[task0]++;
        if(task0 != task1)
          send_count[task1]++;
#else
        int slab_y  = Sp->P[i].IntPos[1] / INTCELL;
        int slab_yy = slab_y + 1;

        if(slab_yy >= GRIDY)
          slab_yy = 0;

        int column0 = slab_x * GRIDY + slab_y;
        int column1 = slab_x * GRIDY + slab_yy;
        int column2 = slab_xx * GRIDY + slab_y;
        int column3 = slab_xx * GRIDY + slab_yy;

        int task0, task1, task2, task3;

        if(column0 < pivotcol)
          task0 = column0 / avg;
        else
          task0 = (column0 - pivotcol) / (avg - 1) + tasklastsection;

        if(column1 < pivotcol)
          task1 = column1 / avg;
        else
          task1 = (column1 - pivotcol) / (avg - 1) + tasklastsection;

        if(column2 < pivotcol)
          task2 = column2 / avg;
        else
          task2 = (column2 - pivotcol) / (avg - 1) + tasklastsection;

        if(column3 < pivotcol)
          task3 = column3 / avg;
        else
          task3 = (column3 - pivotcol) / (avg - 1) + tasklastsection;

        send_count[task0]++;
        if(task1 != task0)
          send_count[task1]++;
        if(task2 != task1 && task2 != task0)
          send_count[task2]++;
        if(task3 != task0 && task3 != task1 && task3 != task2)
          send_count[task3]++;
#endif
      }
  }

  /* collect thread-specific offset table and collect the results from the other threads */
  Sndpm_offset[0] = 0;
  for(int i = 1; i < NTask; i++)
    {
      int ind      = i;
      int ind_prev = i - 1;

      Sndpm_offset[ind] = Sndpm_offset[ind_prev] + Sndpm_count[ind_prev];
    }

  MPI_Alltoall(Sndpm_count, sizeof(size_t), MPI_BYTE, Rcvpm_count, sizeof(size_t), MPI_BYTE, Communicator);

  nimport = 0, nexport = 0, Rcvpm_offset[0] = 0, Sndpm_offset[0] = 0;
  for(int j = 0; j < NTask; j++)
    {
      nexport += Sndpm_count[j];
      nimport += Rcvpm_count[j];

      if(j > 0)
        {
          Sndpm_offset[j] = Sndpm_offset[j - 1] + Sndpm_count[j - 1];
          Rcvpm_offset[j] = Rcvpm_offset[j - 1] + Rcvpm_count[j - 1];
        }
    }

  /* allocate import and export buffer */
  partin  = (partbuf *)Mem.mymalloc("partin", nimport * sizeof(partbuf));
  partout = (partbuf *)Mem.mymalloc("partout", nexport * sizeof(partbuf));

  {
    size_t *send_count  = Sndpm_count;
    size_t *send_offset = Sndpm_offset;

    for(int j = 0; j < NTask; j++)
      send_count[j] = 0;

    /* fill export buffer */
    for(int i = 0; i < Sp->NumPart; i++)
      {
        int slab_x  = Sp->P[i].IntPos[0] / INTCELL;
        int slab_xx = slab_x + 1;

        if(slab_xx >= GRIDX)
          slab_xx = 0;

#ifndef FFT_COLUMN_BASED
        int task0 = myplan.slab_to_task[slab_x];
        int task1 = myplan.slab_to_task[slab_xx];

        size_t ind0 = send_offset[task0] + send_count[task0]++;
        for(int j = 0; j < 3; j++)
          partout[ind0].IntPos[j] = Sp->P[i].IntPos[j];

        if(task0 != task1)
          {
            size_t ind1 = send_offset[task1] + send_count[task1]++;
            for(int j = 0; j < 3; j++)
              partout[ind1].IntPos[j] = Sp->P[i].IntPos[j];
          }
#else
        int slab_y  = Sp->P[i].IntPos[1] / INTCELL;
        int slab_yy = slab_y + 1;

        if(slab_yy >= GRIDY)
          slab_yy = 0;

        int column0 = slab_x * GRIDY + slab_y;
        int column1 = slab_x * GRIDY + slab_yy;
        int column2 = slab_xx * GRIDY + slab_y;
        int column3 = slab_xx * GRIDY + slab_yy;

        int task0, task1, task2, task3;

        if(column0 < pivotcol)
          task0 = column0 / avg;
        else
          task0 = (column0 - pivotcol) / (avg - 1) + tasklastsection;

        if(column1 < pivotcol)
          task1 = column1 / avg;
        else
          task1 = (column1 - pivotcol) / (avg - 1) + tasklastsection;

        if(column2 < pivotcol)
          task2 = column2 / avg;
        else
          task2 = (column2 - pivotcol) / (avg - 1) + tasklastsection;

        if(column3 < pivotcol)
          task3 = column3 / avg;
        else
          task3 = (column3 - pivotcol) / (avg - 1) + tasklastsection;

        size_t ind0 = send_offset[task0] + send_count[task0]++;
        for(int j = 0; j < 3; j++)
          partout[ind0].IntPos[j] = Sp->P[i].IntPos[j];

        if(task1 != task0)
          {
            size_t ind1 = send_offset[task1] + send_count[task1]++;

            for(int j = 0; j < 3; j++)
              partout[ind1].IntPos[j] = Sp->P[i].IntPos[j];
          }
        if(task2 != task1 && task2 != task0)
          {
            size_t ind2 = send_offset[task2] + send_count[task2]++;

            for(int j = 0; j < 3; j++)
              partout[ind2].IntPos[j] = Sp->P[i].IntPos[j];
          }
        if(task3 != task0 && task3 != task1 && task3 != task2)
          {
            size_t ind3 = send_offset[task3] + send_count[task3]++;

            for(int j = 0; j < 3; j++)
              partout[ind3].IntPos[j] = Sp->P[i].IntPos[j];
          }
#endif
      }
  }

  int flag_big = 0, flag_big_all;
  for(int i = 0; i < NTask; i++)
    if(Sndpm_count[i] * sizeof(partbuf) > MPI_MESSAGE_SIZELIMIT_IN_BYTES)
      flag_big = 1;

  /* produce a flag if any of the send sizes is above our transfer limit, in this case we will
   * transfer the data in chunks.
   */
  MPI_Allreduce(&flag_big, &flag_big_all, 1, MPI_INT, MPI_MAX, Communicator);

  /* exchange particle data */
  myMPI_Alltoallv(partout, Sndpm_count, Sndpm_offset, partin, Rcvpm_count, Rcvpm_offset, sizeof(partbuf), flag_big_all, Communicator);

  Mem.myfree(partout);
}

void ngenic::ngenic_compute_transform_of_source_potential(fft_real *pot)
{
  fft_real *workspace = (fft_real *)Mem.mymalloc("workspace", maxfftsize * sizeof(fft_real));

#ifndef FFT_COLUMN_BASED
  my_slab_based_fft(&myplan, &pot[0], &workspace[0], +1);
#else
  my_column_based_fft(&myplan, pot, workspace, +1);  // result is in workspace, not in Phi2
  memcpy(pot, workspace, maxfftsize * sizeof(fft_real));
#endif

  Mem.myfree(workspace);

  double normfac = 1 / (((double)GRIDX) * GRIDY * GRIDZ);

  for(size_t n = 0; n < maxfftsize; n++)
    pot[n] *= normfac;
}

/* this function returns the component 'axes' (0, 1, or 2) of the gradient of a field phi,
 * which is the solution of nabla^2 phi = grid.
 * We input the Fourier transform of grid to the function, and this field is overwritten with
 * the gradient.
 */
void ngenic::ngenic_get_derivate_from_fourier_field(int axes1, int axes2, fft_complex *fft_of_grid)
{
  double kfacx = 2.0 * M_PI / All.BoxSize;
  double kfacy = 2.0 * M_PI / All.BoxSize;
  double kfacz = 2.0 * M_PI / All.BoxSize;

#ifdef FFT_COLUMN_BASED
  for(large_array_offset ip = 0; ip < myplan.second_transposed_ncells; ip += GRIDX)
    {
      large_array_offset ipcell = ip + ((large_array_offset)myplan.second_transposed_firstcol) * GRIDX;
      int y                     = ipcell / (GRIDX * GRIDz);
      int yr                    = ipcell % (GRIDX * GRIDz);
      int z                     = yr / GRIDX;
      if(yr % GRIDX != 0)  // Note: check that x-columns are really complete
        Terminate("x-column seems incomplete. This is not expected");
#else
  for(int y = myplan.slabstart_y; y < myplan.slabstart_y + myplan.nslab_y; y++)
    for(int z = 0; z < GRIDz; z++)
      {
#endif

      for(int x = 0; x < GRIDX; x++)
        {
          int xx, yy, zz;

          if(x >= (GRIDX / 2))
            xx = x - GRIDX;
          else
            xx = x;
          if(y >= (GRIDY / 2))
            yy = y - GRIDY;
          else
            yy = y;
          if(z >= (GRIDZ / 2))
            zz = z - GRIDZ;
          else
            zz = z;

          double kvec[3];
          kvec[0] = kfacx * xx;
          kvec[1] = kfacy * yy;
          kvec[2] = kfacz * zz;

          double kmag2 = kvec[0] * kvec[0] + kvec[1] * kvec[1] + kvec[2] * kvec[2];

          double smth = 1;

#ifdef CORRECT_CIC
          if(axes2 >= 0)
            {
              /* do deconvolution of CIC interpolation */
              double fx = 1, fy = 1, fz = 1;
              if(kvec[0] != 0)
                {
                  fx = (kvec[0] * All.BoxSize / 2) / NGENIC;
                  fx = sin(fx) / fx;
                }
              if(kvec[1] != 0)
                {
                  fy = (kvec[1] * All.BoxSize / 2) / NGENIC;
                  fy = sin(fy) / fy;
                }
              if(kvec[2] != 0)
                {
                  fz = (kvec[2] * All.BoxSize / 2) / NGENIC;
                  fz = sin(fz) / fz;
                }
              double ff = 1 / (fx * fy * fz);
              smth      = ff * ff;
              /* end deconvolution */
            }
#endif

#ifndef FFT_COLUMN_BASED
          large_array_offset elem = ((large_array_offset)GRIDz) * (GRIDX * (y - myplan.slabstart_y) + x) + z;
#else
            large_array_offset elem = ip + x;
#endif

          fft_real re = smth * fft_of_grid[elem][0];
          fft_real im = smth * fft_of_grid[elem][1];

          if(axes2 < 0)
            {
              /* first derivative */
              fft_of_grid[elem][0] = (kmag2 > 0.0 ? -kvec[axes1] / kmag2 * im : 0.0);
              fft_of_grid[elem][1] = (kmag2 > 0.0 ? kvec[axes1] / kmag2 * re : 0.0);
            }
          else
            {
              /* second derivative */
              fft_of_grid[elem][0] = (kmag2 > 0.0 ? kvec[axes1] * kvec[axes2] / kmag2 * re : 0.0);
              fft_of_grid[elem][1] = (kmag2 > 0.0 ? kvec[axes1] * kvec[axes2] / kmag2 * im : 0.0);
            }
        }
    }

#ifdef FFT_COLUMN_BASED
  if(myplan.second_transposed_firstcol == 0)
    fft_of_grid[0][0] = fft_of_grid[0][1] = 0.0;
#else
  if(myplan.slabstart_y == 0)
    fft_of_grid[0][0] = fft_of_grid[0][1] = 0.0;
#endif

  /* Do the inverse FFT to get the displacement field */
  fft_real *workspace = (fft_real *)Mem.mymalloc("workspace", maxfftsize * sizeof(fft_real));

#ifndef FFT_COLUMN_BASED
  my_slab_based_fft(&myplan, &fft_of_grid[0], &workspace[0], -1);
#else
  my_column_based_fft(&myplan, fft_of_grid, workspace, -1);  // result is in workspace
  memcpy(fft_of_grid, workspace, maxfftsize * sizeof(fft_real));
#endif

  Mem.myfree(workspace);
}

void ngenic::ngenic_setup_modes_in_kspace(fft_complex *fft_of_grid)
{
  double fac = pow(2 * M_PI / All.BoxSize, 1.5);

  /* clear local FFT-mesh */
  memset(fft_of_grid, 0, maxfftsize * sizeof(fft_real));

  mpi_printf("NGENIC: setting up modes in kspace...\n");

  double kfacx = 2.0 * M_PI / All.BoxSize;
  double kfacy = 2.0 * M_PI / All.BoxSize;
  double kfacz = 2.0 * M_PI / All.BoxSize;

#ifdef FFT_COLUMN_BASED
  for(large_array_offset ip = 0; ip < myplan.second_transposed_ncells; ip += GRIDX)
    {
      large_array_offset ipcell = ip + ((large_array_offset)myplan.second_transposed_firstcol) * GRIDX;
      int y                     = ipcell / (GRIDX * GRIDz);
      int yr                    = ipcell % (GRIDX * GRIDz);
      int z                     = yr / GRIDX;
      if(yr % GRIDX != 0)  // Note: check that x-columns are really complete
        Terminate("x-column seems incomplete. This is not expected");
#else
  for(int y = myplan.slabstart_y; y < myplan.slabstart_y + myplan.nslab_y; y++)
    for(int z = 0; z < GRIDz; z++)
      {
#endif

      // let's use the y and z plane here, because the x-column is available in full for both FFT schemes
      gsl_rng_set(rnd_generator, seedtable[y * NGENIC + z]);

      // we also create the modes for the conjugate column so that we can fulfill the reality constraint
      // by using the conjugate of the corresponding mode if needed
      int y_conj, z_conj;
      if(y > 0)
        y_conj = GRIDY - y;
      else
        y_conj = 0;

      if(z > 0)
        z_conj = GRIDZ - z;
      else
        z_conj = 0;

      gsl_rng_set(rnd_generator_conjugate, seedtable[y_conj * NGENIC + z_conj]);

#ifndef NGENIC_FIX_MODE_AMPLITUDES
      double mode_ampl[GRIDX], mode_ampl_conj[GRIDX];
#endif
      double mode_phase[GRIDX], mode_phase_conj[GRIDX];

      // in this loop we precompute the modes for both columns, from low-k to high-k,
      // so that after an increase of resolution, one gets the same modes plus new ones
      for(int xoff = 0; xoff < GRIDX / 2; xoff++)
        for(int side = 0; side < 2; side++)
          {
            int x;
            if(side == 0)
              x = xoff;
            else
              x = GRIDX - 1 - xoff;

            double phase      = gsl_rng_uniform(rnd_generator) * 2 * M_PI;
            double phase_conj = gsl_rng_uniform(rnd_generator_conjugate) * 2 * M_PI;

#ifdef NGENIC_MIRROR_PHASES
            phase += M_PI;
            if(phase >= 2 * M_PI)
              phase -= 2 * M_PI;

            phase_conj += M_PI;
            if(phase_conj >= 2 * M_PI)
              phase_conj -= 2 * M_PI;
#endif
            mode_phase[x]      = phase;
            mode_phase_conj[x] = phase_conj;

#ifndef NGENIC_FIX_MODE_AMPLITUDES
            double ampl;
            do
              {
                ampl = gsl_rng_uniform(rnd_generator);
              }
            while(ampl == 0);

            double ampl_conj;
            do
              {
                ampl_conj = gsl_rng_uniform(rnd_generator_conjugate);
              }
            while(ampl_conj == 0);

            mode_ampl[x] = ampl;

            mode_ampl_conj[x] = ampl_conj;
#endif
          }

      // now let's populate the full x-column of modes
      for(int x = 0; x < GRIDX; x++)
        {
          int xx, yy, zz;

          if(x >= (GRIDX / 2))
            xx = x - GRIDX;
          else
            xx = x;
          if(y >= (GRIDY / 2))
            yy = y - GRIDY;
          else
            yy = y;
          if(z >= (GRIDZ / 2))
            zz = z - GRIDZ;
          else
            zz = z;

          double kvec[3];
          kvec[0] = kfacx * xx;
          kvec[1] = kfacy * yy;
          kvec[2] = kfacz * zz;

          double kmag2 = kvec[0] * kvec[0] + kvec[1] * kvec[1] + kvec[2] * kvec[2];
          double kmag  = sqrt(kmag2);

          if(All.SphereMode == 1)
            {
              if(kmag * All.BoxSize / (2 * M_PI) > All.NSample / 2) /* select a sphere in k-space */
                continue;
            }
          else
            {
              if(fabs(kvec[0]) * All.BoxSize / (2 * M_PI) > All.NSample / 2)
                continue;
              if(fabs(kvec[1]) * All.BoxSize / (2 * M_PI) > All.NSample / 2)
                continue;
              if(fabs(kvec[2]) * All.BoxSize / (2 * M_PI) > All.NSample / 2)
                continue;
            }

          double p_of_k = ngenic_power_spec(kmag);

          /* Note: kmag and p_of_k are unaffected by whether or not we use the conjugate mode */

          int conjugate_flag = 0;

          if(z == 0 || z == GRIDZ / 2)
            {
              if(x > GRIDX / 2 && x < GRIDX)
                conjugate_flag = 1;
              else if(x == 0 || x == GRIDX / 2)
                {
                  if(y > GRIDY / 2 && y < GRIDY)
                    conjugate_flag = 1;
                  else if(y == 0 || y == GRIDX / 2)
                    {
                      continue;
                    }
                }
            }

          // determine location of conjugate mode in x column
          int x_conj;

          if(x > 0)
            x_conj = GRIDX - x;
          else
            x_conj = 0;

#ifndef NGENIC_FIX_MODE_AMPLITUDES
          if(conjugate_flag)
            p_of_k *= -log(mode_ampl_conj[x_conj]);
          else
            p_of_k *= -log(mode_ampl[x]);
#endif

          double delta = fac * sqrt(p_of_k) / Dplus; /* scale back to starting redshift */

#ifndef FFT_COLUMN_BASED
          large_array_offset elem = ((large_array_offset)GRIDz) * (GRIDX * (y - myplan.slabstart_y) + x) + z;
#else
            large_array_offset elem = ip + x;
#endif

          if(conjugate_flag)
            {
              fft_of_grid[elem][0] = delta * cos(mode_phase_conj[x_conj]);
              fft_of_grid[elem][1] = -delta * sin(mode_phase_conj[x_conj]);
            }
          else
            {
              fft_of_grid[elem][0] = delta * cos(mode_phase[x]);
              fft_of_grid[elem][1] = delta * sin(mode_phase[x]);
            }
        }
    }
}

void ngenic::ngenic_readout_disp(fft_real *grid, int axis, double pfac, double vfac)
{
#ifdef FFT_COLUMN_BASED
  int columns         = GRIDX * GRIDY;
  int avg             = (columns - 1) / NTask + 1;
  int exc             = NTask * avg - columns;
  int tasklastsection = NTask - exc;
  int pivotcol        = tasklastsection * avg;
#endif

  double *flistin  = (double *)Mem.mymalloc("flistin", nimport * sizeof(double));
  double *flistout = (double *)Mem.mymalloc("flistout", nexport * sizeof(double));

  for(size_t i = 0; i < nimport; i++)
    {
      flistin[i] = 0;

      int slab_x         = partin[i].IntPos[0] / INTCELL;
      MyIntPosType rmd_x = partin[i].IntPos[0] % INTCELL;

      int slab_y         = partin[i].IntPos[1] / INTCELL;
      MyIntPosType rmd_y = partin[i].IntPos[1] % INTCELL;

      int slab_z         = partin[i].IntPos[2] / INTCELL;
      MyIntPosType rmd_z = partin[i].IntPos[2] % INTCELL;

      double dx = rmd_x * (1.0 / INTCELL);
      double dy = rmd_y * (1.0 / INTCELL);
      double dz = rmd_z * (1.0 / INTCELL);

      int slab_xx = slab_x + 1;
      int slab_yy = slab_y + 1;
      int slab_zz = slab_z + 1;

      if(slab_xx >= GRIDX)
        slab_xx = 0;
      if(slab_yy >= GRIDY)
        slab_yy = 0;
      if(slab_zz >= GRIDZ)
        slab_zz = 0;

#ifndef FFT_COLUMN_BASED
      if(myplan.slab_to_task[slab_x] == ThisTask)
        {
          slab_x -= myplan.first_slab_x_of_task[ThisTask];

          flistin[i] += grid[FI(slab_x, slab_y, slab_z)] * (1.0 - dx) * (1.0 - dy) * (1.0 - dz) +
                        grid[FI(slab_x, slab_y, slab_zz)] * (1.0 - dx) * (1.0 - dy) * (dz) +
                        grid[FI(slab_x, slab_yy, slab_z)] * (1.0 - dx) * (dy) * (1.0 - dz) +
                        grid[FI(slab_x, slab_yy, slab_zz)] * (1.0 - dx) * (dy) * (dz);
        }

      if(myplan.slab_to_task[slab_xx] == ThisTask)
        {
          slab_xx -= myplan.first_slab_x_of_task[ThisTask];

          flistin[i] += grid[FI(slab_xx, slab_y, slab_z)] * (dx) * (1.0 - dy) * (1.0 - dz) +
                        grid[FI(slab_xx, slab_y, slab_zz)] * (dx) * (1.0 - dy) * (dz) +
                        grid[FI(slab_xx, slab_yy, slab_z)] * (dx) * (dy) * (1.0 - dz) +
                        grid[FI(slab_xx, slab_yy, slab_zz)] * (dx) * (dy) * (dz);
        }
#else
      int column0 = slab_x * GRIDY + slab_y;
      int column1 = slab_x * GRIDY + slab_yy;
      int column2 = slab_xx * GRIDY + slab_y;
      int column3 = slab_xx * GRIDY + slab_yy;

      if(column0 >= myplan.firstcol_XY && column0 <= myplan.lastcol_XY)
        {
          flistin[i] += grid[FC(column0, slab_z)] * (1.0 - dx) * (1.0 - dy) * (1.0 - dz) +
                        grid[FC(column0, slab_zz)] * (1.0 - dx) * (1.0 - dy) * (dz);
        }
      if(column1 >= myplan.firstcol_XY && column1 <= myplan.lastcol_XY)
        {
          flistin[i] +=
              grid[FC(column1, slab_z)] * (1.0 - dx) * (dy) * (1.0 - dz) + grid[FC(column1, slab_zz)] * (1.0 - dx) * (dy) * (dz);
        }

      if(column2 >= myplan.firstcol_XY && column2 <= myplan.lastcol_XY)
        {
          flistin[i] +=
              grid[FC(column2, slab_z)] * (dx) * (1.0 - dy) * (1.0 - dz) + grid[FC(column2, slab_zz)] * (dx) * (1.0 - dy) * (dz);
        }

      if(column3 >= myplan.firstcol_XY && column3 <= myplan.lastcol_XY)
        {
          flistin[i] += grid[FC(column3, slab_z)] * (dx) * (dy) * (1.0 - dz) + grid[FC(column3, slab_zz)] * (dx) * (dy) * (dz);
        }
#endif
    }

  /* exchange the potential component data */
  int flag_big = 0, flag_big_all;
  for(int i = 0; i < NTask; i++)
    if(Sndpm_count[i] * sizeof(double) > MPI_MESSAGE_SIZELIMIT_IN_BYTES)
      flag_big = 1;

  /* produce a flag if any of the send sizes is above our transfer limit, in this case we will
   * transfer the data in chunks.
   */
  MPI_Allreduce(&flag_big, &flag_big_all, 1, MPI_INT, MPI_MAX, Communicator);

  /* exchange  data */
  myMPI_Alltoallv(flistin, Rcvpm_count, Rcvpm_offset, flistout, Sndpm_count, Sndpm_offset, sizeof(double), flag_big_all, Communicator);

  {
    size_t *send_count  = Sndpm_count;
    size_t *send_offset = Sndpm_offset;

    for(int j = 0; j < NTask; j++)
      send_count[j] = 0;

    for(int i = 0; i < Sp->NumPart; i++)
      {
        int slab_x  = Sp->P[i].IntPos[0] / INTCELL;
        int slab_xx = slab_x + 1;

        if(slab_xx >= GRIDX)
          slab_xx = 0;

#ifndef FFT_COLUMN_BASED
        int task0 = myplan.slab_to_task[slab_x];
        int task1 = myplan.slab_to_task[slab_xx];

        double value = flistout[send_offset[task0] + send_count[task0]++];

        if(task0 != task1)
          value += flistout[send_offset[task1] + send_count[task1]++];
#else
        int slab_y  = Sp->P[i].IntPos[1] / INTCELL;
        int slab_yy = slab_y + 1;

        if(slab_yy >= GRIDY)
          slab_yy = 0;

        int column0 = slab_x * GRIDY + slab_y;
        int column1 = slab_x * GRIDY + slab_yy;
        int column2 = slab_xx * GRIDY + slab_y;
        int column3 = slab_xx * GRIDY + slab_yy;

        int task0, task1, task2, task3;

        if(column0 < pivotcol)
          task0 = column0 / avg;
        else
          task0 = (column0 - pivotcol) / (avg - 1) + tasklastsection;

        if(column1 < pivotcol)
          task1 = column1 / avg;
        else
          task1 = (column1 - pivotcol) / (avg - 1) + tasklastsection;

        if(column2 < pivotcol)
          task2 = column2 / avg;
        else
          task2 = (column2 - pivotcol) / (avg - 1) + tasklastsection;

        if(column3 < pivotcol)
          task3 = column3 / avg;
        else
          task3 = (column3 - pivotcol) / (avg - 1) + tasklastsection;

        double value = flistout[send_offset[task0] + send_count[task0]++];

        if(task1 != task0)
          value += flistout[send_offset[task1] + send_count[task1]++];

        if(task2 != task1 && task2 != task0)
          value += flistout[send_offset[task2] + send_count[task2]++];

        if(task3 != task0 && task3 != task1 && task3 != task2)
          value += flistout[send_offset[task3] + send_count[task3]++];
#endif

        Pdisp[i].deltapos[axis] += pfac * value;
        Sp->P[i].Vel[axis] += vfac * value;
      }
  }

  Mem.myfree(flistout);
  Mem.myfree(flistin);
}

void ngenic::ngenic_initialize_ffts(void)
{
#ifdef LONG_X
  if(LONG_X != (int)(LONG_X))
    Terminate("LONG_X must be an integer if used with PMGRID");
#endif

#ifdef LONG_Y
  if(LONG_Y != (int)(LONG_Y))
    Terminate("LONG_Y must be an integer if used with PMGRID");
#endif

#ifdef LONG_Z
  if(LONG_Z != (int)(LONG_Z))
    Terminate("LONG_Z must be an integer if used with PMGRID");
#endif

  /* Set up the FFTW-3 plan files. */
  int ndimx[1] = {GRIDX}; /* dimension of the 1D transforms */
  int ndimy[1] = {GRIDY}; /* dimension of the 1D transforms */
  int ndimz[1] = {GRIDZ}; /* dimension of the 1D transforms */

  int max_GRID2 = 2 * (std::max<int>(std::max<int>(GRIDX, GRIDY), GRIDZ) / 2 + 1);

  /* temporarily allocate some arrays to make sure that out-of-place plans are created */
  fft_real *DispGrid           = (fft_real *)Mem.mymalloc("DispGrid", max_GRID2 * sizeof(fft_real));
  fft_complex *fft_of_DispGrid = (fft_complex *)Mem.mymalloc("DispGrid", max_GRID2 * sizeof(fft_real));

#ifdef DOUBLEPRECISION_FFTW
  int alignflag = 0;
#else
  /* for single precision, the start of our FFT columns is presently only guaranteed to be 8-byte aligned */
  int alignflag = FFTW_UNALIGNED;
#endif

#ifndef FFT_COLUMN_BASED
  int stride = GRIDz;
#else
  int stride    = 1;
#endif

  myplan.backward_plan_xdir =
      FFTW(plan_many_dft)(1, ndimx, 1, (fft_complex *)DispGrid, 0, stride, GRIDz * GRIDX, fft_of_DispGrid, 0, stride, GRIDz * GRIDX,
                          FFTW_BACKWARD, FFTW_ESTIMATE | FFTW_DESTROY_INPUT | alignflag);

  myplan.backward_plan_ydir =
      FFTW(plan_many_dft)(1, ndimy, 1, (fft_complex *)DispGrid, 0, stride, GRIDz * GRIDY, fft_of_DispGrid, 0, stride, GRIDz * GRIDY,
                          FFTW_BACKWARD, FFTW_ESTIMATE | FFTW_DESTROY_INPUT | alignflag);

  myplan.backward_plan_zdir = FFTW(plan_many_dft_c2r)(1, ndimz, 1, (fft_complex *)DispGrid, 0, 1, GRIDz, (fft_real *)fft_of_DispGrid,
                                                      0, 1, GRID2, FFTW_ESTIMATE | FFTW_DESTROY_INPUT | alignflag);

  myplan.forward_plan_xdir = FFTW(plan_many_dft)(1, ndimx, 1, (fft_complex *)DispGrid, 0, stride, GRIDz * GRIDX, fft_of_DispGrid, 0,
                                                 stride, GRIDz * GRIDX, FFTW_FORWARD, FFTW_ESTIMATE | FFTW_DESTROY_INPUT | alignflag);

  myplan.forward_plan_ydir = FFTW(plan_many_dft)(1, ndimy, 1, (fft_complex *)DispGrid, 0, stride, GRIDz * GRIDY, fft_of_DispGrid, 0,
                                                 stride, GRIDz * GRIDY, FFTW_FORWARD, FFTW_ESTIMATE | FFTW_DESTROY_INPUT | alignflag);

  myplan.forward_plan_zdir = FFTW(plan_many_dft_r2c)(1, ndimz, 1, DispGrid, 0, 1, GRID2, (fft_complex *)fft_of_DispGrid, 0, 1, GRIDz,
                                                     FFTW_ESTIMATE | FFTW_DESTROY_INPUT | alignflag);

  Mem.myfree(fft_of_DispGrid);
  Mem.myfree(DispGrid);

#ifndef FFT_COLUMN_BASED

  my_slab_based_fft_init(&myplan, GRIDX, GRIDY, GRIDZ);

  maxfftsize = std::max<int>(myplan.largest_x_slab * GRIDY, myplan.largest_y_slab * GRIDX) * ((size_t)GRID2);

#else

  my_column_based_fft_init(&myplan, GRIDX, GRIDY, GRIDZ);

  maxfftsize = myplan.max_datasize;

#endif
}

void ngenic::print_spec(void)
{
  if(ThisTask == 0)
    {
      char buf[3 * MAXLEN_PATH];
      sprintf(buf, "%s/inputspec_%s.txt", All.OutputDir, All.SnapshotFileBase);

      FILE *fd = fopen(buf, "w");

      double gf = ngenic_growth_factor(0.001, 1.0) / (1.0 / 0.001);

      double DDD = ngenic_growth_factor(All.cf_atime, 1.0);

      fprintf(fd, "%12g %12g\n", All.cf_redshift, DDD); /* print actual starting redshift and
                                                           linear growth factor for this cosmology */

      double kstart = 2 * M_PI / (1000.0 * (1e6 * PARSEC / All.UnitLength_in_cm));  /* 1000 Mpc/h */
      double kend   = 2 * M_PI / (0.001 * (3.1e6 * PARSEC / All.UnitLength_in_cm)); /* 0.001 Mpc/h */

      for(double k = kstart; k < kend; k *= 1.025)
        {
          double po = ngenic_power_spec(k);
          double dl = 4.0 * M_PI * k * k * k * po;

          double kf = 0.5;

          double po2 = ngenic_power_spec(1.001 * k * kf);
          double po1 = ngenic_power_spec(k * kf);
          double dnl = 0, knl = 0;

          if(po != 0 && po1 != 0 && po2 != 0)
            {
              double neff = (log(po2) - log(po1)) / (log(1.001 * k * kf) - log(k * kf));

              if(1 + neff / 3 > 0)
                {
                  double A     = 0.482 * pow(1 + neff / 3, -0.947);
                  double B     = 0.226 * pow(1 + neff / 3, -1.778);
                  double alpha = 3.310 * pow(1 + neff / 3, -0.244);
                  double beta  = 0.862 * pow(1 + neff / 3, -0.287);
                  double V     = 11.55 * pow(1 + neff / 3, -0.423) * 1.2;

                  dnl = fnl(dl, A, B, alpha, beta, V, gf);
                  knl = k * pow(1 + dnl, 1.0 / 3);
                }
            }

          fprintf(fd, "%12g %12g    %12g %12g\n", k, dl, knl, dnl);
        }
      fclose(fd);
    }
1189
1190
1191
1192
1193
1194
1195
1196
1197
1198
1199
1200
1201
1202
1203
1204
1205
1206
1207
1208
1209
1210
1211
1212
1213
1214
1215
1216
1217
1218
1219
1220
1221
1222
1223
1224
1225
1226
1227
1228
1229
1230
1231
1232
1233
1234
1235
1236
1237
1238
1239
1240
1241

  /* create an output file with the growth factor as a function of redshift, which
   * can be handy in analysis routines
   */
  if(ThisTask == 0)
    {
      if(All.cf_atime < 1.0)
        {
          char buf[3 * MAXLEN_PATH];
          sprintf(buf, "%s/growthfac.txt", All.OutputDir);

          FILE *fd = fopen(buf, "w");

          const int NSTEPS = 100;

          for(int i = 0; i <= NSTEPS; i++)
            {
              double a = exp(log(All.cf_atime) + ((log(1.0) - log(All.cf_atime)) / NSTEPS) * i);

              double d = ngenic_growth_factor(a, 1.0);

              fprintf(fd, "%12g %12g\n", a, 1.0 / d);
            }
          fclose(fd);
        }
    }

  /* also create an output file with the sigma(M), i.e. the variance as a function of the mass scale,
   * for simplifying analytic mass function plots, such as Press-Schechter
   */

  if(ThisTask == 0)
    {
      char buf[3 * MAXLEN_PATH];
      sprintf(buf, "%s/variance.txt", All.OutputDir);

      FILE *fd = fopen(buf, "w");

      double rhoback = All.Omega0 * 3 * All.Hubble * All.Hubble / (8 * M_PI * All.G);

      for(double M = 1.0e5; M <= 1.01e16; M *= pow(10.0, 1.0 / 16))  // mass in solar masses / h
        {
          double Mint = M * (SOLAR_MASS / All.UnitMass_in_g);

          double R = pow(3.0 * Mint / (4 * M_PI * rhoback), 1.0 / 3);

          double sigma2 = ngenic_tophat_sigma2(R);

          fprintf(fd, "%12g %12g %12g %12g\n", M, Mint, sigma2, sqrt(sigma2));
        }

      fclose(fd);
    }
Volker Springel's avatar
Volker Springel committed
1242
1243
1244
}

#endif