GADGET-4
fof_findgroups.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 FOF
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 "../logs/timer.h"
29 #include "../main/simulation.h"
30 #include "../mpi_utils/generic_comm.h"
31 #include "../mpi_utils/mpi_utils.h"
32 #include "../ngbtree/ngbtree.h"
33 #include "../sort/cxxsort.h"
34 #include "../system/system.h"
35 
36 /* local data structure for collecting particle/cell data that is sent to other processors if needed */
37 struct foffind_in : data_in_generic
38 {
39 #if defined(LIGHTCONE_PARTICLES_GROUPS)
40  double DistanceOrigin;
41 #endif
42  MyIntPosType IntPos[3];
43  MyIDStorage MinID;
44  int MinIDTask;
45 };
46 
47 /* local data structure that holds results acquired on remote processors */
48 struct foffind_out
49 {
50  char link_count_flag;
51 };
52 
53 /* routine that fills the relevant particle/cell data into the input structure defined above */
54 
55 template <typename T_tree, typename T_domain, typename T_partset>
56 class foffind_comm : public generic_comm<foffind_in, foffind_out, T_tree, T_domain, T_partset>
57 {
58  public:
60  using gcomm::D;
61  using gcomm::Thread;
62  using gcomm::Tp; // This makes sure that we can access Tp from the base class without having to use "this->Tp"
63  using gcomm::Tree;
64 
65  /* need to call the base class constructor explicitly */
66  foffind_comm(T_domain *dptr, T_tree *tptr, T_partset *pptr) : gcomm(dptr, tptr, pptr) {}
67 
68  void particle2in(foffind_in *in, int i) override
69  {
70  for(int k = 0; k < 3; k++)
71  in->IntPos[k] = Tp->P[i].IntPos[k];
72 
73  in->MinID = Tp->MinID[Tp->Head[i]];
74  in->MinIDTask = Tp->MinIDTask[Tp->Head[i]];
75 
76 #if defined(LIGHTCONE_PARTICLES_GROUPS)
77  in->DistanceOrigin = Tp->DistanceOrigin[Tp->Head[i]];
78 #endif
79  }
80 
81  void out2particle(foffind_out *out, int i, int mode) override
82  {
83  if(mode == MODE_LOCAL_PARTICLES) /* initial store */
84  {
85  /* nothing to be done here */
86  }
87  else /* combine */
88  {
89  if(out->link_count_flag)
90  Tp->Flags[i].Marked = 1;
91  }
92  }
93 
94  int evaluate(int target, int mode, int thread_id, int action, foffind_in *in, int numnodes, node_info *firstnode,
95  foffind_out &out) override
96  {
97  memset(&out, 0, sizeof(out));
98 
99 #if defined(LIGHTCONE_PARTICLES_GROUPS)
100  double target_DistanceOrigin = in->DistanceOrigin;
101 #else
102  double target_DistanceOrigin = 0;
103 #endif
104 
105  int numngb = Tree->treefind_fof_primary(in->IntPos, Tp->LinkL, target, mode, &Thread, numnodes, firstnode, Thread.Ngblist,
106  in->MinID, in->MinIDTask, target_DistanceOrigin);
107 
108  if(mode == MODE_IMPORTED_PARTICLES)
109  {
110  if(numngb > 0)
111  out.link_count_flag = 1;
112  else
113  out.link_count_flag = 0;
114  }
115 
116  return numngb;
117  }
118 };
119 
120 template <typename partset>
121 double fof<partset>::fof_find_groups(void)
122 {
123  double tstart = Logs.second();
124 
125  mpi_printf("FOF: Start linking particles (presently allocated=%g MB)\n", Mem.getAllocatedBytesInMB());
126 
127  Tp->Flags = (typename partset::bit_flags *)Mem.mymalloc_clear("Flags", Tp->NumPart * sizeof(typename partset::bit_flags));
128 
129  FoFNgbTree.FullyLinkedNodePIndex = (int *)Mem.mymalloc("FullyLinkedNodePIndex", FoFNgbTree.NumNodes * sizeof(int));
130  FoFNgbTree.FullyLinkedNodePIndex -= FoFNgbTree.MaxPart;
131 
132  for(int i = 0; i < FoFNgbTree.NumNodes; i++)
133  {
134  int no = i + FoFNgbTree.MaxPart;
135  FoFNgbTree.FullyLinkedNodePIndex[no] = -1;
136  }
137 
138  /* let's link all those primary particles which are in small enough nodes, only process local non-toplevel nodes */
139  double taa = Logs.second();
140  for(int i = 0; i < FoFNgbTree.NumNodes; i++)
141  {
142  int no = i + FoFNgbTree.MaxPart;
143 
144  if(FoFNgbTree.FullyLinkedNodePIndex[no] < 0)
145  {
146  if(FoFNgbTree.get_nodep(no)->level > 0)
147  {
148  double len = (((MyIntPosType)1) << (BITS_FOR_POSITIONS - FoFNgbTree.get_nodep(no)->level)) * Tp->FacIntToCoord;
149  if(FACTSQRT3 * len < Tp->LinkL) // all particles in this node can be linked
150  {
151  int q = FoFNgbTree.treefind_fof_return_a_particle_in_cell_recursive(no);
152 
153  if(q >= 0)
154  FoFNgbTree.fof_link_particles_in_cell_recursive(no, q);
155  }
156  }
157  }
158  }
159  double tbb = Logs.second();
160  mpi_printf("FOF: linking of small cells took %g sec\n", Logs.timediff(taa, tbb));
161 
162  /* first, link only among local particles */
163  int *targetlist = (int *)Mem.mymalloc("targetlist", Tp->NumPart * sizeof(int));
164 
165  int npart = 0;
166  for(int i = 0; i < Tp->NumPart; i++)
167  {
168  if(is_type_primary_link_type(Tp->P[i].getType()))
169  targetlist[npart++] = i;
170  }
171 
172  /* create an object for handling the communication */
173  foffind_comm<foftree<partset>, domain<partset>, partset> commpattern(FoFDomain, &FoFNgbTree, Tp);
174 
175  TIMER_STORE; /* active timer should be CPU_FOF */
176 
177  double t0 = Logs.second();
178 
179  commpattern.execute(npart, targetlist, MODE_LOCAL_NO_EXPORT, logs::CPU_FOFWALK, logs::CPU_FOFWALK, logs::CPU_FOFIMBAL);
180 
181  double t1 = Logs.second();
182 
183  int marked = 0;
184  for(int i = 0; i < Tp->NumPart; i++)
185  {
186  if(is_type_primary_link_type(Tp->P[i].getType()))
187  {
188  if(Tp->Flags[i].Nonlocal)
189  targetlist[marked++] = i;
190  }
191  }
192 
193  double dt = TIMER_DIFF(CPU_FOFWALK), dtmax, dtsum;
194  MPI_Allreduce(&dt, &dtmax, 1, MPI_DOUBLE, MPI_MAX, Communicator);
195  MPI_Allreduce(&dt, &dtsum, 1, MPI_DOUBLE, MPI_SUM, Communicator);
196 
197  long long totmarked, totnpart;
198  sumup_large_ints(1, &marked, &totmarked, Communicator);
199  sumup_large_ints(1, &npart, &totnpart, Communicator);
200  mpi_printf(
201  "FOF: local links done (took %g sec, avg-work=%g, imbalance=%g).\nFOF: Marked=%lld out of the %lld primaries "
202  "which "
203  "are linked\n",
204  Logs.timediff(t0, t1), dtsum / NTask, dtmax / (dtsum / NTask), totmarked, totnpart);
205 
206  npart = marked;
207 
208  mpi_printf("FOF: begin linking across processors (presently allocated=%g MB) \n", Mem.getAllocatedBytesInMB());
209 
210  for(int i = 0; i < Tp->NumPart; i++)
211  Tp->Flags[i].Marked = 1;
212 
213  long long link_across_tot;
214 
215  do
216  {
217  double t0 = Logs.second();
218 
219  npart = 0;
220  for(int i = 0; i < Tp->NumPart; i++)
221  {
222  if(is_type_primary_link_type(Tp->P[i].getType()))
223  {
224  if(Tp->Flags[i].Nonlocal && Tp->Flags[i].Marked)
225  targetlist[npart++] = i;
226  }
227 
228  Tp->Flags[i].MinIDChanged = 0;
229  Tp->Flags[i].Marked = 0;
230  }
231 
232  commpattern.execute(npart, targetlist, MODE_DEFAULT, logs::CPU_FOFWALK, logs::CPU_FOFWALK, logs::CPU_FOFIMBAL);
233 
234  int link_across = commpattern.Thread.Interactions;
235 
236  sumup_large_ints(1, &link_across, &link_across_tot, Communicator);
237 
238  long long ntot;
239  sumup_large_ints(1, &npart, &ntot, Communicator);
240 
241  double t1 = Logs.second();
242 
243  mpi_printf("FOF: have done %15lld cross links (processed %14lld, took %g sec)\n", link_across_tot, ntot, Logs.timediff(t0, t1));
244 
245  /* let's check out which particles have changed their MinID */
246 
247  for(int i = 0; i < Tp->NumPart; i++)
248  {
249  if(Tp->Flags[i].Nonlocal)
250  {
251  if(Tp->Flags[Tp->Head[i]].MinIDChanged)
252  Tp->Flags[i].Marked = 1;
253  }
254  }
255  }
256  while(link_across_tot > 0);
257 
258  Mem.myfree(targetlist);
259  Mem.myfree(FoFNgbTree.FullyLinkedNodePIndex + FoFNgbTree.MaxPart);
260  Mem.myfree(Tp->Flags);
261 
262  mpi_printf("FOF: Local groups found.\n");
263 
264  double tend = Logs.second();
265  return Logs.timediff(tstart, tend);
266 }
267 
273 template <typename partset>
274 int foftree<partset>::treefind_fof_primary(MyIntPosType *searchcenter, MyNgbTreeFloat hsml, int target, int mode, thread_data *thread,
275  int numnodes, node_info *firstnode, int *ngblist, MyIDStorage target_MinID,
276  int target_MinIDTask, double target_DistanceOrigin)
277 {
278  if(mode == MODE_LOCAL_NO_EXPORT)
279  Tp->Flags[target].Nonlocal = 0;
280 
281  MyNgbTreeFloat hsml2 = hsml * hsml;
282 
283  MyIntPosType search_min[3], search_range[3];
284  MyIntPosType inthsml = hsml * Tp->FacCoordToInt;
285 
286  for(int i = 0; i < 3; i++)
287  {
288  search_min[i] = searchcenter[i] - inthsml;
289  search_range[i] = inthsml + inthsml;
290  }
291 
292  int numngb = 0;
293 
294  for(int k = 0; k < numnodes; k++)
295  {
296  int no;
297 
298  if(mode == MODE_LOCAL_PARTICLES || mode == MODE_LOCAL_NO_EXPORT)
299  {
300  no = MaxPart; /* root node */
301  }
302  else
303  {
304  no = firstnode[k].Node;
305  no = get_nodep(no)->nextnode; /* open it */
306  }
307 
308  int shmrank = TreeSharedMem_ThisTask;
309 
310  while(no >= 0)
311  {
312  if(no < MaxPart) /* single particle */
313  {
314  if(shmrank != TreeSharedMem_ThisTask)
315  Terminate("unexpected because in the present algorithm we are only allowed walk local branches");
316 
317  int p = no;
318  auto P = get_Pp(no, shmrank);
319 
320  no = get_nextnodep(shmrank)[no]; /* note: here shmrank cannot change */
321 
322  if(mode == MODE_LOCAL_PARTICLES) /* because we have already linked those in previous phase with MODE_LOCAL_NO_EXPORT */
323  continue;
324 
325  double dx = ((MySignedIntPosType)(P->IntPos[0] - searchcenter[0])) * Tp->FacIntToCoord;
326  double dd = dx * dx;
327  if(dd > hsml2)
328  continue;
329 
330  double dy = ((MySignedIntPosType)(P->IntPos[1] - searchcenter[1])) * Tp->FacIntToCoord;
331  dd += dy * dy;
332  if(dd > hsml2)
333  continue;
334 
335  double dz = ((MySignedIntPosType)(P->IntPos[2] - searchcenter[2])) * Tp->FacIntToCoord;
336  dd += dz * dz;
337  if(dd > hsml2)
338  continue;
339 
340  if(mode == MODE_IMPORTED_PARTICLES)
341  {
342 #if defined(LIGHTCONE_PARTICLES_GROUPS)
343  if(Tp->DistanceOrigin[Tp->Head[p]] > target_DistanceOrigin)
344  {
345  Tp->DistanceOrigin[Tp->Head[p]] = target_DistanceOrigin;
346  Tp->Flags[Tp->Head[p]].MinIDChanged = 1;
347  numngb++;
348  }
349 #endif
350  if(Tp->MinID[Tp->Head[p]].get() > target_MinID.get())
351  {
352  Tp->MinID[Tp->Head[p]] = target_MinID;
353  Tp->MinIDTask[Tp->Head[p]] = target_MinIDTask;
354  Tp->Flags[Tp->Head[p]].MinIDChanged = 1;
355  numngb++;
356  }
357  }
358  else if(mode == MODE_LOCAL_NO_EXPORT)
359  {
360  if(Tp->Head[target] != Tp->Head[p])
361  {
362  int no_target = Father[target];
363  int no_p = Father[p];
364 
365  Tp->link_two_particles(target, p);
366 
367  if(no_target >= 0)
368  {
369  if(FullyLinkedNodePIndex[no_target] >= 0)
370  {
371  if(Tp->Head[FullyLinkedNodePIndex[no_target]] != Tp->Head[target])
372  Terminate("how come that the node is fully linked, but not to us");
373  }
374  else
375  {
376  // this parent node was not yet fully linked to the final group: Is this the case now?
377  // If needed, check also all parent nodes
378  while(treefind_fof_check_single_node_for_full_linking(no_target))
379  no_target = get_nodep(no_target)->father;
380  }
381  }
382 
383  if(no_p >= 0)
384  {
385  if(FullyLinkedNodePIndex[no_p] >= 0)
386  {
387  if(Tp->Head[FullyLinkedNodePIndex[no_p]] != Tp->Head[p])
388  Terminate("how come that the node is fully linked, but not to us");
389  }
390  else
391  {
392  // this parent node was not yet fully linked to the final group: Is this the case now?
393  // If needed, check also all parent nodes
394  while(treefind_fof_check_single_node_for_full_linking(no_p))
395  no_p = get_nodep(no_p)->father;
396  }
397  }
398  }
399  }
400  else
401  Terminate("strange mode");
402  }
403  else if(no < MaxPart + MaxNodes) /* internal node */
404  {
405  fofnode *current = get_nodep(no, shmrank);
406 
407  if(current->level == 0)
408  {
409  /* we always open the root node (its full node length couldn't be stored in the integer type */
410  no = current->nextnode; /* no change in shmrank expected here */
411  shmrank = current->nextnode_shmrank;
412  continue;
413  }
414  /* check whether the node lies outside our search range */
415 
416  if(mode == MODE_IMPORTED_PARTICLES)
417  {
418  if(no < FirstNonTopLevelNode) /* we reached a top-level node again, which means that we are done with the branch */
419  break;
420 
421  if(FullyLinkedNodePIndex[no] >= 0)
422  {
423  int head = Tp->Head[FullyLinkedNodePIndex[no]];
424 
425  if(head >= 0)
426  if(Tp->MinID[head].get() <= target_MinID.get())
427  {
428 #if defined(LIGHTCONE_PARTICLES_GROUPS)
429  if(Tp->DistanceOrigin[FullyLinkedNodePIndex[no]] <= target_DistanceOrigin)
430 #endif
431  {
432  no = current->sibling; /* the node can be discarded */
433  shmrank = current->sibling_shmrank;
434  continue;
435  }
436  }
437  }
438  }
439  else if(mode == MODE_LOCAL_PARTICLES)
440  {
441  int p = current->nextnode;
442 
443  /* in case the next node after opening is not a top-level node, we have either reached a leaf node or are in a local
444  * branch we need to do nothing if we would end up on different shared memory thread */
445  if(p < MaxPart || (p >= FirstNonTopLevelNode && p < MaxPart + MaxNodes))
446  {
447  if(current->nextnode_shmrank != TreeSharedMem_ThisTask)
448  {
449  int task = D->ThisTask + current->nextnode_shmrank - TreeSharedMem_ThisTask;
450 
451  if(target >= 0) /* export */
452  tree_export_node_threads_by_task_and_node(task, no, target, thread);
453 
454  no = current->sibling; /* in case the node can be discarded */
455  shmrank = current->sibling_shmrank;
456  continue;
457  }
458  }
459 
460  if(no >= FirstNonTopLevelNode)
461  {
462  /* we have a node with only local particles, hence we can skip it for mode == 0 */
463  no = current->sibling; /* in case the node can be discarded */
464  shmrank = current->sibling_shmrank;
465  continue;
466  }
467  }
468  else if(mode == MODE_LOCAL_NO_EXPORT)
469  {
470  int p = current->nextnode;
471 
472  /* in case the next node after opening is not a top-level node, we have either reached a leaf node or are in a local
473  * branch we need to do nothing if we would end up on different shared memory thread */
474  if(p < MaxPart || (p >= FirstNonTopLevelNode && p < MaxPart + MaxNodes))
475  {
476  if(current->nextnode_shmrank != TreeSharedMem_ThisTask)
477  {
478  no = current->sibling; /* in case the node can be discarded */
479  shmrank = current->sibling_shmrank;
480 
481  MyIntPosType left[3], right[3];
482 
483  left[0] = current->range_min[0] - search_min[0];
484  right[0] = current->range_max[0] - search_min[0];
485 
486  /* check whether we can stop walking along this branch */
487  if(left[0] > search_range[0] && right[0] > left[0])
488  continue;
489 
490  left[1] = current->range_min[1] - search_min[1];
491  right[1] = current->range_max[1] - search_min[1];
492 
493  /* check whether we can stop walking along this branch */
494  if(left[1] > search_range[1] && right[1] > left[1])
495  continue;
496 
497  left[2] = current->range_min[2] - search_min[2];
498  right[2] = current->range_max[2] - search_min[2];
499 
500  /* check whether we can stop walking along this branch */
501  if(left[2] > search_range[2] && right[2] > left[2])
502  continue;
503 
504  Tp->Flags[target].Nonlocal = 1;
505  continue;
506  }
507  }
508 
509  if(FullyLinkedNodePIndex[no] >= 0)
510  if(Tp->Head[target] == Tp->Head[FullyLinkedNodePIndex[no]]) // all particles in the node are linked to us anyhow
511  {
512  no = current->sibling; /* in case the node can be discarded */
513  shmrank = current->sibling_shmrank;
514  continue;
515  }
516  }
517 
518  MyIntPosType left[3], right[3];
519 
520  left[0] = current->range_min[0] - search_min[0];
521  right[0] = current->range_max[0] - search_min[0];
522 
523  /* check whether we can stop walking along this branch */
524  if(left[0] > search_range[0] && right[0] > left[0])
525  {
526  no = current->sibling; /* in case the node can be discarded */
527  shmrank = current->sibling_shmrank;
528  continue;
529  }
530 
531  left[1] = current->range_min[1] - search_min[1];
532  right[1] = current->range_max[1] - search_min[1];
533 
534  /* check whether we can stop walking along this branch */
535  if(left[1] > search_range[1] && right[1] > left[1])
536  {
537  no = current->sibling; /* in case the node can be discarded */
538  shmrank = current->sibling_shmrank;
539  continue;
540  }
541 
542  left[2] = current->range_min[2] - search_min[2];
543  right[2] = current->range_max[2] - search_min[2];
544 
545  /* check whether we can stop walking along this branch */
546  if(left[2] > search_range[2] && right[2] > left[2])
547  {
548  no = current->sibling; /* in case the node can be discarded */
549  shmrank = current->sibling_shmrank;
550  continue;
551  }
552 
553  no = current->nextnode; /* ok, we need to open the node */
554  shmrank = current->nextnode_shmrank; /* ok, we need to open the node */
555  }
556  else /* pseudo particle */
557  {
558  if(mode == MODE_LOCAL_PARTICLES)
559  {
560  if(target >= 0) /* if no target is given, export will not occur */
561  tree_export_node_threads(no, target, thread);
562  }
563  else if(mode == MODE_LOCAL_NO_EXPORT)
564  {
565  Tp->Flags[target].Nonlocal = 1;
566  }
567 
568  no = get_nextnodep(shmrank)[no - MaxNodes];
569  /* note: here shmrank does not need to change */
570 
571  continue;
572  }
573  }
574  }
575 
576  return numngb;
577 }
578 
579 template <typename partset>
581 {
582  if(no >= MaxPart && no < MaxPart + MaxNodes) /* internal node */
583  {
584  FullyLinkedNodePIndex[no] = q;
585 
586  int p = get_nodep(no)->nextnode;
587 
588  /* in case the next node after opening is not a top-level node, we have either reached a leaf node or are in a local
589  * branch. We need to do nothing if we would end up on different shared memory thread */
590  if(p < MaxPart || (p >= FirstNonTopLevelNode && p < MaxPart + MaxNodes))
591  {
592  if(get_nodep(no)->nextnode_shmrank != TreeSharedMem_ThisTask)
593  return;
594  }
595 
596  while(p != get_nodep(no)->sibling)
597  {
598  if(p < MaxPart) /* a particle */
599  {
600  if(p != q)
601  {
602  Tp->link_two_particles(p, q); // link them if not already linked
603  }
604 
605  p = Nextnode[p];
606  }
607  else if(p < MaxPart + MaxNodes) /* an internal node */
608  {
609  fof_link_particles_in_cell_recursive(p, q);
610 
611  p = get_nodep(p)->sibling;
612  }
613  else /* a pseudo particle */
614  p = Nextnode[p - MaxNodes];
615  }
616  }
617 }
618 
619 template <typename partset>
621 {
622  if(no >= MaxPart && no < MaxPart + MaxNodes) /* internal node */
623  {
624  int p = get_nodep(no)->nextnode;
625 
626  /* in case the next node after opening is not a top-level node, we have either reached a leaf node or are in a local
627  * branch. We need to do nothing if we would end up on different shared memory thread */
628  if(p < MaxPart || (p >= FirstNonTopLevelNode && p < MaxPart + MaxNodes))
629  {
630  if(get_nodep(no)->nextnode_shmrank != TreeSharedMem_ThisTask)
631  return -1;
632  }
633 
634  while(p != get_nodep(no)->sibling)
635  {
636  if(p < MaxPart) /* a particle */
637  {
638  return p;
639 
640  p = Nextnode[p];
641  }
642  else if(p < MaxPart + MaxNodes) /* an internal node */
643  {
644  int ret = treefind_fof_return_a_particle_in_cell_recursive(p);
645 
646  if(ret >= 0)
647  return ret;
648 
649  p = get_nodep(p)->sibling;
650  }
651  else /* a pseudo particle */
652  p = Nextnode[p - MaxNodes];
653  }
654  }
655 
656  return -1;
657 }
658 
659 template <typename partset>
661 {
662  if(no >= MaxPart && no < MaxPart + MaxNodes) /* internal node */
663  {
664  if(FullyLinkedNodePIndex[no] >= 0) // already linked
665  return 0;
666 
667  int head = -1; /* no particle yet */
668 
669  int p = get_nodep(no)->nextnode;
670 
671  /* in case the next node after opening is not a top-level node, we have either reached a leaf node or are in a local
672  * branch. We need to do nothing if we would end up on different shared memory thread */
673  if(p < MaxPart || (p >= FirstNonTopLevelNode && p < MaxPart + MaxNodes))
674  {
675  if(get_nodep(no)->nextnode_shmrank != TreeSharedMem_ThisTask)
676  return 0;
677  }
678 
679  while(p != get_nodep(no)->sibling)
680  {
681  if(p < MaxPart) /* a particle */
682  {
683  if(head == -1)
684  head = Tp->Head[p];
685  else if(head >= 0)
686  {
687  if(head != Tp->Head[p])
688  {
689  head = -2;
690  break;
691  }
692  }
693 
694  p = Nextnode[p];
695  }
696  else if(p < MaxPart + MaxNodes) /* an internal node */
697  {
698  if(FullyLinkedNodePIndex[p] >= 0)
699  {
700  if(head == -1)
701  head = Tp->Head[FullyLinkedNodePIndex[p]];
702  else if(head >= 0)
703  {
704  if(head != Tp->Head[FullyLinkedNodePIndex[p]])
705  {
706  head = -2;
707  break;
708  }
709  }
710  }
711  else
712  {
713  if(treefind_fof_return_a_particle_in_cell_recursive(no) >= 0)
714  {
715  head = -2;
716  break;
717  }
718  }
719 
720  p = get_nodep(p)->sibling;
721  }
722  else /* a pseudo particle */
723  p = Nextnode[p - MaxNodes];
724  }
725 
726  if(head >= 0)
727  {
728  FullyLinkedNodePIndex[no] = treefind_fof_return_a_particle_in_cell_recursive(no);
729 
730  if(Tp->Head[FullyLinkedNodePIndex[no]] != head)
731  Terminate("no=%d FullyLinkedNodePIndex[no]=%d Tp->Head[FullyLinkedNodePIndex[no]]=%d head=%d \n", no,
732  FullyLinkedNodePIndex[no], Tp->Head[FullyLinkedNodePIndex[no]], head);
733 
734  return 1;
735  }
736  }
737 
738  return 0;
739 }
740 
741 /* now make sure that the following classes are really instantiated, otherwise we may get a linking problem */
742 #include "../data/simparticles.h"
743 template class fof<simparticles>;
744 
745 #if defined(LIGHTCONE) && defined(LIGHTCONE_PARTICLES_GROUPS)
746 #include "../data/lcparticles.h"
747 template class fof<lcparticles>;
748 #endif
749 
750 #endif
MyIDType get(void) const
Definition: idstorage.h:87
Definition: domain.h:31
void fof_link_particles_in_cell_recursive(int no, int q)
int treefind_fof_return_a_particle_in_cell_recursive(int no)
int treefind_fof_check_single_node_for_full_linking(int no)
int treefind_fof_primary(MyIntPosType *searchcenter, MyNgbTreeFloat hsml, int target, int mode, thread_data *thread, int numnodes, node_info *firstnode, int *ngblist, MyIDStorage target_MinID, int target_MinIDTask, double target_DistanceOrigin)
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
double getAllocatedBytesInMB(void)
Definition: mymalloc.h:68
#define FACTSQRT3
Definition: constants.h:437
#define MODE_DEFAULT
Definition: constants.h:23
#define MODE_IMPORTED_PARTICLES
Definition: constants.h:22
#define MODE_LOCAL_NO_EXPORT
Definition: constants.h:24
#define MODE_LOCAL_PARTICLES
Definition: constants.h:21
uint32_t MyIntPosType
Definition: dtypes.h:35
float MyNgbTreeFloat
Definition: dtypes.h:88
int32_t MySignedIntPosType
Definition: dtypes.h:36
#define BITS_FOR_POSITIONS
Definition: dtypes.h:37
logs Logs
Definition: main.cc:43
#define TIMER_STORE
Copies the current value of CPU times to a stored variable, such that differences with respect to thi...
Definition: logs.h:257
#define TIMER_DIFF(counter)
Returns amount elapsed for the timer since last save with TIMER_STORE.
Definition: logs.h:250
#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
unsigned char level
Definition: tree.h:64
unsigned char nextnode_shmrank
Definition: tree.h:66
int nextnode
Definition: tree.h:57
unsigned char sibling_shmrank
Definition: tree.h:65
int sibling
Definition: tree.h:53
MyIntPosType range_max[3]
Definition: foftree.h:33
MyIntPosType range_min[3]
Definition: foftree.h:32
Definition: tree.h:77
int Node
Definition: tree.h:78