-
Notifications
You must be signed in to change notification settings - Fork 2
/
forcetree.c
3945 lines (3337 loc) · 99 KB
/
forcetree.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 <stdlib.h>
#include <string.h>
#include <math.h>
#include <time.h>
#include <mpi.h>
#include "allvars.h"
#include "proto.h"
#include "ngravs.h"
/*! \file forcetree.c
* \brief gravitational tree and code for Ewald correction
*
* This file contains the computation of the gravitational force by means
* of a tree. The type of tree implemented is a geometrical oct-tree,
* starting from a cube encompassing all particles. This cube is
* automatically found in the domain decomposition, which also splits up
* the global "top-level" tree along node boundaries, moving the particles
* of different parts of the tree to separate processors. Tree nodes can
* be dynamically updated in drift/kick operations to avoid having to
* reconstruct the tree every timestep.
*/
/*! auxialiary variable used to set-up non-recursive walk */
static int last;
#ifdef PMGRID
// KC 2/8/16
// We make these double because we are now doing subtractions with them, instead of multiplications
// so we need the additional precision.
static double shortrange_fourier_pot[N_GRAVS][N_GRAVS][NTAB], shortrange_fourier_force[N_GRAVS][N_GRAVS][NTAB];
extern struct ngravsInterpolant *ngravsPeriodicTable;
#endif
/*! toggles after first tree-memory allocation, has only influence on log-files */
static int first_flag = 0;
#ifdef PERIODIC
/*! Macro that maps a distance to the nearest periodic neighbour */
#define NEAREST(x) (((x)>boxhalf)?((x)-boxsize):(((x)<-boxhalf)?((x)+boxsize):(x)))
/*! Size of 3D lock-up table for Lattice sum correction force */
/*! 3D lock-up table for Lattice sum correction to force and potential. Only one
* octant is stored, the rest constructed by using the symmetry
*/
static double fcorrx[N_GRAVS][N_GRAVS][NGRAVS_EN + 1][NGRAVS_EN + 1][NGRAVS_EN + 1];
static double fcorry[N_GRAVS][N_GRAVS][NGRAVS_EN + 1][NGRAVS_EN + 1][NGRAVS_EN + 1];
static double fcorrz[N_GRAVS][N_GRAVS][NGRAVS_EN + 1][NGRAVS_EN + 1][NGRAVS_EN + 1];
static double potcorr[N_GRAVS][N_GRAVS][NGRAVS_EN + 1][NGRAVS_EN + 1][NGRAVS_EN + 1];
static double fac_intp;
#endif
/*! This function is a driver routine for constructing the gravitational
* oct-tree, which is done by calling a small number of other functions.
*/
int force_treebuild(int npart)
{
Numnodestree = force_treebuild_single(npart);
force_update_pseudoparticles();
force_flag_localnodes();
TimeOfLastTreeConstruction = All.Time;
return Numnodestree;
}
/*! Constructs the gravitational oct-tree.
*
* The index convention for accessing tree nodes is the following: the
* indices 0...NumPart-1 reference single particles, the indices
* All.MaxPart.... All.MaxPart+nodes-1 reference tree nodes. `Nodes_base'
* points to the first tree node, while `nodes' is shifted such that
* nodes[All.MaxPart] gives the first tree node. Finally, node indices
* with values 'All.MaxPart + MaxNodes' and larger indicate "pseudo
* particles", i.e. multipole moments of top-level nodes that lie on
* different CPUs. If such a node needs to be opened, the corresponding
* particle must be exported to that CPU. The 'Extnodes' structure
* parallels that of 'Nodes'. Its information is only needed for the SPH
* part of the computation. (The data is split onto these two structures
* as a tuning measure. If it is merged into 'Nodes' a somewhat bigger
* size of the nodes also for gravity would result, which would reduce
* cache utilization slightly.
*/
int force_treebuild_single(int npart)
{
int i, j, subnode = 0, parent, numnodes;
int nfree, th, nn, no;
struct NODE *nfreep;
double lenhalf, epsilon;
peanokey key;
/* create an empty root node */
nfree = All.MaxPart; /* index of first free node */
nfreep = &Nodes[nfree]; /* select first node */
nfreep->len = DomainLen;
for(j = 0; j < 3; j++)
nfreep->center[j] = DomainCenter[j];
for(j = 0; j < 8; j++)
nfreep->u.suns[j] = -1;
numnodes = 1;
nfreep++;
nfree++;
/* create a set of empty nodes corresponding to the top-level domain
* grid. We need to generate these nodes first to make sure that we have a
* complete top-level tree which allows the easy insertion of the
* pseudo-particles at the right place
*/
force_create_empty_nodes(All.MaxPart, 0, 1, 0, 0, 0, &numnodes, &nfree);
/* if a high-resolution region in a global tree is used, we need to generate
* an additional set empty nodes to make sure that we have a complete
* top-level tree for the high-resolution inset
*/
nfreep = &Nodes[nfree];
parent = -1; /* note: will not be used below before it is changed */
/* now we insert all particles */
for(i = 0; i < npart; i++)
{
/* the softening is only used to check whether particles are so close
* that the tree needs not to be refined further
*/
epsilon = All.ForceSoftening[P[i].Type];
key = peano_hilbert_key((P[i].Pos[0] - DomainCorner[0]) * DomainFac,
(P[i].Pos[1] - DomainCorner[1]) * DomainFac,
(P[i].Pos[2] - DomainCorner[2]) * DomainFac, BITS_PER_DIMENSION);
no = 0;
while(TopNodes[no].Daughter >= 0)
no = TopNodes[no].Daughter + (key - TopNodes[no].StartKey) / (TopNodes[no].Size / 8);
no = TopNodes[no].Leaf;
th = DomainNodeIndex[no];
while(1)
{
if(th >= All.MaxPart) /* we are dealing with an internal node */
{
subnode = 0;
if(P[i].Pos[0] > Nodes[th].center[0])
subnode += 1;
if(P[i].Pos[1] > Nodes[th].center[1])
subnode += 2;
if(P[i].Pos[2] > Nodes[th].center[2])
subnode += 4;
nn = Nodes[th].u.suns[subnode];
if(nn >= 0) /* ok, something is in the daughter slot already, need to continue */
{
parent = th;
th = nn;
}
else
{
/* here we have found an empty slot where we can attach
* the new particle as a leaf.
*/
Nodes[th].u.suns[subnode] = i;
break; /* done for this particle */
}
}
else
{
/* We try to insert into a leaf with a single particle. Need
* to generate a new internal node at this point.
*/
Nodes[parent].u.suns[subnode] = nfree;
nfreep->len = 0.5 * Nodes[parent].len;
lenhalf = 0.25 * Nodes[parent].len;
if(subnode & 1)
nfreep->center[0] = Nodes[parent].center[0] + lenhalf;
else
nfreep->center[0] = Nodes[parent].center[0] - lenhalf;
if(subnode & 2)
nfreep->center[1] = Nodes[parent].center[1] + lenhalf;
else
nfreep->center[1] = Nodes[parent].center[1] - lenhalf;
if(subnode & 4)
nfreep->center[2] = Nodes[parent].center[2] + lenhalf;
else
nfreep->center[2] = Nodes[parent].center[2] - lenhalf;
nfreep->u.suns[0] = -1;
nfreep->u.suns[1] = -1;
nfreep->u.suns[2] = -1;
nfreep->u.suns[3] = -1;
nfreep->u.suns[4] = -1;
nfreep->u.suns[5] = -1;
nfreep->u.suns[6] = -1;
nfreep->u.suns[7] = -1;
subnode = 0;
if(P[th].Pos[0] > nfreep->center[0])
subnode += 1;
if(P[th].Pos[1] > nfreep->center[1])
subnode += 2;
if(P[th].Pos[2] > nfreep->center[2])
subnode += 4;
#ifndef NOTREERND
if(nfreep->len < 1.0e-3 * epsilon)
{
/* seems like we're dealing with particles at identical (or extremely close)
* locations. Randomize subnode index to allow tree construction. Note: Multipole moments
* of tree are still correct, but this will only happen well below gravitational softening
* length-scale anyway.
*/
subnode = (int) (8.0 * get_random_number((0xffff & P[i].ID) + P[i].GravCost));
P[i].GravCost += 1;
if(subnode >= 8)
subnode = 7;
}
#endif
nfreep->u.suns[subnode] = th;
th = nfree; /* resume trying to insert the new particle at
* the newly created internal node
*/
numnodes++;
nfree++;
nfreep++;
if((numnodes) >= MaxNodes)
{
printf("task %d: maximum number %d of tree-nodes reached.\n", ThisTask, MaxNodes);
printf("for particle %d\n", i);
dump_particles();
endrun(1);
}
}
}
}
/* insert the pseudo particles that represent the mass distribution of other domains */
force_insert_pseudo_particles();
/* now compute the multipole moments recursively */
last = -1;
force_update_node_recursive(All.MaxPart, -1, -1);
if(last >= All.MaxPart)
{
if(last >= All.MaxPart + MaxNodes) /* a pseudo-particle */
Nextnode[last - MaxNodes] = -1;
else
Nodes[last].u.d.nextnode = -1;
}
else
Nextnode[last] = -1;
return numnodes;
}
/*! This function recursively creates a set of empty tree nodes which
* corresponds to the top-level tree for the domain grid. This is done to
* ensure that this top-level tree is always "complete" so that we can
* easily associate the pseudo-particles of other CPUs with tree-nodes at
* a given level in the tree, even when the particle population is so
* sparse that some of these nodes are actually empty.
*/
void force_create_empty_nodes(int no, int topnode, int bits, int x, int y, int z, int *nodecount,
int *nextfree)
{
int i, j, k, n, sub, count;
if(TopNodes[topnode].Daughter >= 0)
{
for(i = 0; i < 2; i++)
for(j = 0; j < 2; j++)
for(k = 0; k < 2; k++)
{
sub = 7 & peano_hilbert_key((x << 1) + i, (y << 1) + j, (z << 1) + k, bits);
count = i + 2 * j + 4 * k;
Nodes[no].u.suns[count] = *nextfree;
Nodes[*nextfree].len = 0.5 * Nodes[no].len;
Nodes[*nextfree].center[0] = Nodes[no].center[0] + (2 * i - 1) * 0.25 * Nodes[no].len;
Nodes[*nextfree].center[1] = Nodes[no].center[1] + (2 * j - 1) * 0.25 * Nodes[no].len;
Nodes[*nextfree].center[2] = Nodes[no].center[2] + (2 * k - 1) * 0.25 * Nodes[no].len;
for(n = 0; n < 8; n++)
Nodes[*nextfree].u.suns[n] = -1;
if(TopNodes[TopNodes[topnode].Daughter + sub].Daughter == -1)
DomainNodeIndex[TopNodes[TopNodes[topnode].Daughter + sub].Leaf] = *nextfree;
*nextfree = *nextfree + 1;
*nodecount = *nodecount + 1;
if((*nodecount) >= MaxNodes)
{
printf("task %d: maximum number %d of tree-nodes reached.\n", ThisTask, MaxNodes);
printf("in create empty nodes\n");
dump_particles();
endrun(11);
}
force_create_empty_nodes(*nextfree - 1, TopNodes[topnode].Daughter + sub,
bits + 1, 2 * x + i, 2 * y + j, 2 * z + k, nodecount, nextfree);
}
}
}
/*! this function inserts pseudo-particles which will represent the mass
* distribution of the other CPUs. Initially, the mass of the
* pseudo-particles is set to zero, and their coordinate is set to the
* center of the domain-cell they correspond to. These quantities will be
* updated later on.
*/
void force_insert_pseudo_particles(void)
{
int i, index, subnode, nn, th;
int j;
for(i = 0; i < NTopleaves; i++)
{
index = DomainNodeIndex[i];
for(j = 0; j < N_GRAVS; ++j) {
DomainMoment[i].s[0][j] = Nodes[index].center[0];
DomainMoment[i].s[1][j] = Nodes[index].center[1];
DomainMoment[i].s[2][j] = Nodes[index].center[2];
DomainMoment[i].mass[j] = 0;
#ifdef N_GRAVS_ACCUMULATOR
DomainMoment[i].Nparticles[j] = 0;
#endif
}
}
for(i = 0; i < NTopleaves; i++)
{
if(i < DomainMyStart || i > DomainMyLast)
{
th = All.MaxPart; /* select index of first node in tree */
while(1)
{
if(th >= All.MaxPart) /* we are dealing with an internal node */
{
if(th >= All.MaxPart + MaxNodes)
endrun(888); /* this can't be */
subnode = 0;
// KC 8/11/14 This looks like a bisection check, so the generalization should be that
// if any of the mass center coordinates exceeds the center, than we subnode up
for(j = 0; j < N_GRAVS; ++j) {
if(DomainMoment[i].s[0][j] > Nodes[th].center[0]) {
subnode += 1;
break;
}
}
for(j = 0; j < N_GRAVS; ++j) {
if(DomainMoment[i].s[1][j] > Nodes[th].center[1]) {
subnode += 2;
break;
}
}
for(j = 0; j < N_GRAVS; ++j) {
if(DomainMoment[i].s[2][j] > Nodes[th].center[2]) {
subnode += 4;
break;
}
}
//fprintf(stderr, "subnode structure: [%d, %d, %d]\n", i, th, subnode);
nn = Nodes[th].u.suns[subnode];
if(nn >= 0) /* ok, something is in the daughter slot already, need to continue */
{
th = nn;
}
else
{
/* here we have found an empty slot where we can
* attach the pseudo particle as a leaf
*/
Nodes[th].u.suns[subnode] = All.MaxPart + MaxNodes + i;
break; /* done for this pseudo particle */
}
}
else
{
endrun(889); /* this can't be */
}
}
}
}
}
// KC 8/10/14 Now we update this function to correctly determine the (monopole) moments
// of all the different types of interactions we can have
/*! this routine determines the multipole moments for a given internal node
* and all its subnodes using a recursive computation. The result is
* stored in the Nodes[] structure in the sequence of this tree-walk.
*
* Note that the bitflags-variable for each node is used to store in the
* lowest bits some special information: Bit 0 flags whether the node
* belongs to the top-level tree corresponding to the domain
* decomposition, while Bit 1 signals whether the top-level node is
* dependent on local mass.
*
* If UNEQUALSOFTENINGS is set, bits 2-4 give the particle type with
* the maximum softening among the particles in the node, and bit 5
* flags whether the node contains any particles with lower softening
* than that.
*/
void force_update_node_recursive(int no, int sib, int father)
{
int j, jj, p, pp, nextsib, suns[8];
FLOAT hmax;
#ifdef UNEQUALSOFTENINGS
#ifndef ADAPTIVE_GRAVSOFT_FORGAS
int maxsofttype, diffsoftflag;
#else
FLOAT maxsoft;
#endif
#endif
struct particle_data *pa;
int i;
// KC 8/11/14 Updated to reflect the different types
// ??? Why is this a double explicitly and not a FLOAT?
// Perhaps it is because we are adding many things, and so we don't
// want to lose precision during the addition...
double s[3][N_GRAVS], vs[3][N_GRAVS], mass[N_GRAVS];
#ifdef NGRAVS_ACCUMULATOR
long Nparticles[N_GRAVS];
#endif
if(no >= All.MaxPart && no < All.MaxPart + MaxNodes) /* internal node */
{
for(j = 0; j < 8; j++)
suns[j] = Nodes[no].u.suns[j]; /* this "backup" is necessary because the nextnode entry will
overwrite one element (union!) */
if(last >= 0)
{
if(last >= All.MaxPart)
{
if(last >= All.MaxPart + MaxNodes) /* a pseudo-particle */
Nextnode[last - MaxNodes] = no;
else
Nodes[last].u.d.nextnode = no;
}
else
Nextnode[last] = no;
}
last = no;
// KC 3/11/14
// Replace memsets, memcpys with explicit looping. 3-9% performance gains!
// Unroll inner loop: 1-2% performance gains! Holy shit.
for(i = 0; i < N_GRAVS; ++i) {
s[0][i] = 0;
s[1][i] = 0;
s[2][i] = 0;
vs[0][i] = 0;
vs[1][i] = 0;
vs[2][i] = 0;
mass[i] = 0;
#ifdef NGRAVS_ACCUMULATOR
Nparticles[i] = 0;
#endif
}
hmax = 0;
#ifdef UNEQUALSOFTENINGS
#ifndef ADAPTIVE_GRAVSOFT_FORGAS
maxsofttype = 7;
diffsoftflag = 0;
#else
maxsoft = 0;
#endif
#endif
for(j = 0; j < 8; j++)
{
if((p = suns[j]) >= 0)
{
/* check if we have a sibling on the same level */
for(jj = j + 1; jj < 8; jj++)
if((pp = suns[jj]) >= 0)
break;
if(jj < 8) /* yes, we do */
nextsib = pp;
else
nextsib = sib;
force_update_node_recursive(p, nextsib, no);
if(p >= All.MaxPart) /* an internal node or pseudo particle */
{
if(p >= All.MaxPart + MaxNodes) /* a pseudo particle */
{
/* nothing to be done here because the mass of the
* pseudo-particle is still zero. This will be changed
* later.
*/
}
else
{
// KC 8/11/14 Loop over interaction types
// Note we cannot just memcpy due to type differences
for(i = 0; i < N_GRAVS; ++i) {
#ifdef NGRAVS_ACCUMULATOR
Nparticles[i] += Nodes[p].u.d.Nparticles[i];
#endif
mass[i] += Nodes[p].u.d.mass[i];
s[0][i] += Nodes[p].u.d.mass[i] * Nodes[p].u.d.s[0][i];
s[1][i] += Nodes[p].u.d.mass[i] * Nodes[p].u.d.s[1][i];
s[2][i] += Nodes[p].u.d.mass[i] * Nodes[p].u.d.s[2][i];
vs[0][i] += Nodes[p].u.d.mass[i] * Extnodes[p].vs[0][i];
vs[1][i] += Nodes[p].u.d.mass[i] * Extnodes[p].vs[1][i];
vs[2][i] += Nodes[p].u.d.mass[i] * Extnodes[p].vs[2][i];
}
if(Extnodes[p].hmax > hmax)
hmax = Extnodes[p].hmax;
#ifdef UNEQUALSOFTENINGS
#ifndef ADAPTIVE_GRAVSOFT_FORGAS
diffsoftflag |= (Nodes[p].u.d.bitflags >> 5) & 1;
if(maxsofttype == 7)
{
maxsofttype = (Nodes[p].u.d.bitflags >> 2) & 7;
}
else
{
if(((Nodes[p].u.d.bitflags >> 2) & 7) != 7)
{
if(All.ForceSoftening[((Nodes[p].u.d.bitflags >> 2) & 7)] >
All.ForceSoftening[maxsofttype])
{
maxsofttype = ((Nodes[p].u.d.bitflags >> 2) & 7);
diffsoftflag = 1;
}
else
{
if(All.ForceSoftening[((Nodes[p].u.d.bitflags >> 2) & 7)] <
All.ForceSoftening[maxsofttype])
diffsoftflag = 1;
}
}
}
#else
if(Nodes[p].maxsoft > maxsoft)
maxsoft = Nodes[p].maxsoft;
#endif
#endif
}
}
else /* a particle */
{
pa = &P[p];
i = TypeToGrav[pa->Type];
mass[i] += pa->Mass;
s[0][i] += pa->Mass * pa->Pos[0];
s[1][i] += pa->Mass * pa->Pos[1];
s[2][i] += pa->Mass * pa->Pos[2];
vs[0][i] += pa->Mass * pa->Vel[0];
vs[1][i] += pa->Mass * pa->Vel[1];
vs[2][i] += pa->Mass * pa->Vel[2];
#ifdef NGRAVS_ACCUMULATOR
Nparticles[i]++;
#endif
#ifdef UNEQUALSOFTENINGS
#ifndef ADAPTIVE_GRAVSOFT_FORGAS
if(maxsofttype == 7)
{
maxsofttype = pa->Type;
}
else
{
if(All.ForceSoftening[pa->Type] > All.ForceSoftening[maxsofttype])
{
maxsofttype = pa->Type;
diffsoftflag = 1;
}
else
{
if(All.ForceSoftening[pa->Type] < All.ForceSoftening[maxsofttype])
diffsoftflag = 1;
}
}
#else
if(pa->Type == 0)
{
if(SphP[p].Hsml > maxsoft)
maxsoft = SphP[p].Hsml;
}
else
{
if(All.ForceSoftening[pa->Type] > maxsoft)
maxsoft = All.ForceSoftening[pa->Type];
}
#endif
#endif
if(pa->Type == 0)
if(SphP[p].Hsml > hmax)
hmax = SphP[p].Hsml;
}
}
}
// KC 8/11/14 Looks like we should be looping over this stuff
// KC 2/7/15 PPP (perform a single division, then the others are multuplications)
// but Volker does this too, so this can't be so bad?
for(i = 0; i < N_GRAVS; ++i) {
if(mass[i] > 0)
{
s[0][i] /= mass[i];
s[1][i] /= mass[i];
s[2][i] /= mass[i];
vs[0][i] /= mass[i];
vs[1][i] /= mass[i];
vs[2][i] /= mass[i];
}
else
{
s[0][i] = Nodes[no].center[0];
s[1][i] = Nodes[no].center[1];
s[2][i] = Nodes[no].center[2];
}
Nodes[no].u.d.s[0][i] = s[0][i];
Nodes[no].u.d.s[1][i] = s[1][i];
Nodes[no].u.d.s[2][i] = s[2][i];
Nodes[no].u.d.mass[i] = mass[i];
#ifdef NGRAVS_ACCUMULATOR
Nodes[no].u.d.Nparticles[i] = Nparticles[i];
#endif
// KC 8/11/14 We move the vs up here as we need to loop to set all these things
Extnodes[no].vs[0][i] = vs[0][i];
Extnodes[no].vs[1][i] = vs[1][i];
Extnodes[no].vs[2][i] = vs[2][i];
#ifdef NGRAVS_ACCUMULATOR_DEBUG
printf("%d: {%g, %g, %g} --> %g = %d\n", i, Nodes[no].center[0], Nodes[no].center[1], Nodes[no].center[2], Nodes[no].len, Nparticles[i]);
#endif
}
#ifdef UNEQUALSOFTENINGS
#ifndef ADAPTIVE_GRAVSOFT_FORGAS
Nodes[no].u.d.bitflags = 4 * maxsofttype + 32 * diffsoftflag;
#else
Nodes[no].u.d.bitflags = 0;
Nodes[no].maxsoft = maxsoft;
#endif
#else
Nodes[no].u.d.bitflags = 0;
#endif
// KC 8/11/14 This is where the vs used to be
Extnodes[no].hmax = hmax;
Nodes[no].u.d.sibling = sib;
Nodes[no].u.d.father = father;
}
else /* single particle or pseudo particle */
{
if(last >= 0)
{
if(last >= All.MaxPart)
{
if(last >= All.MaxPart + MaxNodes) /* a pseudo-particle */
Nextnode[last - MaxNodes] = no;
else
Nodes[last].u.d.nextnode = no;
}
else
Nextnode[last] = no;
}
last = no;
if(no < All.MaxPart) /* only set it for single particles */
Father[no] = father;
}
}
/*! This function updates the multipole moments of the pseudo-particles
* that represent the mass distribution on different CPUs. For that
* purpose, it first exchanges the necessary data, and then updates the
* top-level tree accordingly. The detailed implementation of these two
* tasks is done in separate functions.
*/
void force_update_pseudoparticles(void)
{
force_exchange_pseudodata();
force_treeupdate_pseudos();
}
/*! This function communicates the values of the multipole moments of the
* top-level tree-nodes of the domain grid. This data can then be used to
* update the pseudo-particles on each CPU accordingly.
*/
void force_exchange_pseudodata(void)
{
int i, no;
MPI_Status status;
int level, sendTask, recvTask;
int n, k;
for(i = DomainMyStart; i <= DomainMyLast; i++)
{
no = DomainNodeIndex[i];
/* read out the multipole moments from the local base cells */
// KC 11/2/14
// Unroll like a beast
for(n = 0; n < N_GRAVS; ++n) {
DomainMoment[i].s[0][n] = Nodes[no].u.d.s[0][n];
DomainMoment[i].s[1][n] = Nodes[no].u.d.s[1][n];
DomainMoment[i].s[2][n] = Nodes[no].u.d.s[2][n];
DomainMoment[i].mass[n] = Nodes[no].u.d.mass[n];
#ifdef NGRAVS_ACCUMULATOR
DomainMoment[i].Nparticles[n] = Nodes[no].u.d.Nparticles[n];
#endif
DomainMoment[i].vs[0][n] = Extnodes[no].vs[0][n];
DomainMoment[i].vs[1][n] = Extnodes[no].vs[1][n];
DomainMoment[i].vs[2][n] = Extnodes[no].vs[2][n];
}
#ifdef UNEQUALSOFTENINGS
#ifndef ADAPTIVE_GRAVSOFT_FORGAS
DomainMoment[i].bitflags = Nodes[no].u.d.bitflags;
#else
DomainMoment[i].maxsoft = Nodes[no].maxsoft;
#endif
#endif
}
/* share the pseudo-particle data accross CPUs */
for(level = 1; level < (1 << PTask); level++)
{
sendTask = ThisTask;
recvTask = ThisTask ^ level;
if(recvTask < NTask)
MPI_Sendrecv(&DomainMoment[DomainStartList[sendTask]],
(DomainEndList[sendTask] - DomainStartList[sendTask] + 1) * sizeof(struct DomainNODE),
MPI_BYTE, recvTask, TAG_DMOM,
&DomainMoment[DomainStartList[recvTask]],
(DomainEndList[recvTask] - DomainStartList[recvTask] + 1) * sizeof(struct DomainNODE),
MPI_BYTE, recvTask, TAG_DMOM, MPI_COMM_WORLD, &status);
}
}
/*! This function updates the top-level tree after the multipole moments of
* the pseudo-particles have been updated.
*/
void force_treeupdate_pseudos(void)
{
int i, k, no;
// KC 8/11/14 Extension of the variables
FLOAT sold[3][N_GRAVS], vsold[3][N_GRAVS], snew[3][N_GRAVS], vsnew[3][N_GRAVS], massold[N_GRAVS], massnew[N_GRAVS], mm[N_GRAVS];
#ifdef UNEQUALSOFTENINGS
#ifndef ADAPTIVE_GRAVSOFT_FORGAS
int maxsofttype, diffsoftflag;
#else
FLOAT maxsoft;
#endif
#endif
int j;
#ifdef NGRAVS_ACCUMULATOR
int partnew[N_GRAVS], partold[N_GRAVS];
#endif
for(i = 0; i < NTopleaves; i++)
if(i < DomainMyStart || i > DomainMyLast)
{
no = DomainNodeIndex[i];
// KC 11/4/14
// No need to unroll here because we keep hitting different variables each time
// so little cache gain?
//
// He does keep it so that DomainMoment occurs near in the code
for(j = 0; j < N_GRAVS; j++) {
for(k = 0; k < 3; ++k) {
vsold[k][j] = Extnodes[no].vs[k][j];
sold[k][j] = Nodes[no].u.d.s[k][j];
}
massold[j] = Nodes[no].u.d.mass[j];
#ifdef NGRAVS_ACCUMULATOR
partold[j] = Nodes[no].u.d.Nparticles[j];
#endif
for(k = 0; k < 3; ++k) {
snew[k][j] = DomainMoment[i].s[k][j];
vsnew[k][j] = DomainMoment[i].vs[k][j];
}
massnew[j] = DomainMoment[i].mass[j];
#ifdef NGRAVS_ACCUMULATOR
partnew[j] = DomainMoment[i].Nparticles[j];
#endif
}
#ifdef UNEQUALSOFTENINGS
#ifndef ADAPTIVE_GRAVSOFT_FORGAS
maxsofttype = (DomainMoment[i].bitflags >> 2) & 7;
diffsoftflag = (DomainMoment[i].bitflags >> 5) & 1;
#else
maxsoft = DomainMoment[i].maxsoft;
#endif
#endif
do
{
// KC 11/4/14
// This should correctly compute mass centers, as long as the values accumulated
// in mass[i] represent asymptotic masses of the sourcing entities
for(j = 0; j < N_GRAVS; ++j) {
#ifdef NGRAVS_ACCUMULATOR
Nodes[no].u.d.Nparticles[j] += partnew[j] - partold[j];
#endif
mm[j] = Nodes[no].u.d.mass[j] + massnew[j] - massold[j];
for(k = 0; k < 3; k++)
{
if(mm[j] > 0)
{
Nodes[no].u.d.s[k][j] =
(Nodes[no].u.d.mass[j] * Nodes[no].u.d.s[k][j] +
massnew[j] * snew[k][j] -
massold[j] * sold[k][j]) / mm[j];
Extnodes[no].vs[k][j] =
(Nodes[no].u.d.mass[j] * Extnodes[no].vs[k][j] +
massnew[j] * vsnew[k][j] -
massold[j] * vsold[k][j]) / mm[j];
}
}
Nodes[no].u.d.mass[j] = mm[j];
}
#ifdef UNEQUALSOFTENINGS
#ifndef ADAPTIVE_GRAVSOFT_FORGAS
diffsoftflag |= (Nodes[no].u.d.bitflags >> 5) & 1;
if(maxsofttype == 7)
maxsofttype = (Nodes[no].u.d.bitflags >> 2) & 7;
else
{
if(((Nodes[no].u.d.bitflags >> 2) & 7) != 7)
{
if(All.ForceSoftening[((Nodes[no].u.d.bitflags >> 2) & 7)] >
All.ForceSoftening[maxsofttype])
{
maxsofttype = ((Nodes[no].u.d.bitflags >> 2) & 7);
diffsoftflag = 1;
}
else
{
if(All.ForceSoftening[((Nodes[no].u.d.bitflags >> 2) & 7)] <
All.ForceSoftening[maxsofttype])
diffsoftflag = 1;
}
}
}
Nodes[no].u.d.bitflags = (Nodes[no].u.d.bitflags & 3) + 4 * maxsofttype + 32 * diffsoftflag;
#else
if(Nodes[no].maxsoft < maxsoft)
Nodes[no].maxsoft = maxsoft;
maxsoft = Nodes[no].maxsoft;
#endif
#endif
no = Nodes[no].u.d.father;
}
while(no >= 0);
}
}
/*! This function flags nodes in the top-level tree that are dependent on
* local particle data.
*/
void force_flag_localnodes(void)
{
int no, i;
/* mark all top-level nodes */
for(i = 0; i < NTopleaves; i++)
{
no = DomainNodeIndex[i];
while(no >= 0)
{
if((Nodes[no].u.d.bitflags & 1))
break;
Nodes[no].u.d.bitflags |= 1;
no = Nodes[no].u.d.father;
}
}
/* mark top-level nodes that contain local particles */
for(i = DomainMyStart; i <= DomainMyLast; i++)
{
/*
if(DomainMoment[i].mass > 0)
*/
{
no = DomainNodeIndex[i];
while(no >= 0)
{
if((Nodes[no].u.d.bitflags & 2))
break;
Nodes[no].u.d.bitflags |= 2;
no = Nodes[no].u.d.father;
}
}
}
}
/*! This function updates the side-length of tree nodes in case the tree is