GADGET-4
domain_toplevel.cc
Go to the documentation of this file.
1 /*******************************************************************************
2  * \copyright This file is part of the GADGET4 N-body/SPH code developed
3  * \copyright by Volker Springel. Copyright (C) 2014-2020 by Volker Springel
4  * \copyright (vspringel@mpa-garching.mpg.de) and all contributing authors.
5  *******************************************************************************/
6 
12 #include "gadgetconfig.h"
13 
14 #include <mpi.h>
15 #include <algorithm>
16 #include <cmath>
17 #include <cstdio>
18 #include <cstdlib>
19 #include <cstring>
20 
21 #include "../data/allvars.h"
22 #include "../data/dtypes.h"
23 #include "../data/mymalloc.h"
24 #include "../domain/domain.h"
25 #include "../logs/timer.h"
26 #include "../main/simulation.h"
27 #include "../mpi_utils/mpi_utils.h"
28 #include "../sort/cxxsort.h"
29 #include "../sort/peano.h"
30 #include "../system/system.h"
31 
32 template <typename partset>
33 void domain<partset>::domain_do_local_refine(int n, int *list)
34 {
35  double *worklist = (double *)Mem.mymalloc("worklist", 8 * n * sizeof(double));
36  long long *countlist = (long long *)Mem.mymalloc("countlist", 8 * n * sizeof(long long));
37 
38  /* create the new nodes */
39  for(int k = 0; k < n; k++)
40  {
41  int i = list[k];
42  TopNodes[i].Daughter = NTopnodes;
43  NTopnodes += 8;
44  NTopleaves += 7;
45 
46  for(int j = 0; j < 8; j++)
47  {
48  int sub = TopNodes[i].Daughter + j;
49 
50  TopNodes[sub].Daughter = -1;
51  topNodes[sub].Level = topNodes[i].Level + 1;
52  topNodes[sub].StartKey = topNodes[i].StartKey + get_peanokey_offset(j, (3 * (BITS_FOR_POSITIONS - topNodes[sub].Level)));
53  topNodes[sub].PIndex = topNodes[i].PIndex;
54  topNodes[sub].Count = 0;
55  topNodes[sub].Cost = 0;
56  }
57 
58  int sub = TopNodes[i].Daughter;
59 
60  for(int p = topNodes[i].PIndex, j = 0; p < topNodes[i].PIndex + topNodes[i].Count; p++)
61  {
62  if(j < 7)
63  while(mp[p].key >= topNodes[sub + 1].StartKey)
64  {
65  j++;
66  sub++;
67  topNodes[sub].PIndex = p;
68  if(j >= 7)
69  break;
70  }
71 
72  topNodes[sub].Cost += mp[p].cost;
73  topNodes[sub].Count++;
74  }
75 
76  for(int j = 0; j < 8; j++)
77  {
78  int sub = TopNodes[i].Daughter + j;
79  worklist[k * 8 + j] = topNodes[sub].Cost;
80  countlist[k * 8 + j] = topNodes[sub].Count;
81  }
82  }
83 
84  allreduce_sum<double>(worklist, 8 * n, Communicator);
85  allreduce_sum<long long>(countlist, 8 * n, Communicator);
86 
87  /* store the results in the corresponding top nodes */
88  for(int k = 0; k < n; k++)
89  {
90  int i = list[k];
91 
92  for(int j = 0; j < 8; j++)
93  {
94  int sub = TopNodes[i].Daughter + j;
95  topNodes[sub].Cost = worklist[k * 8 + j];
96  topNodes[sub].CountTot = countlist[k * 8 + j];
97  }
98  }
99 
100  Mem.myfree(countlist);
101  Mem.myfree(worklist);
102 }
103 
109 template <typename partset>
111 {
112  if(TopNodes[no].Daughter == -1)
113  {
114  TopNodes[no].Leaf = NTopleaves;
115 
116  domain_leaf_cost[TopNodes[no].Leaf].Cost = topNodes[no].Cost;
117 
118  NTopleaves++;
119  }
120  else
121  {
122  for(int i = 0; i < 8; i++)
123  domain_walktoptree(TopNodes[no].Daughter + i);
124  }
125 }
126 
127 template <>
129 {
130  double cost = 0;
131 
132  for(int n = 0; n < NumTimeBinsToBeBalanced; n++)
133  {
134  int bin = ListOfTimeBinsToBeBalanced[n];
135 
136  if(bin >= Tp->P[i].TimeBinGrav)
137  cost += GravCostNormFactors[n] * Tp->P[i].GravCost;
138 
139  if(Tp->P[i].getType() == 0)
140  {
141 #ifdef SUBFIND
142  if(Mode == COLL_SUBFIND)
143  {
144  if(Tp->PS[i].DomainFlag)
145  cost += HydroCostNormFactors[n];
146  }
147  else
148 #endif
149  {
150  if(bin >= Tp->P[i].getTimeBinHydro())
151  cost += HydroCostNormFactors[n];
152  }
153  }
154  }
155 
156  return cost;
157 }
158 
159 #ifdef LIGHTCONE_PARTICLES
160 template <>
162 {
163  return 0;
164 }
165 #endif
166 
173 template <typename partset>
175 {
176  double t0 = Logs.second();
177  int message_printed = 0;
178 
179  mp = (domain_peano_hilbert_data *)Mem.mymalloc_movable(&mp, "mp", sizeof(domain_peano_hilbert_data) * Tp->NumPart);
180 
181  int count = 0;
182  double sum = 0.0;
183 
184  for(int i = 0; i < Tp->NumPart; i++)
185  {
186  mp[i].key = domain_key[i] = peano_hilbert_key(Tp->P[i].IntPos[0], Tp->P[i].IntPos[1], Tp->P[i].IntPos[2], BITS_FOR_POSITIONS);
187  mp[i].cost = 0;
188  count++;
189 
190  mp[i].cost += domain_get_cost_summed_over_timebins(i);
191 
192  mp[i].cost += NormFactorLoad;
193 
194  if(Tp->P[i].getType() == 0)
195  mp[i].cost += NormFactorLoadSph;
196 
197  sum += mp[i].cost;
198  }
199 
200  MPI_Allreduce(MPI_IN_PLACE, &sum, 1, MPI_DOUBLE, MPI_SUM, Communicator);
201  domain_printf("DOMAIN: Sum=%g TotalCost=%g NumTimeBinsToBeBalanced=%d MultipleDomains=%d\n", sum, TotalCost,
202  NumTimeBinsToBeBalanced, MultipleDomains);
203 
204  mycxxsort(mp, mp + Tp->NumPart, domain_compare_key);
205 
206  NTopnodes = 1;
207  NTopleaves = 1;
208  TopNodes[0].Daughter = -1;
209  topNodes[0].Level = 0;
210  topNodes[0].StartKey = {0, 0, 0};
211  topNodes[0].PIndex = 0;
212  topNodes[0].Count = count; /* this is the local number */
213  topNodes[0].CountTot = Tp->TotNumPart;
214  topNodes[0].Cost = sum;
215 
216  /* in list[], we store the node indices hat should be refined */
217  int *list = (int *)Mem.mymalloc_movable(&list, "list", MaxTopNodes * sizeof(int));
218 
219  double limit = 1.0 / (All.TopNodeFactor * NTask);
220 
221  int iter = 0;
222 
223  do
224  {
225  count = 0;
226 
227  for(int n = 0; n < NTopnodes; n++)
228  if(TopNodes[n].Daughter == -1) // consider only leaf nodes
229  {
230  if(topNodes[n].CountTot > 1) // only split nodes with at list 1 particles
231  if(topNodes[n].Cost > limit)
232  {
233  while(NTopnodes + 8 * (count + 1) > MaxTopNodes)
234  {
235  domain_printf("DOMAIN: Increasing TopNodeAllocFactor=%g ", All.TopNodeAllocFactor);
236  All.TopNodeAllocFactor *= 1.3;
237  domain_printf("new value=%g\n", All.TopNodeAllocFactor);
238  if(All.TopNodeAllocFactor > 1000)
239  Terminate("something seems to be going seriously wrong here. Stopping.\n");
240 
241  MaxTopNodes = All.TopNodeAllocFactor * std::max<int>(All.TopNodeFactor * MultipleDomains * NTask, BASENUMBER);
242 
243  topNodes = (local_topnode_data *)Mem.myrealloc_movable(topNodes, (MaxTopNodes * sizeof(local_topnode_data)));
244  TopNodes = (topnode_data *)Mem.myrealloc_movable(TopNodes, (MaxTopNodes * sizeof(topnode_data)));
245  TaskOfLeaf = (int *)Mem.myrealloc_movable(TaskOfLeaf, (MaxTopNodes * sizeof(int)));
246  domain_leaf_cost =
247  (domain_cost_data *)Mem.myrealloc_movable(domain_leaf_cost, (MaxTopNodes * sizeof(domain_cost_data)));
248  list = (int *)Mem.myrealloc_movable(list, MaxTopNodes * sizeof(int));
249  }
250 
251  if(topNodes[n].Level >= BITS_FOR_POSITIONS - 2)
252  {
253  if(message_printed == 0)
254  {
255  domain_printf(
256  "DOMAIN: Note: we would like to refine top-tree beyond the level allowed by the selected positional "
257  "accuracy.\n");
258  message_printed = 1;
259  }
260  }
261  else
262  {
263  list[count] = n;
264  count++;
265  }
266  }
267  }
268 
269  if(count > 0)
270  {
271  domain_do_local_refine(count, list);
272  iter++;
273  }
274  }
275  while(count > 0);
276 
277  Mem.myfree(list);
278  Mem.myfree(mp);
279 
280  /* count the number of top leaves */
281  NTopleaves = 0;
282  domain_walktoptree(0);
283 
284  double t1 = Logs.second();
285 
286  domain_printf("DOMAIN: NTopleaves=%d, determination of top-level tree involved %d iterations and took %g sec\n", NTopleaves, iter,
287  Logs.timediff(t0, t1));
288 }
289 
290 #include "../data/simparticles.h"
291 template class domain<simparticles>;
292 
293 #ifdef LIGHTCONE_PARTICLES
294 #include "../data/lcparticles.h"
295 template class domain<lcparticles>;
296 #endif
global_data_all_processes All
Definition: main.cc:40
Definition: domain.h:31
domain_options Mode
Definition: domain.h:50
double timediff(double t0, double t1)
Definition: logs.cc:488
double second(void)
Definition: logs.cc:471
#define BASENUMBER
Definition: constants.h:303
double mycxxsort(T *begin, T *end, Tcomp comp)
Definition: cxxsort.h:39
@ COLL_SUBFIND
Definition: domain.h:24
peanokey get_peanokey_offset(unsigned int j, int bits)
Definition: dtypes.h:270
#define BITS_FOR_POSITIONS
Definition: dtypes.h:37
logs Logs
Definition: main.cc:43
#define Terminate(...)
Definition: macros.h:19
memory Mem
Definition: main.cc:44
peanokey peano_hilbert_key(MyIntPosType x, MyIntPosType y, MyIntPosType z, int bits)
Definition: peano.cc:83