-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathpoisson.cc
1197 lines (955 loc) · 34.6 KB
/
poisson.cc
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
/*
poisson.cc - This file is part of MUSIC -
a code to generate multi-scale initial conditions
for cosmological simulations
Copyright (C) 2010 Oliver Hahn
*/
/****** ABSTRACT FACTORY PATTERN IMPLEMENTATION *******/
#include "poisson.hh"
#include "Numerics.hh"
std::map< std::string, poisson_plugin_creator *>&
get_poisson_plugin_map()
{
static std::map< std::string, poisson_plugin_creator* > poisson_plugin_map;
return poisson_plugin_map;
}
void print_poisson_plugins()
{
std::map< std::string, poisson_plugin_creator *>& m = get_poisson_plugin_map();
std::map< std::string, poisson_plugin_creator *>::iterator it;
it = m.begin();
std::cout << " - Available poisson solver plug-ins:\n";
while( it!=m.end() )
{
if( (*it).second )
std::cout << "\t\'" << (*it).first << "\'\n";
LOGINFO("Poisson plug-in :: %s",std::string((*it).first).c_str());
++it;
}
}
/****** CALL IMPLEMENTATIONS OF POISSON SOLVER CLASSES ******/
#include "mg_solver.hh"
#include "fd_schemes.hh"
#ifdef SINGLE_PRECISION
typedef multigrid::solver< stencil_7P<float>, interp_O3_fluxcorr, mg_straight, float > poisson_solver_O2;
typedef multigrid::solver< stencil_13P<float>, interp_O5_fluxcorr, mg_straight, float > poisson_solver_O4;
typedef multigrid::solver< stencil_19P<float>, interp_O7_fluxcorr, mg_straight, float > poisson_solver_O6;
#else
typedef multigrid::solver< stencil_7P<double>, interp_O3_fluxcorr, mg_straight, double > poisson_solver_O2;
typedef multigrid::solver< stencil_13P<double>, interp_O5_fluxcorr, mg_straight, double > poisson_solver_O4;
typedef multigrid::solver< stencil_19P<double>, interp_O7_fluxcorr, mg_straight, double > poisson_solver_O6;
#endif
/**************************************************************************************/
/**************************************************************************************/
#pragma mark -
double multigrid_poisson_plugin::solve( grid_hierarchy& f, grid_hierarchy& u )
{
LOGUSER("Initializing multi-grid Poisson solver...");
unsigned verbosity = cf_.getValueSafe<unsigned>("setup","verbosity",2);
if( verbosity > 0 )
{
std::cout << "-------------------------------------------------------------\n";
std::cout << " - Invoking multi-grid Poisson solver..." << std::endl;
}
double acc = 1e-5, err;
std::string ps_smoother_name;
unsigned ps_presmooth, ps_postsmooth, order;
acc = cf_.getValueSafe<double>("poisson","accuracy",acc);
ps_presmooth = cf_.getValueSafe<unsigned>("poisson","pre_smooth",3);
ps_postsmooth = cf_.getValueSafe<unsigned>("poisson","post_smooth",3);
ps_smoother_name = cf_.getValueSafe<std::string>("poisson","smoother","gs");
order = cf_.getValueSafe<unsigned>( "poisson", "laplace_order", 4 );
multigrid::opt::smtype ps_smtype = multigrid::opt::sm_gauss_seidel;
if ( ps_smoother_name == std::string("gs") )
{
ps_smtype = multigrid::opt::sm_gauss_seidel;
LOGUSER("Selected Gauss-Seidel multigrid smoother");
}
else if ( ps_smoother_name == std::string("jacobi") )
{
ps_smtype = multigrid::opt::sm_jacobi;
LOGUSER("Selected Jacobi multigrid smoother");
}
else if ( ps_smoother_name == std::string("sor") )
{
ps_smtype = multigrid::opt::sm_sor;
LOGUSER("Selected SOR multigrid smoother");
}
else
{
LOGWARN("Unknown multigrid smoother \'%s\' specified. Reverting to Gauss-Seidel.",ps_smoother_name.c_str());
std::cerr << " - Warning: unknown smoother \'" << ps_smoother_name << "\' for multigrid solver!\n"
<< " reverting to \'gs\' (Gauss-Seidel)" << std::endl;
}
double tstart, tend;
#ifdef _OPENMP
tstart = omp_get_wtime();
#else
tstart = (double)clock() / CLOCKS_PER_SEC;
#endif
//----- run Poisson solver -----//
if( order == 2 )
{
LOGUSER("Running multigrid solver with 2nd order Laplacian...");
poisson_solver_O2 ps( f, ps_smtype, ps_presmooth, ps_postsmooth );
err = ps.solve( u, acc, true );
}
else if( order == 4 )
{
LOGUSER("Running multigrid solver with 4th order Laplacian...");
poisson_solver_O4 ps( f, ps_smtype, ps_presmooth, ps_postsmooth );
err = ps.solve( u, acc, true );
}
else if( order == 6 )
{
LOGUSER("Running multigrid solver with 6th order Laplacian..");
poisson_solver_O6 ps( f, ps_smtype, ps_presmooth, ps_postsmooth );
err = ps.solve( u, acc, true );
}
else
{
LOGERR("Invalid order specified for Laplace operator");
throw std::runtime_error("Invalid order specified for Laplace operator");
}
//------------------------------//
#ifdef _OPENMP
tend = omp_get_wtime();
if( verbosity > 1 )
std::cout << " - Poisson solver took " << tend-tstart << "s with " << omp_get_max_threads() << " threads." << std::endl;
#else
tend = (double)clock() / CLOCKS_PER_SEC;
if( verbosity > 1 )
std::cout << " - Poisson solver took " << tend-tstart << "s." << std::endl;
#endif
return err;
}
double multigrid_poisson_plugin::gradient( int dir, grid_hierarchy& u, grid_hierarchy& Du )
{
Du = u;
unsigned order = cf_.getValueSafe<unsigned>( "poisson", "grad_order", 4 );
switch( order )
{
case 2:
implementation().gradient_O2( dir, u, Du );
break;
case 4:
implementation().gradient_O4( dir, u, Du );
break;
case 6:
implementation().gradient_O6( dir, u, Du );
break;
default:
LOGERR("Invalid order %d specified for gradient operator", order);
throw std::runtime_error("Invalid order specified for gradient operator!");
}
return 0.0;
}
double multigrid_poisson_plugin::gradient_add( int dir, grid_hierarchy& u, grid_hierarchy& Du )
{
//Du = u;
unsigned order = cf_.getValueSafe<unsigned>( "poisson", "grad_order", 4 );
switch( order )
{
case 2:
implementation().gradient_add_O2( dir, u, Du );
break;
case 4:
implementation().gradient_add_O4( dir, u, Du );
break;
case 6:
implementation().gradient_add_O6( dir, u, Du );
break;
default:
LOGERR("Invalid order %d specified for gradient operator!",order);
throw std::runtime_error("Invalid order specified for gradient operator!");
}
return 0.0;
}
void multigrid_poisson_plugin::implementation::gradient_O2( int dir, grid_hierarchy& u, grid_hierarchy& Du )
{
LOGUSER("Computing a 2nd order finite difference gradient...");
for( unsigned ilevel=u.levelmin(); ilevel<=u.levelmax(); ++ilevel )
{
double h = pow(2.0,ilevel);
meshvar_bnd *pvar = Du.get_grid(ilevel);
if( dir == 0 )
#pragma omp parallel for
for( int ix = 0; ix < (int)(*u.get_grid(ilevel)).size(0); ++ix )
for( int iy = 0; iy < (int)(*u.get_grid(ilevel)).size(1); ++iy )
for( int iz = 0; iz < (int)(*u.get_grid(ilevel)).size(2); ++iz )
(*pvar)(ix,iy,iz) = 0.5*((*u.get_grid(ilevel))(ix+1,iy,iz)-(*u.get_grid(ilevel))(ix-1,iy,iz))*h;
else if( dir == 1 )
#pragma omp parallel for
for( int ix = 0; ix < (int)(*u.get_grid(ilevel)).size(0); ++ix )
for( int iy = 0; iy < (int)(*u.get_grid(ilevel)).size(1); ++iy )
for( int iz = 0; iz < (int)(*u.get_grid(ilevel)).size(2); ++iz )
(*pvar)(ix,iy,iz) = 0.5*((*u.get_grid(ilevel))(ix,iy+1,iz)-(*u.get_grid(ilevel))(ix,iy-1,iz))*h;
else if( dir == 2 )
#pragma omp parallel for
for( int ix = 0; ix < (int)(*u.get_grid(ilevel)).size(0); ++ix )
for( int iy = 0; iy < (int)(*u.get_grid(ilevel)).size(1); ++iy )
for( int iz = 0; iz < (int)(*u.get_grid(ilevel)).size(2); ++iz )
(*pvar)(ix,iy,iz) = 0.5*((*u.get_grid(ilevel))(ix,iy,iz+1)-(*u.get_grid(ilevel))(ix,iy,iz-1))*h;
}
LOGUSER("Done computing a 2nd order finite difference gradient.");
}
void multigrid_poisson_plugin::implementation::gradient_add_O2( int dir, grid_hierarchy& u, grid_hierarchy& Du )
{
LOGUSER("Computing a 2nd order finite difference gradient...");
for( unsigned ilevel=u.levelmin(); ilevel<=u.levelmax(); ++ilevel )
{
double h = pow(2.0,ilevel);
meshvar_bnd *pvar = Du.get_grid(ilevel);
if( dir == 0 )
#pragma omp parallel for
for( int ix = 0; ix < (int)(*Du.get_grid(ilevel)).size(0); ++ix )
for( int iy = 0; iy < (int)(*Du.get_grid(ilevel)).size(1); ++iy )
for( int iz = 0; iz < (int)(*Du.get_grid(ilevel)).size(2); ++iz )
(*pvar)(ix,iy,iz) += 0.5*((*u.get_grid(ilevel))(ix+1,iy,iz)-(*u.get_grid(ilevel))(ix-1,iy,iz))*h;
else if( dir == 1 )
#pragma omp parallel for
for( int ix = 0; ix < (int)(*Du.get_grid(ilevel)).size(0); ++ix )
for( int iy = 0; iy < (int)(*Du.get_grid(ilevel)).size(1); ++iy )
for( int iz = 0; iz < (int)(*Du.get_grid(ilevel)).size(2); ++iz )
(*pvar)(ix,iy,iz) += 0.5*((*u.get_grid(ilevel))(ix,iy+1,iz)-(*u.get_grid(ilevel))(ix,iy-1,iz))*h;
else if( dir == 2 )
#pragma omp parallel for
for( int ix = 0; ix < (int)(*Du.get_grid(ilevel)).size(0); ++ix )
for( int iy = 0; iy < (int)(*Du.get_grid(ilevel)).size(1); ++iy )
for( int iz = 0; iz < (int)(*Du.get_grid(ilevel)).size(2); ++iz )
(*pvar)(ix,iy,iz) += 0.5*((*u.get_grid(ilevel))(ix,iy,iz+1)-(*u.get_grid(ilevel))(ix,iy,iz-1))*h;
}
LOGUSER("Done computing a 4th order finite difference gradient.");
}
void multigrid_poisson_plugin::implementation::gradient_O4( int dir, grid_hierarchy& u, grid_hierarchy& Du )
{
LOGUSER("Computing a 4th order finite difference gradient...");
for( unsigned ilevel=u.levelmin(); ilevel<=u.levelmax(); ++ilevel )
{
double h = pow(2.0,ilevel);
meshvar_bnd *pvar = Du.get_grid(ilevel);
h /= 12.0;
if( dir == 0 )
#pragma omp parallel for
for( int ix = 0; ix < (int)(*u.get_grid(ilevel)).size(0); ++ix )
for( int iy = 0; iy < (int)(*u.get_grid(ilevel)).size(1); ++iy )
for( int iz = 0; iz < (int)(*u.get_grid(ilevel)).size(2); ++iz )
(*pvar)(ix,iy,iz) = ((*u.get_grid(ilevel))(ix-2,iy,iz)
-8.0*(*u.get_grid(ilevel))(ix-1,iy,iz)
+8.0*(*u.get_grid(ilevel))(ix+1,iy,iz)
-(*u.get_grid(ilevel))(ix+2,iy,iz))*h;
else if( dir == 1 )
#pragma omp parallel for
for( int ix = 0; ix < (int)(*u.get_grid(ilevel)).size(0); ++ix )
for( int iy = 0; iy < (int)(*u.get_grid(ilevel)).size(1); ++iy )
for( int iz = 0; iz < (int)(*u.get_grid(ilevel)).size(2); ++iz )
(*pvar)(ix,iy,iz) = ((*u.get_grid(ilevel))(ix,iy-2,iz)
-8.0*(*u.get_grid(ilevel))(ix,iy-1,iz)
+8.0*(*u.get_grid(ilevel))(ix,iy+1,iz)
-(*u.get_grid(ilevel))(ix,iy+2,iz))*h;
else if( dir == 2 )
#pragma omp parallel for
for( int ix = 0; ix < (int)(*u.get_grid(ilevel)).size(0); ++ix )
for( int iy = 0; iy < (int)(*u.get_grid(ilevel)).size(1); ++iy )
for( int iz = 0; iz < (int)(*u.get_grid(ilevel)).size(2); ++iz )
(*pvar)(ix,iy,iz) = ((*u.get_grid(ilevel))(ix,iy,iz-2)
-8.0*(*u.get_grid(ilevel))(ix,iy,iz-1)
+8.0*(*u.get_grid(ilevel))(ix,iy,iz+1)
-(*u.get_grid(ilevel))(ix,iy,iz+2))*h;
}
LOGUSER("Done computing a 4th order finite difference gradient.");
}
void multigrid_poisson_plugin::implementation::gradient_add_O4( int dir, grid_hierarchy& u, grid_hierarchy& Du )
{
LOGUSER("Computing a 4th order finite difference gradient...");
for( unsigned ilevel=u.levelmin(); ilevel<=u.levelmax(); ++ilevel )
{
double h = pow(2.0,ilevel);
meshvar_bnd *pvar = Du.get_grid(ilevel);
h /= 12.0;
if( dir == 0 )
#pragma omp parallel for
for( int ix = 0; ix < (int)(*u.get_grid(ilevel)).size(0); ++ix )
for( int iy = 0; iy < (int)(*u.get_grid(ilevel)).size(1); ++iy )
for( int iz = 0; iz < (int)(*u.get_grid(ilevel)).size(2); ++iz )
(*pvar)(ix,iy,iz) += ((*u.get_grid(ilevel))(ix-2,iy,iz)
-8.0*(*u.get_grid(ilevel))(ix-1,iy,iz)
+8.0*(*u.get_grid(ilevel))(ix+1,iy,iz)
-(*u.get_grid(ilevel))(ix+2,iy,iz))*h;
else if( dir == 1 )
#pragma omp parallel for
for( int ix = 0; ix < (int)(*u.get_grid(ilevel)).size(0); ++ix )
for( int iy = 0; iy < (int)(*u.get_grid(ilevel)).size(1); ++iy )
for( int iz = 0; iz < (int)(*u.get_grid(ilevel)).size(2); ++iz )
(*pvar)(ix,iy,iz) += ((*u.get_grid(ilevel))(ix,iy-2,iz)
-8.0*(*u.get_grid(ilevel))(ix,iy-1,iz)
+8.0*(*u.get_grid(ilevel))(ix,iy+1,iz)
-(*u.get_grid(ilevel))(ix,iy+2,iz))*h;
else if( dir == 2 )
#pragma omp parallel for
for( int ix = 0; ix < (int)(*u.get_grid(ilevel)).size(0); ++ix )
for( int iy = 0; iy < (int)(*u.get_grid(ilevel)).size(1); ++iy )
for( int iz = 0; iz < (int)(*u.get_grid(ilevel)).size(2); ++iz )
(*pvar)(ix,iy,iz) += ((*u.get_grid(ilevel))(ix,iy,iz-2)
-8.0*(*u.get_grid(ilevel))(ix,iy,iz-1)
+8.0*(*u.get_grid(ilevel))(ix,iy,iz+1)
-(*u.get_grid(ilevel))(ix,iy,iz+2))*h;
}
LOGUSER("Done computing a 4th order finite difference gradient.");
}
void multigrid_poisson_plugin::implementation::gradient_O6( int dir, grid_hierarchy& u, grid_hierarchy& Du )
{
LOGUSER("Computing a 6th order finite difference gradient...");
for( unsigned ilevel=u.levelmin(); ilevel<=u.levelmax(); ++ilevel )
{
double h = pow(2.0,ilevel);
meshvar_bnd *pvar = Du.get_grid(ilevel);
h /= 60.;
if( dir == 0 )
#pragma omp parallel for
for( int ix = 0; ix < (int)(*u.get_grid(ilevel)).size(0); ++ix )
for( int iy = 0; iy < (int)(*u.get_grid(ilevel)).size(1); ++iy )
for( int iz = 0; iz < (int)(*u.get_grid(ilevel)).size(2); ++iz )
(*pvar)(ix,iy,iz) =
(-(*u.get_grid(ilevel))(ix-3,iy,iz)
+9.0*(*u.get_grid(ilevel))(ix-2,iy,iz)
-45.0*(*u.get_grid(ilevel))(ix-1,iy,iz)
+45.0*(*u.get_grid(ilevel))(ix+1,iy,iz)
-9.0*(*u.get_grid(ilevel))(ix+2,iy,iz)
+(*u.get_grid(ilevel))(ix+3,iy,iz))*h;
else if( dir == 1 )
#pragma omp parallel for
for( int ix = 0; ix < (int)(*u.get_grid(ilevel)).size(0); ++ix )
for( int iy = 0; iy < (int)(*u.get_grid(ilevel)).size(1); ++iy )
for( int iz = 0; iz < (int)(*u.get_grid(ilevel)).size(2); ++iz )
(*pvar)(ix,iy,iz) =
(-(*u.get_grid(ilevel))(ix,iy-3,iz)
+9.0*(*u.get_grid(ilevel))(ix,iy-2,iz)
-45.0*(*u.get_grid(ilevel))(ix,iy-1,iz)
+45.0*(*u.get_grid(ilevel))(ix,iy+1,iz)
-9.0*(*u.get_grid(ilevel))(ix,iy+2,iz)
+(*u.get_grid(ilevel))(ix,iy+3,iz))*h;
else if( dir == 2 )
#pragma omp parallel for
for( int ix = 0; ix < (int)(*u.get_grid(ilevel)).size(0); ++ix )
for( int iy = 0; iy < (int)(*u.get_grid(ilevel)).size(1); ++iy )
for( int iz = 0; iz < (int)(*u.get_grid(ilevel)).size(2); ++iz )
(*pvar)(ix,iy,iz) =
(-(*u.get_grid(ilevel))(ix,iy,iz-3)
+9.0*(*u.get_grid(ilevel))(ix,iy,iz-2)
-45.0*(*u.get_grid(ilevel))(ix,iy,iz-1)
+45.0*(*u.get_grid(ilevel))(ix,iy,iz+1)
-9.0*(*u.get_grid(ilevel))(ix,iy,iz+2)
+(*u.get_grid(ilevel))(ix,iy,iz+3))*h;
}
LOGUSER("Done computing a 6th order finite difference gradient.");
}
void multigrid_poisson_plugin::implementation::gradient_add_O6( int dir, grid_hierarchy& u, grid_hierarchy& Du )
{
LOGUSER("Computing a 6th order finite difference gradient...");
for( unsigned ilevel=u.levelmin(); ilevel<=u.levelmax(); ++ilevel )
{
double h = pow(2.0,ilevel);
meshvar_bnd *pvar = Du.get_grid(ilevel);
h /= 60.;
if( dir == 0 )
#pragma omp parallel for
for( int ix = 0; ix < (int)(*u.get_grid(ilevel)).size(0); ++ix )
for( int iy = 0; iy < (int)(*u.get_grid(ilevel)).size(1); ++iy )
for( int iz = 0; iz < (int)(*u.get_grid(ilevel)).size(2); ++iz )
(*pvar)(ix,iy,iz) +=
(-(*u.get_grid(ilevel))(ix-3,iy,iz)
+9.0*(*u.get_grid(ilevel))(ix-2,iy,iz)
-45.0*(*u.get_grid(ilevel))(ix-1,iy,iz)
+45.0*(*u.get_grid(ilevel))(ix+1,iy,iz)
-9.0*(*u.get_grid(ilevel))(ix+2,iy,iz)
+(*u.get_grid(ilevel))(ix+3,iy,iz))*h;
else if( dir == 1 )
#pragma omp parallel for
for( int ix = 0; ix < (int)(*u.get_grid(ilevel)).size(0); ++ix )
for( int iy = 0; iy < (int)(*u.get_grid(ilevel)).size(1); ++iy )
for( int iz = 0; iz < (int)(*u.get_grid(ilevel)).size(2); ++iz )
(*pvar)(ix,iy,iz) +=
(-(*u.get_grid(ilevel))(ix,iy-3,iz)
+9.0*(*u.get_grid(ilevel))(ix,iy-2,iz)
-45.0*(*u.get_grid(ilevel))(ix,iy-1,iz)
+45.0*(*u.get_grid(ilevel))(ix,iy+1,iz)
-9.0*(*u.get_grid(ilevel))(ix,iy+2,iz)
+(*u.get_grid(ilevel))(ix,iy+3,iz))*h;
else if( dir == 2 )
#pragma omp parallel for
for( int ix = 0; ix < (int)(*u.get_grid(ilevel)).size(0); ++ix )
for( int iy = 0; iy < (int)(*u.get_grid(ilevel)).size(1); ++iy )
for( int iz = 0; iz < (int)(*u.get_grid(ilevel)).size(2); ++iz )
(*pvar)(ix,iy,iz) +=
(-(*u.get_grid(ilevel))(ix,iy,iz-3)
+9.0*(*u.get_grid(ilevel))(ix,iy,iz-2)
-45.0*(*u.get_grid(ilevel))(ix,iy,iz-1)
+45.0*(*u.get_grid(ilevel))(ix,iy,iz+1)
-9.0*(*u.get_grid(ilevel))(ix,iy,iz+2)
+(*u.get_grid(ilevel))(ix,iy,iz+3))*h;
}
LOGUSER("Done computing a 6th order finite difference gradient.");
}
/**************************************************************************************/
/**************************************************************************************/
#pragma mark -
#include "general.hh"
double fft_poisson_plugin::solve( grid_hierarchy& f, grid_hierarchy& u )
{
LOGUSER("Entering k-space Poisson solver...");
unsigned verbosity = cf_.getValueSafe<unsigned>("setup","verbosity",2);
if( f.levelmin() != f.levelmax() )
{
LOGERR("Attempt to run k-space Poisson solver on non unigrid mesh.");
throw std::runtime_error("fft_poisson_plugin::solve : k-space method can only be used in unigrid mode (levelmin=levelmax)");
}
if( verbosity > 0 )
{
std::cout << "-------------------------------------------------------------\n";
std::cout << " - Invoking unigrid FFT Poisson solver..." << std::endl;
}
int nx,ny,nz,nzp;
nx = f.get_grid(f.levelmax())->size(0);
ny = f.get_grid(f.levelmax())->size(1);
nz = f.get_grid(f.levelmax())->size(2);
nzp = 2*(nz/2+1);
//... copy data ..................................................
fftw_real *data = new fftw_real[(size_t)nx*(size_t)ny*(size_t)nzp];
fftw_complex *cdata = reinterpret_cast<fftw_complex*> (data);
#pragma omp parallel for
for( int i=0; i<nx; ++i )
for( int j=0; j<ny; ++j )
for( int k=0; k<nz; ++k )
{
size_t idx = (size_t)(i*ny+j)*(size_t)nzp+(size_t)k;
data[idx] = (*f.get_grid(f.levelmax()))(i,j,k);
}
//... perform FFT and Poisson solve................................
LOGUSER("Performing forward transform.");
#ifdef FFTW3
#ifdef SINGLE_PRECISION
fftwf_plan
plan = fftwf_plan_dft_r2c_3d( nx, ny, nz, data, cdata, FFTW_ESTIMATE ),
iplan = fftwf_plan_dft_c2r_3d( nx, ny, nz, cdata, data, FFTW_ESTIMATE );
fftwf_execute(plan);
#else
fftw_plan
plan = fftw_plan_dft_r2c_3d( nx, ny, nz, data, cdata, FFTW_ESTIMATE ),
iplan = fftw_plan_dft_c2r_3d( nx, ny, nz, cdata, data, FFTW_ESTIMATE );
fftw_execute(plan);
#endif
#else
rfftwnd_plan
plan = rfftw3d_create_plan( nx,ny,nz,
FFTW_REAL_TO_COMPLEX, FFTW_ESTIMATE|FFTW_IN_PLACE),
iplan = rfftw3d_create_plan( nx,ny,nz,
FFTW_COMPLEX_TO_REAL, FFTW_ESTIMATE|FFTW_IN_PLACE);
#ifndef SINGLETHREAD_FFTW
rfftwnd_threads_one_real_to_complex( omp_get_max_threads(), plan, data, NULL );
#else
rfftwnd_one_real_to_complex( plan, data, NULL );
#endif
#endif
double kfac = 2.0*M_PI;
double fac = -1.0/(double)((size_t)nx*(size_t)ny*(size_t)nz);
#pragma omp parallel for
for( int i=0; i<nx; ++i )
for( int j=0; j<ny; ++j )
for( int k=0; k<nz/2+1; ++k )
{
int ii = i; if(ii>nx/2) ii-=nx;
int jj = j; if(jj>ny/2) jj-=ny;
double ki = (double)ii;
double kj = (double)jj;
double kk = (double)k;
double kk2 = kfac*kfac*(ki*ki+kj*kj+kk*kk);
size_t idx = (size_t)(i*ny+j)*(size_t)(nzp/2)+(size_t)k;
RE(cdata[idx]) *= -1.0/kk2*fac;
IM(cdata[idx]) *= -1.0/kk2*fac;
}
RE(cdata[0]) = 0.0;
IM(cdata[0]) = 0.0;
LOGUSER("Performing backward transform.");
#ifdef FFTW3
#ifdef SINGLE_PRECISION
fftwf_execute(iplan);
fftwf_destroy_plan(plan);
fftwf_destroy_plan(iplan);
#else
fftw_execute(iplan);
fftw_destroy_plan(plan);
fftw_destroy_plan(iplan);
#endif
#else
#ifndef SINGLETHREAD_FFTW
rfftwnd_threads_one_complex_to_real( omp_get_max_threads(), iplan, cdata, NULL );
#else
rfftwnd_one_complex_to_real( iplan, cdata, NULL );
#endif
rfftwnd_destroy_plan(plan);
rfftwnd_destroy_plan(iplan);
#endif
//... copy data ..........................................
#pragma omp parallel for
for( int i=0; i<nx; ++i )
for( int j=0; j<ny; ++j )
for( int k=0; k<nz; ++k )
{
size_t idx = (size_t)(i*ny+j)*(size_t)nzp+(size_t)k;
(*u.get_grid(u.levelmax()))(i,j,k) = data[idx];
}
delete[] data;
//... set boundary values ................................
int nb = u.get_grid(u.levelmax())->m_nbnd;
for( int iy=-nb; iy<ny+nb; ++iy )
for( int iz=-nb; iz<nz+nb; ++iz )
{
int iiy( (iy+ny)%ny ), iiz( (iz+nz)%nz );
for( int i=-nb; i<0; ++i )
{
(*u.get_grid(u.levelmax()))(i,iy,iz) = (*u.get_grid(u.levelmax()))(nx+i,iiy,iiz);
(*u.get_grid(u.levelmax()))(nx-1-i,iy,iz) = (*u.get_grid(u.levelmax()))(-1-i,iiy,iiz);
}
}
for( int ix=-nb; ix<nx+nb; ++ix )
for( int iz=-nb; iz<nz+nb; ++iz )
{
int iix( (ix+nx)%nx ), iiz( (iz+nz)%nz );
for( int i=-nb; i<0; ++i )
{
(*u.get_grid(u.levelmax()))(ix,i,iz) = (*u.get_grid(u.levelmax()))(iix,ny+i,iiz);
(*u.get_grid(u.levelmax()))(ix,ny-1-i,iz) = (*u.get_grid(u.levelmax()))(iix,-1-i,iiz);
}
}
for( int ix=-nb; ix<nx+nb; ++ix )
for( int iy=-nb; iy<ny+nb; ++iy )
{
int iix( (ix+nx)%nx ), iiy( (iy+ny)%ny );
for( int i=-nb; i<0; ++i )
{
(*u.get_grid(u.levelmax()))(ix,iy,i) = (*u.get_grid(u.levelmax()))(iix,iiy,nz+i);
(*u.get_grid(u.levelmax()))(ix,iy,nz-1-i) = (*u.get_grid(u.levelmax()))(iix,iiy,-1-i);
}
}
LOGUSER("Done with k-space Poisson solver.");
return 0.0;
}
double fft_poisson_plugin::gradient( int dir, grid_hierarchy& u, grid_hierarchy& Du )
{
LOGUSER("Computing a gradient in k-space...\n");
if( u.levelmin() != u.levelmax() )
throw std::runtime_error("fft_poisson_plugin::gradient : k-space method can only be used in unigrid mode (levelmin=levelmax)");
Du = u;
int nx,ny,nz,nzp;
nx = u.get_grid(u.levelmax())->size(0);
ny = u.get_grid(u.levelmax())->size(1);
nz = u.get_grid(u.levelmax())->size(2);
nzp = 2*(nz/2+1);
//... copy data ..................................................
fftw_real *data = new fftw_real[(size_t)nx*(size_t)ny*(size_t)nzp];
fftw_complex *cdata = reinterpret_cast<fftw_complex*> (data);
#pragma omp parallel for
for( int i=0; i<nx; ++i )
for( int j=0; j<ny; ++j )
for( int k=0; k<nz; ++k )
{
size_t idx = (size_t)(i*ny+j)*(size_t)nzp+(size_t)k;
data[idx] = (*u.get_grid(u.levelmax()))(i,j,k);
}
//... perform FFT and Poisson solve................................
#ifdef FFTW3
#ifdef SINGLE_PRECISION
fftwf_plan
plan = fftwf_plan_dft_r2c_3d(nx, ny, nz, data, cdata, FFTW_ESTIMATE),
iplan = fftwf_plan_dft_c2r_3d(nx, ny, nz, cdata, data, FFTW_ESTIMATE);
fftwf_execute(plan);
#else
fftw_plan
plan = fftw_plan_dft_r2c_3d(nx, ny, nz, data, cdata, FFTW_ESTIMATE),
iplan = fftw_plan_dft_c2r_3d(nx, ny, nz, cdata, data, FFTW_ESTIMATE);
fftw_execute(plan);
#endif
#else
rfftwnd_plan
plan = rfftw3d_create_plan( nx,ny,nz,
FFTW_REAL_TO_COMPLEX, FFTW_ESTIMATE|FFTW_IN_PLACE),
iplan = rfftw3d_create_plan( nx,ny,nz,
FFTW_COMPLEX_TO_REAL, FFTW_ESTIMATE|FFTW_IN_PLACE);
#ifndef SINGLETHREAD_FFTW
rfftwnd_threads_one_real_to_complex( omp_get_max_threads(), plan, data, NULL );
#else
rfftwnd_one_real_to_complex( plan, data, NULL );
#endif
#endif
double fac = -1.0/(double)((size_t)nx*(size_t)ny*(size_t)nz);
double kfac = 2.0*M_PI;
bool do_glass = cf_.getValueSafe<bool>("output","glass",false);
bool deconvolve_cic = do_glass | cf_.getValueSafe<bool>("output","glass_cicdeconvolve",false);
if( deconvolve_cic )
LOGINFO("CIC deconvolution is enabled for kernel!");
#pragma omp parallel for
for( int i=0; i<nx; ++i )
for( int j=0; j<ny; ++j )
for( int k=0; k<nz/2+1; ++k )
{
size_t idx = (size_t)(i*ny+j)*(size_t)(nzp/2)+(size_t)k;
int ii = i; if(ii>nx/2) ii-=nx;
int jj = j; if(jj>ny/2) jj-=ny;
const double ki = (double)ii;
const double kj = (double)jj;
const double kk = (double)k;
const double kkdir[3] = {kfac*ki,kfac*kj,kfac*kk};
const double kdir = kkdir[dir];
double re = RE(cdata[idx]);
double im = IM(cdata[idx]);
RE(cdata[idx]) = fac*im*kdir;
IM(cdata[idx]) = -fac*re*kdir;
#ifdef FFTW3
if( deconvolve_cic )
{
double dfx, dfy, dfz;
dfx = M_PI*ki/(double)nx; dfx = (i!=0)? sin(dfx)/dfx : 1.0;
dfy = M_PI*kj/(double)ny; dfy = (j!=0)? sin(dfy)/dfy : 1.0;
dfz = M_PI*kk/(double)nz; dfz = (k!=0)? sin(dfz)/dfz : 1.0;
dfx = 1.0/(dfx*dfy*dfz); dfx = dfx*dfx;
cdata[idx][0] *= dfx;
cdata[idx][1] *= dfx;
}
#else
if( deconvolve_cic )
{
double dfx, dfy, dfz;
dfx = M_PI*ki/(double)nx; dfx = (i!=0)? sin(dfx)/dfx : 1.0;
dfy = M_PI*kj/(double)ny; dfy = (j!=0)? sin(dfy)/dfy : 1.0;
dfz = M_PI*kk/(double)nz; dfz = (k!=0)? sin(dfz)/dfz : 1.0;
dfx = 1.0/(dfx*dfy*dfz); dfx = dfx*dfx;
cdata[idx].re *= dfx;
cdata[idx].im *= dfx;
}
#endif
/*double ktot = sqrt(ii*ii+jj*jj+k*k);
if( ktot >= nx/2 )//dir == 0 && i==nx/2 || dir == 1 && j==ny/2 || dir == 2 && k==nz/2 )
{
cdata[idx].re = 0.0;
cdata[idx].im = 0.0;
}*/
}
RE(cdata[0]) = 0.0;
IM(cdata[0]) = 0.0;
#ifdef FFTW3
#ifdef SINGLE_PRECISION
fftwf_execute(iplan);
fftwf_destroy_plan(plan);
fftwf_destroy_plan(iplan);
#else
fftw_execute(iplan);
fftw_destroy_plan(plan);
fftw_destroy_plan(iplan);
#endif
#else
#ifndef SINGLETHREAD_FFTW
rfftwnd_threads_one_complex_to_real( omp_get_max_threads(), iplan, cdata, NULL );
#else
rfftwnd_one_complex_to_real( iplan, cdata, NULL );
#endif
rfftwnd_destroy_plan(plan);
rfftwnd_destroy_plan(iplan);
#endif
//... copy data ..........................................
double dmax = 0.0;
for( int i=0; i<nx; ++i )
for( int j=0; j<ny; ++j )
for( int k=0; k<nz; ++k )
{
size_t idx = ((size_t)i*ny+(size_t)j)*nzp+(size_t)k;
(*Du.get_grid(u.levelmax()))(i,j,k) = data[idx];
if(fabs(data[idx])>dmax)
dmax = fabs(data[idx]);
}
delete[] data;
LOGUSER("Done with k-space gradient.\n");
return 0.0;
}
/**************************************************************************************/
/**************************************************************************************/
#pragma mark -
template<int order>
double poisson_hybrid_kernel( int idir, int i, int j, int k, int n )
{
return 1.0;
}
template<>
inline double poisson_hybrid_kernel<2>(int idir, int i, int j, int k, int n )
{
if(i==0&&j==0&&k==0)
return 0.0;
double
ki(M_PI*(double)i/(double)n),
kj(M_PI*(double)j/(double)n),
kk(M_PI*(double)k/(double)n),
kr(sqrt(ki*ki+kj*kj+kk*kk));
double grad = 1.0, laplace = 1.0;
if( idir==0 )
grad = sin(ki);
else if( idir==1 )
grad = sin(kj);
else
grad = sin(kk);
laplace = 2.0*((-cos(ki)+1.0)+(-cos(kj)+1.0)+(-cos(kk)+1.0));
double kgrad = 1.0;
if( idir==0 )
kgrad = ki;
else if( idir ==1)
kgrad = kj;
else if( idir ==2)
kgrad = kk;
return kgrad/kr/kr-grad/laplace;
}
template<>
inline double poisson_hybrid_kernel<4>(int idir, int i, int j, int k, int n )
{
if(i==0&&j==0&&k==0)
return 0.0;
double
ki(M_PI*(double)i/(double)n),
kj(M_PI*(double)j/(double)n),
kk(M_PI*(double)k/(double)n),
kr(sqrt(ki*ki+kj*kj+kk*kk));
double grad = 1.0, laplace = 1.0;
if( idir==0 )
grad = 0.166666666667*(-sin(2.*ki)+8.*sin(ki));
else if( idir==1 )
grad = 0.166666666667*(-sin(2.*kj)+8.*sin(kj));
else if( idir==2 )
grad = 0.166666666667*(-sin(2.*kk)+8.*sin(kk));
laplace = 0.1666666667*((cos(2*ki)-16.*cos(ki)+15.)
+(cos(2*kj)-16.*cos(kj)+15.)
+(cos(2*kk)-16.*cos(kk)+15.));
double kgrad = 1.0;
if( idir==0 )
kgrad = ki;
else if( idir ==1)
kgrad = kj;
else if( idir ==2)
kgrad = kk;
return kgrad/kr/kr-grad/laplace;
}
template<>
inline double poisson_hybrid_kernel<6>(int idir, int i, int j, int k, int n )
{
double
ki(M_PI*(double)i/(double)n),
kj(M_PI*(double)j/(double)n),
kk(M_PI*(double)k/(double)n),
kr(sqrt(ki*ki+kj*kj+kk*kk));
if(i==0&&j==0&&k==0)
return 0.0;
double grad = 1.0, laplace = 1.0;
if( idir==0 )
grad = 0.0333333333333*(sin(3.*ki)-9.*sin(2.*ki)+45.*sin(ki));
else if( idir==1 )
grad = 0.0333333333333*(sin(3.*kj)-9.*sin(2.*kj)+45.*sin(kj));
else if( idir==2 )
grad = 0.0333333333333*(sin(3.*kk)-9.*sin(2.*kk)+45.*sin(kk));
laplace = 0.01111111111111*(
(-2.*cos(3.0*ki)+27.*cos(2.*ki)-270.*cos(ki)+245.)
+(-2.*cos(3.0*kj)+27.*cos(2.*kj)-270.*cos(kj)+245.)
+(-2.*cos(3.0*kk)+27.*cos(2.*kk)-270.*cos(kk)+245.));
double kgrad = 1.0;
if( idir==0 )
kgrad = ki;
else if( idir ==1)
kgrad = kj;
else if( idir ==2)
kgrad = kk;
// if( i*i+j*j+k*k >= n*n )
// kgrad = 0.0;
return kgrad/kr/kr-grad/laplace;
}
template<int order>
void do_poisson_hybrid( fftw_real* data, int idir, int nxp, int nyp, int nzp, bool periodic, bool deconvolve_cic )
{
double fftnorm = 1.0/((double)nxp*(double)nyp*(double)nzp);
fftw_complex *cdata = reinterpret_cast<fftw_complex*>(data);
if( deconvolve_cic )
LOGINFO("CIC deconvolution step is enabled.");
#ifdef FFTW3
#ifdef SINGLE_PRECISION
fftwf_plan iplan, plan;
plan = fftwf_plan_dft_r2c_3d(nxp, nyp, nzp, data, cdata, FFTW_ESTIMATE);
iplan = fftwf_plan_dft_c2r_3d(nxp, nyp, nzp, cdata, data, FFTW_ESTIMATE);
fftwf_execute(plan);
#else
fftw_plan iplan, plan;
plan = fftw_plan_dft_r2c_3d(nxp, nyp, nzp, data, cdata, FFTW_ESTIMATE);
iplan = fftw_plan_dft_c2r_3d(nxp, nyp, nzp, cdata, data, FFTW_ESTIMATE);
fftw_execute(plan);
#endif
#else
rfftwnd_plan iplan, plan;
plan = rfftw3d_create_plan( nxp, nyp, nzp,
FFTW_REAL_TO_COMPLEX, FFTW_ESTIMATE|FFTW_IN_PLACE);
iplan = rfftw3d_create_plan( nxp, nyp, nzp,
FFTW_COMPLEX_TO_REAL, FFTW_ESTIMATE|FFTW_IN_PLACE);
#ifndef SINGLETHREAD_FFTW
rfftwnd_threads_one_real_to_complex( omp_get_max_threads(), plan, data, NULL );
#else
rfftwnd_one_real_to_complex( plan, data, NULL );
#endif
#endif