GADGET-4
subfind_so.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 <math.h>
17 #include <mpi.h>
18 #include <stdio.h>
19 #include <stdlib.h>
20 #include <string.h>
21 
22 #include "../data/allvars.h"
23 #include "../data/dtypes.h"
24 #include "../data/intposconvert.h"
25 #include "../data/mymalloc.h"
26 #include "../domain/domain.h"
27 #include "../fof/fof.h"
28 #include "../gravtree/gravtree.h"
29 #include "../logs/timer.h"
30 #include "../main/simulation.h"
31 #include "../mpi_utils/generic_comm.h"
32 #include "../subfind/subfind.h"
33 #include "../system/system.h"
34 
35 static double *R200, *M200;
36 
37 /* local data structure for collecting particle/cell data that is sent to other processors if needed */
38 struct sodata_in : data_in_generic
39 {
40  MyIntPosType IntPos[3];
41  double R200;
42 };
43 
44 struct sodata_out
45 {
46  double Mass;
47 };
48 
49 template <typename T_tree, typename T_domain, typename T_partset>
50 class sodata_comm : public generic_comm<sodata_in, sodata_out, T_tree, T_domain, T_partset>
51 {
52  public:
54  using gcomm::D;
55  using gcomm::Thread;
56  using gcomm::Tp; // This makes sure that we can access Tp from the base class without having to use "this->Tp"
57  using gcomm::Tree;
58 
59  typename fof<T_partset>::group_properties *Group;
60 
61  /* need to call the base class constructor explicitly */
62  sodata_comm(T_domain *dptr, T_tree *tptr, T_partset *pptr, typename fof<T_partset>::group_properties *grpptr)
63  : gcomm(dptr, tptr, pptr)
64  {
65  Group = grpptr;
66  }
67 
68  /* routine that fills the relevant particle/cell data into the input structure defined above */
69  void particle2in(sodata_in *in, int i) override
70  {
71  in->IntPos[0] = Group[i].IntPos[0];
72  in->IntPos[1] = Group[i].IntPos[1];
73  in->IntPos[2] = Group[i].IntPos[2];
74 
75  in->R200 = R200[i];
76  }
77 
78  /* routine to store or combine result data */
79  void out2particle(sodata_out *out, int i, int mode) override
80  {
81  if(mode == MODE_LOCAL_PARTICLES) /* initial store */
82  M200[i] = out->Mass;
83  else /* combine */
84  M200[i] += out->Mass;
85  }
86 
91  int evaluate(int target, int mode, int thread_id, int action, sodata_in *in, int numnodes, node_info *firstnode,
92  sodata_out &out) override
93  {
94  MyIntPosType *intpos = in->IntPos;
95  double hsml = in->R200;
96 
97  double mass = 0;
98 
99  for(int k = 0; k < numnodes; k++)
100  {
101  int no;
102 
103  if(mode == MODE_LOCAL_PARTICLES)
104  {
105  no = Tree->MaxPart; /* root node */
106  }
107  else
108  {
109  no = firstnode[k].Node;
110  no = Tree->get_nodep(no)->nextnode; /* open it */
111  }
112 
113  unsigned int shmrank = Tree->TreeSharedMem_ThisTask;
114 
115  while(no >= 0)
116  {
117  if(no < Tree->MaxPart) /* single particle */
118  {
119  auto P = Tree->get_Pp(no, shmrank);
120 
121  no = Tree->get_nextnodep(shmrank)[no]; /* note: here shmrank cannot change */
122 
123  double dxyz[3];
124  Tp->nearest_image_intpos_to_pos(P->IntPos, intpos, dxyz); /* converts the integer distance to floating point */
125 
126  double h2 = hsml * hsml;
127 
128  double r2 = dxyz[0] * dxyz[0];
129  if(r2 > h2)
130  continue;
131 
132  r2 += dxyz[1] * dxyz[1];
133  if(r2 > h2)
134  continue;
135 
136  r2 += dxyz[2] * dxyz[2];
137  if(r2 > h2)
138  continue;
139 
140  mass += P->getMass();
141  }
142  else if(no < Tree->MaxPart + Tree->MaxNodes) /* internal node */
143  {
144  if(mode == MODE_IMPORTED_PARTICLES)
145  {
146  /* we reached a top-level node again, which means that we are done with the branch */
147  if(no < Tree->FirstNonTopLevelNode)
148  break;
149  }
150 
151  gravnode *current = Tree->get_nodep(no, shmrank);
152 
153  no = current->sibling; /* in case the node can be discarded */
154  shmrank = current->sibling_shmrank;
155 
156  /* converts the integer distance to floating point */
157  double dxyz[3];
158  Tp->nearest_image_intpos_to_pos(current->center.da, intpos, dxyz);
159 
160  double lenhalf = (((MyIntPosType)1) << (BITS_FOR_POSITIONS - 1 - current->level)) * Tp->FacIntToCoord;
161 
162  double dist = hsml + lenhalf;
163 
164  if(fabs(dxyz[0]) > dist)
165  continue;
166  if(fabs(dxyz[1]) > dist)
167  continue;
168  if(fabs(dxyz[2]) > dist)
169  continue;
170 
171  /* now test against the minimal sphere enclosing everything */
172  dist += FACT1 * 2.0 * lenhalf;
173 
174  double r2 = dxyz[0] * dxyz[0] + dxyz[1] * dxyz[1] + dxyz[2] * dxyz[2];
175 
176  if(r2 > dist * dist)
177  continue;
178 
179  if(no >= Tree->FirstNonTopLevelNode) /* only do this for fully local nodes */
180  {
181  double lenhalf = (((MyIntPosType)1) << (BITS_FOR_POSITIONS - 1 - current->level)) * Tp->FacIntToCoord;
182 
183  /* test whether the node is contained within the sphere */
184  dist = hsml - FACTSQRT3 * lenhalf;
185  if(dist > 0)
186  if(r2 < dist * dist)
187  {
188  mass += current->mass;
189  continue;
190  }
191  }
192 
193  no = current->nextnode; /* ok, we need to open the node */
194  shmrank = current->nextnode_shmrank;
195  }
196  else if(no >= Tree->ImportedNodeOffset) /* point from imported nodelist */
197  {
198  int n = no - Tree->ImportedNodeOffset;
199  no = Tree->Nextnode[no - Tree->MaxNodes];
200  /* note: here shmrank cannot change */
201 
202  double dxyz[3];
203  Tp->nearest_image_intpos_to_pos(Tree->Points[n].IntPos, intpos,
204  dxyz); /* converts the integer distance to floating point */
205 
206  double h2 = hsml * hsml;
207 
208  double r2 = dxyz[0] * dxyz[0];
209  if(r2 > h2)
210  continue;
211 
212  r2 += dxyz[1] * dxyz[1];
213  if(r2 > h2)
214  continue;
215 
216  r2 += dxyz[2] * dxyz[2];
217  if(r2 > h2)
218  continue;
219 
220  mass += Tree->Points[n].Mass;
221  }
222  else /* pseudo particle */
223  {
224  if(mode == MODE_LOCAL_PARTICLES)
225  Tree->tree_export_node_threads(no, target, &Thread);
226 
227  no = Tree->Nextnode[no - Tree->MaxNodes];
228  }
229  }
230  }
231 
232  out.Mass = mass;
233 
234  return 0;
235  }
236 };
237 
238 template <typename partset>
239 double fof<partset>::subfind_get_overdensity_value(int type, double ascale)
240 {
241  double z = 1 / ascale - 1;
242  double omegaz =
243  All.Omega0 * pow(1 + z, 3) / (All.Omega0 * pow(1 + z, 3) + (1 - All.Omega0 - All.OmegaLambda) * pow(1 + z, 2) + All.OmegaLambda);
244  double x = omegaz - 1;
245 
246  if(type == 0)
247  {
248  return 200.0; // Mean200
249  }
250  else if(type == 1)
251  { // Generalized Tophat overdensity
252  return (18 * M_PI * M_PI + 82 * x - 39 * x * x) / omegaz;
253  }
254  else if(type == 2)
255  {
256  return 200.0 / omegaz; // DeltaCrit200
257  }
258  else if(type == 3)
259  {
260  return 500.0 / omegaz; // DeltaCrit500
261  }
262  else
263  Terminate("can't be");
264 
265  return 0;
266 }
267 
268 template <typename partset>
269 double fof<partset>::subfind_overdensity(void)
270 {
271  int *TargetList = (int *)Mem.mymalloc("TargetList", Ngroups * sizeof(int));
272 
273  double *Left = (double *)Mem.mymalloc("Left", sizeof(double) * Ngroups);
274  double *Right = (double *)Mem.mymalloc("Right", sizeof(double) * Ngroups);
275 
276  R200 = (double *)Mem.mymalloc("R200", sizeof(double) * Ngroups);
277  M200 = (double *)Mem.mymalloc("M200", sizeof(double) * Ngroups);
278 
279  double rhoback = 3 * All.Omega0 * All.Hubble * All.Hubble / (8 * M_PI * All.G);
280 
281  double tstart = Logs.second();
282 
283  sodata_comm<gravtree<partset>, domain<partset>, partset> commpattern{FoFDomain, &FoFGravTree, Tp, Group};
284 
285  for(int rep = 0; rep < 4; rep++) /* repeat for all four overdensity values */
286  {
287  int Nso = 0;
288 
289  for(int i = 0; i < Ngroups; i++)
290  {
291  if(Group[i].Nsubs > 0)
292  {
293  double rguess = pow(All.G * Group[i].Mass / (100 * All.Hubble * All.Hubble), 1.0 / 3);
294 
295  TargetList[Nso++] = i;
296 
297  Right[i] = 3 * rguess;
298  Left[i] = 0;
299 
300  R200[i] = 0.5 * (Left[i] + Right[i]);
301  }
302  }
303 
304  int iter = 0;
305  long long ntot;
306 
307  /* we will repeat the whole thing for those groups where we didn't converge to a SO radius yet */
308  do
309  {
310  double t0 = Logs.second();
311 
312  commpattern.execute(Nso, TargetList, MODE_DEFAULT);
313 
314  /* do final operations on results */
315  int npleft = 0;
316  for(int n = 0; n < Nso; n++)
317  {
318  int i = TargetList[n];
319 
320  double overdensity = M200[i] / (4.0 * M_PI / 3.0 * R200[i] * R200[i] * R200[i]) / rhoback;
321 
322  if((Right[i] - Left[i]) > 1.0e-4 * Left[i])
323  {
324  /* need to redo this group */
325  TargetList[npleft++] = i;
326 
327  double delta = subfind_get_overdensity_value(rep, Group[i].Ascale);
328  if(overdensity > delta)
329  Left[i] = R200[i];
330  else
331  Right[i] = R200[i];
332 
333  R200[i] = 0.5 * (Left[i] + Right[i]);
334 
335  if(iter >= MAXITER - 10)
336  {
337  printf("gr=%d task=%d R200=%g Left=%g Right=%g Menclosed=%g Right-Left=%g\n pos=(%g|%g|%g)\n", i, ThisTask,
338  R200[i], Left[i], Right[i], M200[i], Right[i] - Left[i], Group[i].Pos[0], Group[i].Pos[1],
339  Group[i].Pos[2]);
340  myflush(stdout);
341  }
342  }
343  }
344 
345  Nso = npleft;
346 
347  sumup_large_ints(1, &npleft, &ntot, Communicator);
348 
349  double t1 = Logs.second();
350 
351  if(ntot > 0)
352  {
353  iter++;
354 
355  if(iter > 0)
356  mpi_printf("SUBFIND: SO iteration %2d: need to repeat for %12lld halo centers. (took %g sec)\n", iter, ntot,
357  Logs.timediff(t0, t1));
358 
359  if(iter > MAXITER)
360  Terminate("failed to converge in SO iteration");
361  }
362  }
363  while(ntot > 0);
364 
365  for(int i = 0; i < Ngroups; i++)
366  {
367  if(Group[i].Nsubs > 0)
368  {
369  double overdensity = M200[i] / (4.0 * M_PI / 3.0 * R200[i] * R200[i] * R200[i]) / rhoback;
370 
371  double delta = subfind_get_overdensity_value(rep, Group[i].Ascale);
372 
373  if((overdensity - delta) > 0.1 * delta)
374  {
375  R200[i] = M200[i] = 0;
376  }
377  else if(M200[i] < 5 * Group[i].Mass / Group[i].Len)
378  {
379  R200[i] = M200[i] = 0;
380  }
381  }
382  else
383  R200[i] = M200[i] = 0;
384 
385  switch(rep)
386  {
387  case 0:
388  Group[i].M_Mean200 = M200[i];
389  Group[i].R_Mean200 = R200[i];
390  break;
391  case 1:
392  Group[i].M_TopHat200 = M200[i];
393  Group[i].R_TopHat200 = R200[i];
394  break;
395  case 2:
396  Group[i].M_Crit200 = M200[i];
397  Group[i].R_Crit200 = R200[i];
398  break;
399  case 3:
400  Group[i].M_Crit500 = M200[i];
401  Group[i].R_Crit500 = R200[i];
402  break;
403  }
404  }
405  }
406 
407  Mem.myfree(M200);
408  Mem.myfree(R200);
409  Mem.myfree(Right);
410  Mem.myfree(Left);
411  Mem.myfree(TargetList);
412 
413  double tend = Logs.second();
414  return Logs.timediff(tstart, tend);
415 }
416 
417 /* now make sure that the following classes are really instantiated, otherwise we may get a linking problem */
418 #include "../data/simparticles.h"
419 template class fof<simparticles>;
420 
421 #if defined(LIGHTCONE) && defined(LIGHTCONE_PARTICLES_GROUPS)
422 #include "../data/lcparticles.h"
423 template class fof<lcparticles>;
424 #endif
425 
426 #endif
global_data_all_processes All
Definition: main.cc:40
Definition: domain.h:31
virtual void out2particle(T_out *out, int i, int mode)=0
virtual void particle2in(T_in *in, int i)=0
virtual int evaluate(int target, int mode, int thread_id, int action, T_in *in, int numnodes, node_info *firstnode, T_out &out)=0
double timediff(double t0, double t1)
Definition: logs.cc:488
double second(void)
Definition: logs.cc:471
T da[3]
Definition: symtensors.h:106
#define FACTSQRT3
Definition: constants.h:437
#define MAXITER
Definition: constants.h:305
#define MODE_DEFAULT
Definition: constants.h:23
#define MODE_IMPORTED_PARTICLES
Definition: constants.h:22
#define FACT1
Definition: constants.h:435
#define MODE_LOCAL_PARTICLES
Definition: constants.h:21
#define M_PI
Definition: constants.h:56
uint32_t MyIntPosType
Definition: dtypes.h:35
#define BITS_FOR_POSITIONS
Definition: dtypes.h:37
logs Logs
Definition: main.cc:43
#define Terminate(...)
Definition: macros.h:19
void sumup_large_ints(int n, int *src, long long *res, MPI_Comm comm)
memory Mem
Definition: main.cc:44
half fabs(half arg)
Definition: half.hpp:2627
expr pow(half base, half exp)
Definition: half.hpp:2803
unsigned char level
Definition: tree.h:64
unsigned char nextnode_shmrank
Definition: tree.h:66
int nextnode
Definition: tree.h:57
vector< MyIntPosType > center
Definition: tree.h:51
unsigned char sibling_shmrank
Definition: tree.h:65
int sibling
Definition: tree.h:53
MyDouble mass
Definition: gravtree.h:65
Definition: tree.h:77
int Node
Definition: tree.h:78
void myflush(FILE *fstream)
Definition: system.cc:125