GADGET-4
subfind_orphanids.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_ORPHAN_TREATMENT
15 
16 #include <errno.h>
17 #include <hdf5.h>
18 #include <math.h>
19 #include <mpi.h>
20 #include <stdio.h>
21 #include <stdlib.h>
22 #include <string.h>
23 #include <sys/stat.h>
24 
25 #include "../cooling_sfr/cooling.h"
26 #include "../data/allvars.h"
27 #include "../data/dtypes.h"
28 #include "../data/mymalloc.h"
29 #include "../fof/fof.h"
30 #include "../io/hdf5_util.h"
31 #include "../io/io.h"
32 #include "../logs/timer.h"
33 #include "../main/main.h"
34 #include "../main/simulation.h"
35 #include "../mergertree/io_readsnap.h"
36 #include "../mergertree/mergertree.h"
37 #include "../mpi_utils/mpi_utils.h"
38 #include "../sort/parallel_sort.h"
39 #include "../system/system.h"
40 
41 template <>
42 void fof<simparticles>::subfind_match_ids_of_previously_most_bound_ids(simparticles *Sp)
43 {
44  int *Send_count = (int *)Mem.mymalloc("Send_count", sizeof(int) * NTask);
45  int *Send_offset = (int *)Mem.mymalloc("Send_offset", sizeof(int) * NTask);
46  int *Recv_count = (int *)Mem.mymalloc("Recv_count", sizeof(int) * NTask);
47  int *Recv_offset = (int *)Mem.mymalloc("Recv_offset", sizeof(int) * NTask);
48 
49  mycxxsort_parallel(Sp->IdStore.ID, Sp->IdStore.ID + Sp->IdStore.NumPart, Sp->compare_IDs, Communicator);
50  mycxxsort_parallel(Sp->P, Sp->P + Sp->NumPart, Sp->compare_SpP_ID, Communicator);
51 
52  MyIDType *list_min_id = (MyIDType *)Mem.mymalloc("list_min_id", NTask * sizeof(MyIDType));
53  MyIDType *list_max_id = (MyIDType *)Mem.mymalloc("list_max_id", NTask * sizeof(MyIDType));
54 
55  MyIDType idmin = Sp->P[0].ID.get();
56  MyIDType idmax = Sp->P[Sp->NumPart - 1].ID.get();
57 
58  MPI_Allgather(&idmin, sizeof(MyIDType), MPI_BYTE, list_min_id, sizeof(MyIDType), MPI_BYTE, Communicator);
59  MPI_Allgather(&idmax, sizeof(MyIDType), MPI_BYTE, list_max_id, sizeof(MyIDType), MPI_BYTE, Communicator);
60 
61  int *num_list = (int *)Mem.mymalloc("num_list", NTask * sizeof(int));
62  MPI_Allgather(&Sp->NumPart, 1, MPI_INT, num_list, 1, MPI_INT, Communicator);
63 
64  int nexport = 0, nimport = 0;
65  MyIDType *import_data = NULL, *export_data = NULL;
66 
67  for(int mode = 0; mode < 2; mode++)
68  {
69  for(int i = 0; i < NTask; i++)
70  Send_count[i] = 0;
71 
72  int target = 0;
73 
74  for(int i = 0; i < Sp->IdStore.NumPart; i++)
75  {
76  while(target < NTask - 1 && (num_list[target] == 0 || Sp->IdStore.ID[i] > list_max_id[target]))
77  target++;
78 
79  if(num_list[target] == 0)
80  Terminate("How can this be? target=%d", target);
81 
82  if(Sp->IdStore.ID[i] >= list_min_id[target] && Sp->IdStore.ID[i] <= list_max_id[target])
83  {
84  if(mode == 0)
85  Send_count[target]++;
86  else
87  {
88  int off = Send_offset[target] + Send_count[target]++;
89  export_data[off] = Sp->IdStore.ID[i];
90  }
91  }
92  }
93 
94  if(mode == 0)
95  {
96  MPI_Alltoall(Send_count, 1, MPI_INT, Recv_count, 1, MPI_INT, Communicator);
97  Recv_offset[0] = Send_offset[0] = 0;
98  for(int j = 0; j < NTask; j++)
99  {
100  nimport += Recv_count[j];
101  nexport += Send_count[j];
102  if(j > 0)
103  {
104  Send_offset[j] = Send_offset[j - 1] + Send_count[j - 1];
105  Recv_offset[j] = Recv_offset[j - 1] + Recv_count[j - 1];
106  }
107  }
108 
109  export_data = (MyIDType *)Mem.mymalloc("export_data", nexport * sizeof(MyIDType));
110  import_data = (MyIDType *)Mem.mymalloc("import_data", nimport * sizeof(MyIDType));
111  }
112  }
113 
114  for(int ngrp = 0; ngrp < (1 << PTask); ngrp++) /* note: here we also have a transfer from each task to itself (for ngrp=0) */
115  {
116  int recvTask = ThisTask ^ ngrp;
117  if(recvTask < NTask)
118  if(Send_count[recvTask] > 0 || Recv_count[recvTask] > 0)
119  MPI_Sendrecv(&export_data[Send_offset[recvTask]], Send_count[recvTask] * sizeof(MyIDType), MPI_BYTE, recvTask, TAG_DENS_B,
120  &import_data[Recv_offset[recvTask]], Recv_count[recvTask] * sizeof(MyIDType), MPI_BYTE, recvTask, TAG_DENS_B,
121  Communicator, MPI_STATUS_IGNORE);
122  }
123 
124  /* incoming data should already be sorted, so now do the match */
125 
126  int nmarked = 0;
127  for(int i = 0, j = 0; i < Sp->NumPart && j < nimport;)
128  {
129  if(Sp->P[i].ID.get() < import_data[j])
130  i++;
131  else if(Sp->P[i].ID.get() > import_data[j])
132  j++;
133  else
134  {
135  if(!Sp->P[i].ID.is_previously_most_bound())
136  {
138  nmarked++;
139  }
140  i++;
141  j++;
142  }
143  }
144 
145  Mem.myfree(import_data);
146  Mem.myfree(export_data);
147 
148  Mem.myfree(num_list);
149  Mem.myfree(list_max_id);
150  Mem.myfree(list_min_id);
151 
152  Mem.myfree(Recv_offset);
153  Mem.myfree(Recv_count);
154  Mem.myfree(Send_offset);
155  Mem.myfree(Send_count);
156 
157  long long tot_ncheck, tot_nmarked;
158  sumup_large_ints(1, &Sp->IdStore.NumPart, &tot_ncheck, Communicator);
159  sumup_large_ints(1, &nmarked, &tot_nmarked, Communicator);
160 
161  mpi_printf("SUBFIND_ORPHAN_TREATMENT: Got %lld particles from previous snapshot, led to %lld additionally marked particles\n",
162  tot_ncheck, tot_nmarked);
163 }
164 
165 #endif
void mark_as_formerly_most_bound(void)
Definition: idstorage.h:107
MyIDType get(void) const
Definition: idstorage.h:87
bool is_previously_most_bound(void)
Definition: idstorage.h:117
particle_data * P
Definition: simparticles.h:54
static bool compare_IDs(const MyIDType &a, const MyIDType &b)
Definition: simparticles.h:72
unsigned int MyIDType
Definition: dtypes.h:68
#define Terminate(...)
Definition: macros.h:19
void sumup_large_ints(int n, int *src, long long *res, MPI_Comm comm)
#define TAG_DENS_B
Definition: mpi_utils.h:51
memory Mem
Definition: main.cc:44
double mycxxsort_parallel(T *begin, T *end, Comp comp, MPI_Comm comm)
MyIDStorage ID
Definition: particle_data.h:70