GADGET-4
subfind_unbind.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 #ifdef SUBFIND
15 
16 #include <gsl/gsl_math.h>
17 #include <mpi.h>
18 #include <algorithm>
19 #include <cmath>
20 #include <cstdio>
21 #include <cstdlib>
22 #include <cstring>
23 
24 #include "../data/allvars.h"
25 #include "../data/dtypes.h"
26 #include "../data/intposconvert.h"
27 #include "../data/mymalloc.h"
28 #include "../domain/domain.h"
29 #include "../fof/fof.h"
30 #include "../gravtree/gravtree.h"
31 #include "../logs/timer.h"
32 #include "../main/simulation.h"
33 #include "../mpi_utils/mpi_utils.h"
34 #include "../sort/cxxsort.h"
35 #include "../sort/parallel_sort.h"
36 #include "../sort/peano.h"
37 #include "../subfind/subfind.h"
38 #include "../system/system.h"
39 
40 #define MAX_UNBOUND_FRAC_BEFORE_BULK_VELOCITY_UPDATE 0.02
41 #define MAX_UNBOUND_FRAC_BEFORE_POTENTIAL_UPDATE 0.20
42 
43 /* this function takes a list of particles given via their indices in d[] and subjects them
44  * to a gravitational unbinding procedure. The number of bound particles is returned,
45  * and the array d[] is updated accordingly.
46  */
47 template <typename partset>
48 int fof<partset>::subfind_unbind(domain<partset> *D, MPI_Comm Communicator, int *d, int num)
49 {
50  double fac_vel_to_phys, fac_hubbleflow, fac_comov_to_phys;
51  subfind_get_factors(fac_vel_to_phys, fac_hubbleflow, fac_comov_to_phys);
52 
53  /* get the local communication context */
54  int commNTask, commThisTask;
55  MPI_Comm_size(Communicator, &commNTask);
56  MPI_Comm_rank(Communicator, &commThisTask);
57 
58  typename partset::pdata *P = Tp->P;
59  subfind_data *PS = Tp->PS;
60 
61  /* we will start out by recomputing the potential for all particles based on all particles */
62  int phaseflag = RECOMPUTE_ALL;
63 
64  int iter = 0;
65 
66  long long totnum = num;
67  MPI_Allreduce(MPI_IN_PLACE, &totnum, 1, MPI_LONG_LONG, MPI_SUM, Communicator);
68 
69  int num_removed = 0;
70  long long totremoved;
71 
72  int *dremoved = (int *)Mem.mymalloc("dremoved", num * sizeof(int));
73  double *potold = (double *)Mem.mymalloc("potold", num * sizeof(double));
74 
75  do
76  {
77  FoFGravTree.treeallocate(Tp->NumPart, Tp, D);
78 
79  if(phaseflag == RECOMPUTE_ALL)
80  {
81  FoFGravTree.treebuild(num, d);
82  }
83  else
84  {
85  FoFGravTree.treebuild(num_removed, dremoved);
86 
87  for(int i = 0; i < num; i++)
88  potold[i] = PS[d[i]].u.s.u.DM_Potential;
89  }
90 
91  /* let's compute the potential energy */
92  subfind_potential_compute(D, num, d);
93 
94  FoFGravTree.treefree();
95 
96  if(phaseflag == RECOMPUTE_ALL)
97  {
98  /* subtract self-potential and convert to physical potential */
99  for(int i = 0; i < num; i++)
100  {
101  int softtype = P[d[i]].getSofteningClass();
102 
103  PS[d[i]].u.s.u.DM_Potential += Tp->P[d[i]].getMass() / (All.ForceSoftening[softtype] / 2.8);
104  PS[d[i]].u.s.u.DM_Potential *= All.G / fac_comov_to_phys;
105  }
106  }
107  else
108  {
109  /* do not correct for self-potential and instead use calculated potential as correction to previous potential */
110  for(int i = 0; i < num; i++)
111  {
112  PS[d[i]].u.s.u.DM_Potential *= All.G / fac_comov_to_phys;
113  PS[d[i]].u.s.u.DM_Potential = potold[i] - PS[d[i]].u.s.u.DM_Potential;
114  }
115  }
116 
117  /* At this point, we have in num/d[] the list of still considered particles, and the potential is current */
118 
119  /* we will not unbind all particles with positive energy, until none are left. Because the kinetic energy depends on
120  * the velocity frame, we recompute the bulk velocity and the center of mass after 2% of the particles have been removed
121  */
122 
123  /* Determine in intpos[] the potential minimum among the particles,
124  * which we take that as the center of the halo
125  */
126  int minindex = -1;
127  struct
128  {
129  double pot;
130  int rank;
131  } local = {MAX_DOUBLE_NUMBER, commThisTask}, global;
132 
133  for(int i = 0; i < num; i++)
134  if(PS[d[i]].u.s.u.DM_Potential < local.pot)
135  {
136  local.pot = PS[d[i]].u.s.u.DM_Potential;
137  minindex = d[i];
138  }
139 
140  MPI_Allreduce(&local, &global, 1, MPI_DOUBLE_INT, MPI_MINLOC, Communicator);
141 
142  MyIntPosType intpos[3]; /* potential minimum */
143 
144  if(commThisTask == global.rank)
145  for(int j = 0; j < 3; j++)
146  intpos[j] = P[minindex].IntPos[j];
147 
148  MPI_Bcast(intpos, 3 * sizeof(MyIntPosType), MPI_BYTE, global.rank, Communicator);
149 
150  /* we start with zero removed particles */
151  num_removed = 0;
152 
153  long long totunbound;
154  do
155  {
156  /* let's get bulk velocity and the center-of-mass */
157  double massloc = 0, sloc[3] = {0, 0, 0}, vloc[3] = {0, 0, 0};
158 
159  for(int i = 0; i < num; i++)
160  {
161  int part_index = d[i];
162 
163  double dxyz[3];
164  Tp->nearest_image_intpos_to_pos(P[part_index].IntPos, intpos, dxyz);
165 
166  for(int j = 0; j < 3; j++)
167  sloc[j] += Tp->P[part_index].getMass() * dxyz[j];
168 
169  for(int j = 0; j < 3; j++)
170  vloc[j] += Tp->P[part_index].getMass() * P[part_index].Vel[j];
171 
172  massloc += Tp->P[part_index].getMass();
173  }
174 
175  double s[3], v[3], mass;
176  MPI_Allreduce(sloc, s, 3, MPI_DOUBLE, MPI_SUM, Communicator);
177  MPI_Allreduce(vloc, v, 3, MPI_DOUBLE, MPI_SUM, Communicator);
178  MPI_Allreduce(&massloc, &mass, 1, MPI_DOUBLE, MPI_SUM, Communicator);
179 
180  for(int j = 0; j < 3; j++)
181  {
182  s[j] /= mass; /* center of mass offset relative to minimum potential */
183  v[j] /= mass; /* center of mass velocity */
184  }
185 
186  MySignedIntPosType off[3];
187  Tp->pos_to_signedintpos(s, off);
188 
189  /* get integer version of absolute center of mass position */
190  MyIntPosType int_cm[3];
191  int_cm[0] = off[0] + intpos[0];
192  int_cm[1] = off[1] + intpos[1];
193  int_cm[2] = off[2] + intpos[2];
194 
195  double *bnd_energy = (double *)Mem.mymalloc("bnd_energy", num * sizeof(double));
196 
197  /* calculate the binding energies */
198  for(int i = 0; i < num; i++)
199  {
200  int part_index = d[i];
201 
202  /* distance to center of mass */
203  double dx[3];
204  Tp->nearest_image_intpos_to_pos(P[part_index].IntPos, int_cm, dx);
205 
206  /* get physical velocity relative to center of mass */
207  double dv[3];
208  for(int j = 0; j < 3; j++)
209  {
210  dv[j] = fac_vel_to_phys * (P[part_index].Vel[j] - v[j]);
211  dv[j] += fac_hubbleflow * fac_comov_to_phys * dx[j];
212  }
213 
214  PS[part_index].v.DM_BindingEnergy =
215  PS[part_index].u.s.u.DM_Potential + 0.5 * (dv[0] * dv[0] + dv[1] * dv[1] + dv[2] * dv[2]);
216 #ifndef LEAN
217  if(P[part_index].getType() == 0)
218  PS[part_index].v.DM_BindingEnergy += PS[part_index].Utherm;
219 #endif
220  bnd_energy[i] = PS[part_index].v.DM_BindingEnergy;
221  }
222 
223  /* sort by binding energy, highest energies (num_unbound / most weakly bound) first */
224  mycxxsort_parallel(bnd_energy, bnd_energy + num, subfind_compare_binding_energy, Communicator);
225 
226  int *npart = (int *)Mem.mymalloc("npart", commNTask * sizeof(int));
227  MPI_Allgather(&num, 1, MPI_INT, npart, 1, MPI_INT, Communicator);
228 
229  /* (global) index of limiting energy value for least tightly bound fraction - those we may
230  * remove at most in one iteration */
231  long long j = std::max<long long>(5, (long long)(MAX_UNBOUND_FRAC_BEFORE_BULK_VELOCITY_UPDATE * totnum));
232 
233  /* find the processor where this lies */
234  int task = 0;
235  while(j >= npart[task])
236  {
237  j -= npart[task];
238  task++;
239  }
240 
241  double energy_limit = MAX_DOUBLE_NUMBER;
242 
243  if(commThisTask == task)
244  energy_limit = bnd_energy[j];
245 
246  MPI_Allreduce(MPI_IN_PLACE, &energy_limit, 1, MPI_DOUBLE, MPI_MIN, Communicator);
247 
248  /* now unbind particles */
249  int num_unbound = 0;
250 
251  for(int i = 0; i < num; i++)
252  {
253  int p = d[i];
254 
255  if(PS[p].v.DM_BindingEnergy > 0 && PS[p].v.DM_BindingEnergy > energy_limit)
256  {
257  num_unbound++;
258 
259  dremoved[num_removed++] = d[i];
260 
261  d[i] = d[num - 1];
262  num--;
263  i--;
264  }
265  }
266 
267  Mem.myfree(npart);
268  Mem.myfree(bnd_energy);
269 
270  totunbound = num_unbound;
271  totremoved = num_removed;
272  totnum = num;
273 
274  MPI_Allreduce(MPI_IN_PLACE, &totunbound, 1, MPI_LONG_LONG, MPI_SUM, Communicator);
275  MPI_Allreduce(MPI_IN_PLACE, &totremoved, 1, MPI_LONG_LONG, MPI_SUM, Communicator);
276  MPI_Allreduce(MPI_IN_PLACE, &totnum, 1, MPI_LONG_LONG, MPI_SUM, Communicator);
277  }
278  while(totunbound > 0 && totnum >= All.DesLinkNgb && totremoved < MAX_UNBOUND_FRAC_BEFORE_POTENTIAL_UPDATE * totnum);
279 
280  iter++;
281 
282  if(iter > MAX_ITER_UNBIND)
283  Terminate("too many iterations");
284 
285  if(phaseflag == RECOMPUTE_ALL)
286  {
287  if(totremoved > 0)
288  phaseflag = UPDATE_ALL;
289  }
290  else
291  {
292  if(totremoved == 0)
293  {
294  phaseflag = RECOMPUTE_ALL; /* this will make us repeat everything once more for all particles */
295  totremoved = 1; /* to make the code check once more all particles */
296  }
297  }
298  }
299  while(totremoved > 0 && totnum >= All.DesLinkNgb);
300 
301  Mem.myfree(potold);
302  Mem.myfree(dremoved);
303 
304  return num;
305 }
306 
307 /* now make sure that the following classes are really instantiated, otherwise we may get a linking problem */
308 #include "../data/simparticles.h"
309 template class fof<simparticles>;
310 
311 #if defined(LIGHTCONE) && defined(LIGHTCONE_PARTICLES_GROUPS)
312 #include "../data/lcparticles.h"
313 template class fof<lcparticles>;
314 #endif
315 
316 #endif
global_data_all_processes All
Definition: main.cc:40
Definition: domain.h:31
#define MAX_DOUBLE_NUMBER
Definition: constants.h:81
uint32_t MyIntPosType
Definition: dtypes.h:35
int32_t MySignedIntPosType
Definition: dtypes.h:36
#define Terminate(...)
Definition: macros.h:19
memory Mem
Definition: main.cc:44
double mycxxsort_parallel(T *begin, T *end, Comp comp, MPI_Comm comm)
double ForceSoftening[NSOFTCLASSES+NSOFTCLASSES_HYDRO+2]
Definition: allvars.h:258
union subfind_data::@0 u