12#include "gadgetconfig.h"
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"
32template <
typename partset>
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));
39 for(
int k = 0; k < n; k++)
42 TopNodes[i].Daughter = NTopnodes;
46 for(
int j = 0; j < 8; j++)
48 int sub = TopNodes[i].Daughter + j;
50 TopNodes[sub].Daughter = -1;
51 topNodes[sub].Level = topNodes[i].Level + 1;
53 topNodes[sub].PIndex = topNodes[i].PIndex;
54 topNodes[sub].Count = 0;
55 topNodes[sub].Cost = 0;
58 int sub = TopNodes[i].Daughter;
60 for(
int p = topNodes[i].PIndex, j = 0; p < topNodes[i].PIndex + topNodes[i].Count; p++)
63 while(mp[p].key >= topNodes[sub + 1].StartKey)
67 topNodes[sub].PIndex = p;
72 topNodes[sub].Cost += mp[p].cost;
73 topNodes[sub].Count++;
76 for(
int j = 0; j < 8; j++)
78 int sub = TopNodes[i].Daughter + j;
79 worklist[k * 8 + j] = topNodes[sub].Cost;
80 countlist[k * 8 + j] = topNodes[sub].Count;
84 allreduce_sum<double>(worklist, 8 * n, Communicator);
85 allreduce_sum<long long>(countlist, 8 * n, Communicator);
88 for(
int k = 0; k < n; k++)
92 for(
int j = 0; j < 8; j++)
94 int sub = TopNodes[i].Daughter + j;
95 topNodes[sub].Cost = worklist[k * 8 + j];
96 topNodes[sub].CountTot = countlist[k * 8 + j];
100 Mem.myfree(countlist);
101 Mem.myfree(worklist);
109template <
typename partset>
112 if(TopNodes[no].Daughter == -1)
114 TopNodes[no].Leaf = NTopleaves;
116 domain_leaf_cost[TopNodes[no].Leaf].Cost = topNodes[no].Cost;
122 for(
int i = 0; i < 8; i++)
123 domain_walktoptree(TopNodes[no].Daughter + i);
132 for(
int n = 0; n < NumTimeBinsToBeBalanced; n++)
134 int bin = ListOfTimeBinsToBeBalanced[n];
136 if(bin >= Tp->P[i].TimeBinGrav)
137 cost += GravCostNormFactors[n] * Tp->P[i].GravCost;
139 if(Tp->P[i].getType() == 0)
144 if(Tp->PS[i].DomainFlag)
145 cost += HydroCostNormFactors[n];
150 if(bin >= Tp->P[i].getTimeBinHydro())
151 cost += HydroCostNormFactors[n];
159#ifdef LIGHTCONE_PARTICLES
173template <
typename partset>
177 int message_printed = 0;
179 mp = (domain_peano_hilbert_data *)
Mem.mymalloc_movable(&mp,
"mp",
sizeof(domain_peano_hilbert_data) * Tp->NumPart);
184 for(
int i = 0; i < Tp->NumPart; i++)
190 mp[i].cost += domain_get_cost_summed_over_timebins(i);
192 mp[i].cost += NormFactorLoad;
194 if(Tp->P[i].getType() == 0)
195 mp[i].cost += NormFactorLoadSph;
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);
204 mycxxsort(mp, mp + Tp->NumPart, domain_compare_key);
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;
213 topNodes[0].CountTot = Tp->TotNumPart;
214 topNodes[0].Cost = sum;
217 int *list = (
int *)
Mem.mymalloc_movable(&list,
"list", MaxTopNodes *
sizeof(
int));
227 for(
int n = 0; n < NTopnodes; n++)
228 if(TopNodes[n].Daughter == -1)
230 if(topNodes[n].CountTot > 1)
231 if(topNodes[n].Cost > limit)
233 while(NTopnodes + 8 * (count + 1) > MaxTopNodes)
239 Terminate(
"something seems to be going seriously wrong here. Stopping.\n");
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)));
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));
253 if(message_printed == 0)
256 "DOMAIN: Note: we would like to refine top-tree beyond the level allowed by the selected positional "
271 domain_do_local_refine(count, list);
282 domain_walktoptree(0);
286 domain_printf(
"DOMAIN: NTopleaves=%d, determination of top-level tree involved %d iterations and took %g sec\n", NTopleaves, iter,
290#include "../data/simparticles.h"
293#ifdef LIGHTCONE_PARTICLES
294#include "../data/lcparticles.h"
global_data_all_processes All
double timediff(double t0, double t1)
double mycxxsort(T *begin, T *end, Tcomp comp)
peanokey get_peanokey_offset(unsigned int j, int bits)
#define BITS_FOR_POSITIONS
peanokey peano_hilbert_key(MyIntPosType x, MyIntPosType y, MyIntPosType z, int bits)
double TopNodeAllocFactor