GADGET-4
rearrange.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 MERGERTREE
15
16#include <gsl/gsl_rng.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#include <sys/types.h>
25#include <unistd.h>
26
27#include "../data/allvars.h"
28#include "../data/dtypes.h"
29#include "../data/intposconvert.h"
30#include "../data/mymalloc.h"
31#include "../io/hdf5_util.h"
32#include "../io/snap_io.h"
33#include "../lightcone/lightcone.h"
34#include "../lightcone/lightcone_particle_io.h"
35#include "../main/main.h"
36#include "../main/simulation.h"
37#include "../mergertree/io_readtrees_mbound.h"
38#include "../sort/cxxsort.h"
39#include "../sort/parallel_sort.h"
40#include "../system/system.h"
41#include "../time_integration/driftfac.h"
42
43#if defined(REARRANGE_OPTION) && defined(MERGERTREE)
44void sim::rearrange_snapshot(int argc, char **argv)
45{
46 if(argc < 5)
47 Terminate("too few arguments: Snapshot rearrange requires <firstnum> <lastnum>");
48
49 int conenr = 0; // not needed here
50 int firstnum = atoi(argv[3]);
51 int lastnum = atoi(argv[4]);
52
55
56 rearrange_generic<simparticles>(Sp, conenr, firstnum, lastnum);
57}
58
59/* specialization */
60template <>
61void sim::rearrange_read<simparticles>(simparticles &Tp, int num, int conenr)
62{
63 /* read the lightcone */
64 snap_io Snap{&Tp, Communicator, All.SnapFormat}; /* get an I/O object */
65 Snap.read_snapshot(num, MOST_BOUND_PARTICLE_SNAPHOT);
66}
67
68/* specialization */
69template <>
70void sim::rearrange_write<simparticles>(simparticles &Tp, int num, int conenr)
71{
72 /* write the snapshot */
73 snap_io Snap{&Tp, &MergerTree, Communicator, All.SnapFormat}; /* get an I/O object */
74 Snap.write_snapshot(num, MOST_BOUND_PARTICLE_SNAPHOT_REORDERED); // true?
75}
76#endif
77
78#if defined(LIGHTCONE) && defined(LIGHTCONE_PARTICLES) && defined(REARRANGE_OPTION)
79
80void sim::rearrange_lightcone(int argc, char **argv)
81{
82 if(argc < 6)
83 Terminate("too few arguments: Lightcone rearrange requires <conenr> <firstnum> <lastnum>");
84
85 int conenr = atoi(argv[3]);
86 int firstnum = atoi(argv[4]);
87 int lastnum = atoi(argv[5]);
88
91
92 if(LightCone.lightcone_init_times())
93 endrun();
94
95#ifdef LIGHTCONE_MASSMAPS
96 LightCone.lightcone_init_massmaps();
97 if(LightCone.lightcone_massmap_report_boundaries())
98 endrun();
99#endif
100
101 double linklength = 0;
102 LightCone.lightcone_init_intposconverter(linklength);
103
104 rearrange_generic<lcparticles>(Lp, conenr, firstnum, lastnum);
105}
106
107/* specialization */
108template <>
109void sim::rearrange_read<lcparticles>(lcparticles &Tp, int num, int conenr)
110{
111 /* read the lightcone */
112 lightcone_particle_io LcIO{&Tp, &LightCone, &MergerTree, Communicator, All.SnapFormat}; /* get an I/O object */
113 LcIO.lightcone_read(num, conenr);
114}
115
116/* specialization */
117template <>
118void sim::rearrange_write<lcparticles>(lcparticles &Tp, int num, int conenr)
119{
120 /* write the lightcone */
121 lightcone_particle_io LcIO{&Tp, &LightCone, &MergerTree, Communicator, All.SnapFormat}; /* get an I/O object */
122 LcIO.lightcone_save(num, conenr, true);
123}
124
125#endif
126
127template <typename partset>
128void sim::rearrange_generic(partset &Tp, int conenr, int firstnum, int lastnum)
129{
130 /* read the merger tree mostbound halo IDs */
131 {
132 readtrees_mbound_io MtreeIO{&MergerTree, Communicator, All.SnapFormat}; /* get an I/O object */
133 MtreeIO.read_trees_mostbound(); // optimize this here to only read most-bound ID
134 }
135
136 /* let's now sort the tree data according to the IDs */
137 mycxxsort_parallel(MergerTree.HaloIDdata, MergerTree.HaloIDdata + MergerTree.Nhalos, MergerTree.compare_HaloIDdata_ID, Communicator);
138
139 /* establish some lists with the minimum and maximum IDs stored on each processor */
140 long long halominID = (MergerTree.Nhalos > 0) ? MergerTree.HaloIDdata[0].SubMostBoundID : 0;
141 long long halomaxID = (MergerTree.Nhalos > 0) ? MergerTree.HaloIDdata[MergerTree.Nhalos - 1].SubMostBoundID : 0;
142
143 long long *list_halominID = (long long *)Mem.mymalloc("list_halominID", NTask * sizeof(long long));
144 long long *list_halomaxID = (long long *)Mem.mymalloc("list_halomaxID", NTask * sizeof(long long));
145 int *list_Nhalos = (int *)Mem.mymalloc("list_Nhalos", NTask * sizeof(int));
146
147 MPI_Allgather(&halominID, 1, MPI_LONG_LONG, list_halominID, 1, MPI_LONG_LONG, Communicator);
148 MPI_Allgather(&halomaxID, 1, MPI_LONG_LONG, list_halomaxID, 1, MPI_LONG_LONG, Communicator);
149 MPI_Allgather(&MergerTree.Nhalos, 1, MPI_INT, list_Nhalos, 1, MPI_INT, Communicator);
150
151 typedef mergertree::treehalo_ids_type treehalo_ids_type;
152
153 for(int num = firstnum; num <= lastnum; num++)
154 {
155 rearrange_read(Tp, num, conenr);
156
157 mpi_printf("REARRANGE: On Task=%d, %d particles\n", ThisTask, Tp.NumPart);
158
159 /* let's now sort the lightcone_particle_data according to ID */
160 mycxxsort_parallel(Tp.P, Tp.P + Tp.NumPart, Tp.compare_ID, Communicator);
161
162 /* establish some lists with the minimum and maximum IDs stored on each processor for the lightcone particles */
163 long long coneminID = (Tp.NumPart > 0) ? Tp.P[0].ID.get() : 0;
164 long long conemaxID = (Tp.NumPart > 0) ? Tp.P[Tp.NumPart - 1].ID.get() : 0;
165
166 long long *list_coneminID = (long long *)Mem.mymalloc("list_coneminID", NTask * sizeof(long long));
167 long long *list_conemaxID = (long long *)Mem.mymalloc("list_conemaxID", NTask * sizeof(long long));
168 int *list_NumPart = (int *)Mem.mymalloc("list_NumPart", NTask * sizeof(int));
169
170 MPI_Allgather(&coneminID, 1, MPI_LONG_LONG, list_coneminID, 1, MPI_LONG_LONG, Communicator);
171 MPI_Allgather(&conemaxID, 1, MPI_LONG_LONG, list_conemaxID, 1, MPI_LONG_LONG, Communicator);
172 MPI_Allgather(&Tp.NumPart, 1, MPI_INT, list_NumPart, 1, MPI_INT, Communicator);
173
174 // let's assign matchIDs to the degree possible and assign TreeIDs in accordance to that
175
176 // default is no matching Tree
177 for(int n = 0; n < Tp.NumPart; n++)
178 Tp.P[n].TreeID = -1;
179
180 MPI_Request *requests = (MPI_Request *)Mem.mymalloc("requests", NTask * sizeof(MPI_Request));
181
182 int nreq = 0;
183
184 for(int task = 0; task < NTask; task++)
185 {
186 if(MergerTree.Nhalos > 0 && list_NumPart[task] > 0)
187 if(!(halomaxID < list_coneminID[task] || halominID > list_conemaxID[task]))
188 {
189 MPI_Issend(MergerTree.HaloIDdata, MergerTree.Nhalos * sizeof(mergertree::treehalo_ids_type), MPI_BYTE, task, TAG_N,
190 Communicator, &requests[nreq++]);
191 }
192 }
193
194 for(int task = 0; task < NTask; task++)
195 {
196 if(list_Nhalos[task] > 0 && Tp.NumPart > 0)
197 if(!(list_halomaxID[task] < coneminID || list_halominID[task] > conemaxID))
198 {
199 treehalo_ids_type *halos = (treehalo_ids_type *)Mem.mymalloc("halos", list_Nhalos[task] * sizeof(treehalo_ids_type));
200
201 MPI_Recv(halos, list_Nhalos[task] * sizeof(treehalo_ids_type), MPI_BYTE, task, TAG_N, Communicator, MPI_STATUS_IGNORE);
202
203 int i_halo = 0;
204 int i_cone = 0;
205
206 while(i_halo < list_Nhalos[task] && i_cone < Tp.NumPart)
207 {
208 if(halos[i_halo].SubMostBoundID == Tp.P[i_cone].ID.get())
209 {
210 Tp.P[i_cone++].TreeID = halos[i_halo++].TreeID;
211 }
212 else if(halos[i_halo].SubMostBoundID < Tp.P[i_cone].ID.get())
213 i_halo++;
214 else
215 i_cone++;
216 }
217
218 Mem.myfree(halos);
219 }
220 }
221
222 if(nreq)
223 MPI_Waitall(nreq, requests, MPI_STATUSES_IGNORE);
224
225 Mem.myfree(requests);
226 Mem.myfree(list_NumPart);
227 Mem.myfree(list_conemaxID);
228 Mem.myfree(list_coneminID);
229
230 /* now we sort the lightcone_particle_data according to TreeID */
231 mycxxsort_parallel(Tp.P, Tp.P + Tp.NumPart, Tp.compare_TreeID_ID, Communicator);
232
233 rearrange_fill_treetable<partset>(Tp);
234
235 /* write the lightcone */
236 rearrange_write(Tp, num, conenr);
237
238 /* free the storage again */
239 Tp.free_memory();
240 }
241
242 Mem.myfree(list_Nhalos);
243 Mem.myfree(list_halomaxID);
244 Mem.myfree(list_halominID);
245}
246
247template <typename partset>
248void sim::rearrange_fill_treetable(partset &Tp)
249{
250 /*
251 we use HaloCount just like ParticleCount, and FirstHalo becomes "ParticleFirst */
252
253 int *Send_count = (int *)Mem.mymalloc("Send_count", sizeof(int) * NTask);
254 int *Send_offset = (int *)Mem.mymalloc("Send_offset", sizeof(int) * NTask);
255 int *Recv_count = (int *)Mem.mymalloc("Recv_count", sizeof(int) * NTask);
256
257 for(int i = 0; i < NTask; i++)
258 Send_count[i] = 0;
259
260 long long *TreeID_list = (long long *)Mem.mymalloc("TreeID_list", sizeof(long long) * Tp.NumPart);
261
262 long long maxTreeID = (MergerTree.Ntrees > 0) ? MergerTree.TreeTable[MergerTree.Ntrees - 1].TreeID : -1;
263
264 long long *maxTreeID_list = (long long *)Mem.mymalloc("maxTreeID_list", sizeof(long long) * NTask);
265 MPI_Allgather(&maxTreeID, sizeof(long long), MPI_BYTE, maxTreeID_list, sizeof(long long), MPI_BYTE, Communicator);
266
267 int target_task = 0;
268
269 for(int i = 0; i < Tp.NumPart; i++)
270 {
271 TreeID_list[i] = Tp.P[i].TreeID;
272
273 while(target_task < NTask - 1 && TreeID_list[i] > maxTreeID_list[target_task])
274 target_task++;
275
276 Send_count[target_task]++;
277 }
278
279 myMPI_Alltoall(Send_count, 1, MPI_INT, Recv_count, 1, MPI_INT, Communicator);
280
281 int nexport = 0, nimport = 0;
282 Send_offset[0] = 0;
283
284 for(int j = 0; j < NTask; j++)
285 {
286 nexport += Send_count[j];
287 nimport += Recv_count[j];
288
289 if(j > 0)
290 Send_offset[j] = Send_offset[j - 1] + Send_count[j - 1];
291 }
292
293 for(int i = 0; i < MergerTree.Ntrees; i++)
294 MergerTree.TreeTable[i].HaloCount = 0;
295
296 /* exchange data */
297 for(int ngrp = 0; ngrp < (1 << PTask); ngrp++)
298 {
299 int recvTask = ThisTask ^ ngrp;
300
301 if(recvTask < NTask)
302 {
303 if(Send_count[recvTask] > 0 || Recv_count[recvTask] > 0)
304 {
305 long long *treeid_tmp = (long long *)Mem.mymalloc("treeid_tmp", sizeof(long long) * Recv_count[recvTask]);
306
307 myMPI_Sendrecv(&TreeID_list[Send_offset[recvTask]], Send_count[recvTask] * sizeof(long long), MPI_BYTE, recvTask,
308 TAG_DENS_A, treeid_tmp, Recv_count[recvTask] * sizeof(long long), MPI_BYTE, recvTask, TAG_DENS_A,
309 Communicator, MPI_STATUS_IGNORE);
310
311 for(int i = 0; i < Recv_count[recvTask]; i++)
312 {
313 if(treeid_tmp[i] != -1)
314 {
315 int off = treeid_tmp[i] - MergerTree.TreeTable[0].TreeID;
316
317 if(off < 0 || off >= MergerTree.Ntrees)
318 Terminate(
319 "strange: off=%d MergerTree.Ntrees=%d treeid_tmp[i]=%lld MergerTree.TreeTable[0].TreeID=%lld i=%d "
320 "Recv_count[recvTask]=%d ",
321 off, MergerTree.Ntrees, treeid_tmp[i], MergerTree.TreeTable[0].TreeID, i, Recv_count[recvTask]);
322
323 MergerTree.TreeTable[off].HaloCount++;
324 }
325 }
326
327 Mem.myfree(treeid_tmp);
328 }
329 }
330 }
331
332 Mem.myfree(maxTreeID_list);
333 Mem.myfree(TreeID_list);
334
335 Mem.myfree(Recv_count);
336 Mem.myfree(Send_offset);
337 Mem.myfree(Send_count);
338
339 /* now also fix the cumulative count */
340
341 if(MergerTree.Ntrees > 0)
342 MergerTree.TreeTable[0].FirstHalo = 0;
343
344 for(int i = 1; i < MergerTree.Ntrees; i++)
345 MergerTree.TreeTable[i].FirstHalo = MergerTree.TreeTable[i - 1].FirstHalo + MergerTree.TreeTable[i - 1].HaloCount;
346
347 long long cumul = 0;
348 for(int i = 0; i < MergerTree.Ntrees; i++)
349 cumul += MergerTree.TreeTable[i].HaloCount;
350
351 long long *cumul_list = (long long *)Mem.mymalloc("cumul_list", sizeof(long long) * NTask);
352
353 MPI_Allgather(&cumul, sizeof(long long), MPI_BYTE, cumul_list, sizeof(long long), MPI_BYTE, Communicator);
354
355 cumul = 0;
356 for(int i = 0; i < ThisTask; i++)
357 cumul += cumul_list[i];
358
359 for(int i = 0; i < MergerTree.Ntrees; i++)
360 MergerTree.TreeTable[i].FirstHalo += cumul;
361
362 Mem.myfree(cumul_list);
363}
364
365#endif
global_data_all_processes All
Definition: main.cc:40
void init_drift_table(void)
Definition: driftfac.cc:26
void mpi_printf(const char *fmt,...)
Definition: setcomm.h:55
int ThisTask
Definition: setcomm.h:33
int NTask
Definition: setcomm.h:32
int PTask
Definition: setcomm.h:34
MPI_Comm Communicator
Definition: setcomm.h:31
void rearrange_write(partset &Tp, int num, int conenr)
void rearrange_generic(partset &Tp, int conenr, int firstnr, int lastnr)
void rearrange_snapshot(int argc, char **argv)
void rearrange_lightcone(int argc, char **argv)
void rearrange_read(partset &Tp, int num, int conenr)
void rearrange_fill_treetable(partset &Tp)
void endrun(void)
This function aborts the simulations.
Definition: begrun.cc:395
simparticles Sp
Definition: simulation.h:56
#define TIMEBASE
Definition: constants.h:333
@ MOST_BOUND_PARTICLE_SNAPHOT
Definition: dtypes.h:307
@ MOST_BOUND_PARTICLE_SNAPHOT_REORDERED
Definition: dtypes.h:308
#define Terminate(...)
Definition: macros.h:15
driftfac Driftfac
Definition: main.cc:41
#define TAG_DENS_A
Definition: mpi_utils.h:50
#define TAG_N
Definition: mpi_utils.h:25
int myMPI_Alltoall(const void *sendbuf, int sendcount, MPI_Datatype sendtype, void *recvbuf, int recvcount, MPI_Datatype recvtype, MPI_Comm comm)
Definition: myalltoall.cc:301
int myMPI_Sendrecv(void *sendbuf, size_t sendcount, MPI_Datatype sendtype, int dest, int sendtag, void *recvbuf, size_t recvcount, MPI_Datatype recvtype, int source, int recvtag, MPI_Comm comm, MPI_Status *status)
memory Mem
Definition: main.cc:44
expr log(half arg)
Definition: half.hpp:2745
double mycxxsort_parallel(T *begin, T *end, Comp comp, MPI_Comm comm)