GADGET-4
fof_io.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 <hdf5.h>
17#include <mpi.h>
18#include <sys/stat.h>
19#include <algorithm>
20#include <cmath>
21#include <cstdio>
22#include <cstdlib>
23#include <cstring>
24
25#include "../data/allvars.h"
26#include "../data/dtypes.h"
27#include "../data/mymalloc.h"
28#include "../fof/fof.h"
29#include "../fof/fof_io.h"
30#include "../gitversion/version.h"
31#include "../io/hdf5_util.h"
32#include "../io/io.h"
33#include "../logs/timer.h"
34#include "../main/main.h"
35#include "../main/simulation.h"
36#include "../mpi_utils/mpi_utils.h"
37#include "../sort/parallel_sort.h"
38#include "../subfind/subfind.h"
39#include "../system/system.h"
40
45template <typename partset>
46fof_io<partset>::fof_io(fof<partset> *FoF_ptr, MPI_Comm comm, int format) : IO_Def(comm, format)
47{
48 FoF = FoF_ptr;
49
50 this->N_IO_Fields = 0;
51 this->N_DataGroups = 3;
52 this->header_size = sizeof(catalogue_header);
53 this->header_buf = &catalogue_header;
54 this->type_of_file = FILE_IS_GROUPCAT;
55 sprintf(this->info, "FOF/SUBFIND: writing group catalogue");
56
57 init_field("FLEN", "GroupLen", mem_len_type, file_len_type, READ_IF_PRESENT, 1, A_G, &FoF->Group[0].Len, NULL, GROUPS, 0, 0, 0, 0, 0,
58 0, 0, true);
59
60 init_field("FMAS", "GroupMass", MEM_MY_DOUBLE, FILE_MY_IO_FLOAT, READ_IF_PRESENT, 1, A_G, &FoF->Group[0].Mass, NULL, GROUPS, 0, 0, 0,
61 0, 0, 0, 0);
62
63 init_field("FPOS", "GroupPos", MEM_MY_DOUBLE, FILE_MY_IO_FLOAT, READ_IF_PRESENT, 3, A_G, &FoF->Group[0].Pos[0], NULL, GROUPS, 0, 0,
64 0, 0, 0, 0, 0);
65
66 init_field("FVEL", "GroupVel", MEM_MY_FLOAT, FILE_MY_IO_FLOAT, READ_IF_PRESENT, 3, A_G, &FoF->Group[0].Vel[0], NULL, GROUPS, 0, 0, 0,
67 0, 0, 0, 0);
68
69 init_field("FLTY", "GroupLenType", mem_len_type, file_len_type, READ_IF_PRESENT, NTYPES, A_G, &FoF->Group[0].LenType[0], NULL,
70 GROUPS, 0, 0, 0, 0, 0, 0, 0, true);
71
72 init_field("FLOT", "GroupOffsetType", MEM_INT64, FILE_INT64, READ_IF_PRESENT, NTYPES, A_G, &FoF->Group[0].OffsetType[0], NULL,
73 GROUPS, 0, 0, 0, 0, 0, 0, 0, true);
74
75 init_field("FMTY", "GroupMassType", MEM_MY_DOUBLE, FILE_MY_IO_FLOAT, READ_IF_PRESENT, NTYPES, A_G, &FoF->Group[0].MassType[0], NULL,
76 GROUPS, 0, 0, 0, 0, 0, 0, 0);
77
78 init_field("FMAA", "GroupAscale", MEM_MY_DOUBLE, FILE_MY_IO_FLOAT, READ_IF_PRESENT, 1, A_G, &FoF->Group[0].Ascale, NULL, GROUPS, 0,
79 0, 0, 0, 0, 0, 0);
80
81#if defined(SUBFIND_ORPHAN_TREATMENT)
82 init_field("FPEN", "GroupLenPrevMostBnd", MEM_INT, FILE_INT, READ_IF_PRESENT, 1, A_G, &FoF->Group[0].LenPrevMostBnd, NULL, GROUPS, 0,
83 0, 0, 0, 0, 0, 0, true);
84#endif
85
86#ifdef STARFORMATION
87 init_field("FSFR", "GroupSFR", MEM_MY_FLOAT, FILE_MY_IO_FLOAT, SKIP_ON_READ, 1, A_G, &FoF->Group[0].Sfr, NULL, GROUPS, 0, 0, 0, 0, 0,
88 0, 0);
89#endif
90
91#ifdef MERGERTREE
92 init_field("GFLO", "GroupFileOffset", MEM_MY_FILEOFFSET, FILE_NONE, SKIP_ON_READ, 1, A_G, &FoF->Group[0].FileOffset, NULL, GROUPS, 0,
93 0, 0, 0, 0, 0, 0, true);
94#endif
95
96#ifdef SUBFIND
97 // ------------------ additional Group fields -------------------------------------
98
99 init_field("FMM2", "Group_M_Mean200", MEM_MY_FLOAT, FILE_MY_IO_FLOAT, READ_IF_PRESENT, 1, A_G, &FoF->Group[0].M_Mean200, NULL,
100 GROUPS, 0, 0, 0, 0, 0, 0, 0);
101 init_field("FRM2", "Group_R_Mean200", MEM_MY_FLOAT, FILE_MY_IO_FLOAT, READ_IF_PRESENT, 1, A_G, &FoF->Group[0].R_Mean200, NULL,
102 GROUPS, 0, 0, 0, 0, 0, 0, 0);
103
104 init_field("FMC2", "Group_M_Crit200", MEM_MY_FLOAT, FILE_MY_IO_FLOAT, READ_IF_PRESENT, 1, A_G, &FoF->Group[0].M_Crit200, NULL,
105 GROUPS, 0, 0, 0, 0, 0, 0, 0);
106 init_field("FRC2", "Group_R_Crit200", MEM_MY_FLOAT, FILE_MY_IO_FLOAT, READ_IF_PRESENT, 1, A_G, &FoF->Group[0].R_Crit200, NULL,
107 GROUPS, 0, 0, 0, 0, 0, 0, 0);
108
109 init_field("FMC5", "Group_M_Crit500", MEM_MY_FLOAT, FILE_MY_IO_FLOAT, READ_IF_PRESENT, 1, A_G, &FoF->Group[0].M_Crit500, NULL,
110 GROUPS, 0, 0, 0, 0, 0, 0, 0);
111 init_field("FRC5", "Group_R_Crit500", MEM_MY_FLOAT, FILE_MY_IO_FLOAT, READ_IF_PRESENT, 1, A_G, &FoF->Group[0].R_Crit500, NULL,
112 GROUPS, 0, 0, 0, 0, 0, 0, 0);
113
114 init_field("FMT2", "Group_M_TopHat200", MEM_MY_FLOAT, FILE_MY_IO_FLOAT, READ_IF_PRESENT, 1, A_G, &FoF->Group[0].M_TopHat200, NULL,
115 GROUPS, 0, 0, 0, 0, 0, 0, 0);
116 init_field("FMR2", "Group_R_TopHat200", MEM_MY_FLOAT, FILE_MY_IO_FLOAT, READ_IF_PRESENT, 1, A_G, &FoF->Group[0].R_TopHat200, NULL,
117 GROUPS, 0, 0, 0, 0, 0, 0, 0);
118
119 init_field("FNSH", "GroupNsubs", MEM_INT, FILE_INT, READ_IF_PRESENT, 1, A_G, &FoF->Group[0].Nsubs, NULL, GROUPS, 0, 0, 0, 0, 0, 0, 0,
120 true);
121
122 init_field("FFSH", "GroupFirstSub", MEM_INT64, FILE_INT64, READ_IF_PRESENT, 1, A_G, &FoF->Group[0].FirstSub, NULL, GROUPS, 0, 0, 0,
123 0, 0, 0, 0, true);
124
125 // ------------------- genuine Subhalo fields ------------------------------------
126
127 init_field("SGNR", "SubhaloGroupNr", MEM_INT64, FILE_INT64, READ_IF_PRESENT, 1, A_S, &FoF->Subhalo[0].GroupNr, NULL, SUBGROUPS, 0, 0,
128 0, 0, 0, 0, 0, true);
129
130 init_field("SGRG", "SubhaloRankInGr", MEM_INT, FILE_INT, READ_IF_PRESENT, 1, A_S, &FoF->Subhalo[0].SubRankInGr, NULL, SUBGROUPS, 0,
131 0, 0, 0, 0, 0, 0, true);
132
133 init_field("SLEN", "SubhaloLen", mem_len_type, file_len_type, READ_IF_PRESENT, 1, A_S, &FoF->Subhalo[0].Len, NULL, SUBGROUPS, 0, 0,
134 0, 0, 0, 0, 0, true);
135
136 init_field("SMAS", "SubhaloMass", MEM_MY_FLOAT, FILE_MY_IO_FLOAT, READ_IF_PRESENT, 1, A_S, &FoF->Subhalo[0].Mass, NULL, SUBGROUPS, 0,
137 0, 0, 0, 0, 0, 0);
138
139 init_field("SPOS", "SubhaloPos", MEM_MY_FLOAT, FILE_MY_IO_FLOAT, READ_IF_PRESENT, 3, A_S, &FoF->Subhalo[0].Pos[0], NULL, SUBGROUPS,
140 0, 0, 0, 0, 0, 0, 0);
141
142 init_field("SVEL", "SubhaloVel", MEM_MY_FLOAT, FILE_MY_IO_FLOAT, READ_IF_PRESENT, 3, A_S, &FoF->Subhalo[0].Vel[0], NULL, SUBGROUPS,
143 0, 0, 0, 0, 0, 0, 0);
144
145 init_field("SLTY", "SubhaloLenType", mem_len_type, file_len_type, READ_IF_PRESENT, NTYPES, A_S, &FoF->Subhalo[0].LenType[0], NULL,
146 SUBGROUPS, 0, 0, 0, 0, 0, 0, 0, true);
147
148 init_field("SLOT", "SubhaloOffsetType", MEM_INT64, FILE_INT64, READ_IF_PRESENT, NTYPES, A_S, &FoF->Subhalo[0].OffsetType[0], NULL,
149 SUBGROUPS, 0, 0, 0, 0, 0, 0, 0, true);
150
151 init_field("SMTY", "SubhaloMassType", MEM_MY_FLOAT, FILE_MY_IO_FLOAT, READ_IF_PRESENT, NTYPES, A_S, &FoF->Subhalo[0].MassType[0],
152 NULL, SUBGROUPS, 0, 0, 0, 0, 0, 0, 0);
153
154 init_field("SCMP", "SubhaloCM", MEM_MY_FLOAT, FILE_MY_IO_FLOAT, READ_IF_PRESENT, 3, A_S, &FoF->Subhalo[0].CM[0], NULL, SUBGROUPS, 0,
155 0, 0, 0, 0, 0, 0);
156
157 init_field("SSPI", "SubhaloSpin", MEM_MY_FLOAT, FILE_MY_IO_FLOAT, READ_IF_PRESENT, 3, A_S, &FoF->Subhalo[0].Spin[0], NULL, SUBGROUPS,
158 0, 0, 0, 0, 0, 0, 0);
159
160 init_field("SSPI", "SubhaloVelDisp", MEM_MY_FLOAT, FILE_MY_IO_FLOAT, READ_IF_PRESENT, 1, A_S, &FoF->Subhalo[0].SubVelDisp, NULL,
161 SUBGROUPS, 0, 0, 0, 0, 0, 0, 0);
162
163 init_field("SVMX", "SubhaloVmax", MEM_MY_FLOAT, FILE_MY_IO_FLOAT, READ_IF_PRESENT, 1, A_S, &FoF->Subhalo[0].SubVmax, NULL, SUBGROUPS,
164 0, 0, 0, 0, 0, 0, 0);
165
166 init_field("SVRX", "SubhaloVmaxRad", MEM_MY_FLOAT, FILE_MY_IO_FLOAT, READ_IF_PRESENT, 1, A_S, &FoF->Subhalo[0].SubVmaxRad, NULL,
167 SUBGROUPS, 0, 0, 0, 0, 0, 0, 0);
168
169 init_field("SHMR", "SubhaloHalfmassRad", MEM_MY_FLOAT, FILE_MY_IO_FLOAT, READ_IF_PRESENT, 1, A_S, &FoF->Subhalo[0].SubHalfMassRad,
170 NULL, SUBGROUPS, 0, 0, 0, 0, 0, 0, 0);
171 init_field("SHMT", "SubhaloHalfmassRadType", MEM_MY_FLOAT, FILE_MY_IO_FLOAT, READ_IF_PRESENT, NTYPES, A_S,
172 &FoF->Subhalo[0].SubHalfMassRadType[0], NULL, SUBGROUPS, 0, 0, 0, 0, 0, 0, 0);
173
174 init_field("SIDM", "SubhaloIDMostbound", MEM_MY_ID_TYPE, FILE_MY_ID_TYPE, READ_IF_PRESENT, 1, A_S, &FoF->Subhalo[0].SubMostBoundID,
175 NULL, SUBGROUPS, 0, 0, 0, 0, 0, 0, 0, true);
176
177 init_field("SPRT", "SubhaloParentRank", MEM_INT, FILE_INT, READ_IF_PRESENT, 1, A_S, &FoF->Subhalo[0].SubParentRank, NULL, SUBGROUPS,
178 0, 0, 0, 0, 0, 0, 0, true);
179
180#ifdef MERGERTREE
181 init_field("SFLO", "SubhaloFileOffset", MEM_MY_FILEOFFSET, FILE_NONE, SKIP_ON_READ, 1, A_S, &FoF->Subhalo[0].FileOffset, NULL,
182 SUBGROUPS, 0, 0, 0, 0, 0, 0, 0, true);
183#endif
184
185#if defined(SUBFIND_ORPHAN_TREATMENT)
186 init_field("LPMO", "SubhaloLenPrevMostBnd", MEM_INT, FILE_INT, READ_IF_PRESENT, 1, A_S, &FoF->Subhalo[0].SubhaloLenPrevMostBnd, NULL,
187 SUBGROUPS, 0, 0, 0, 0, 0, 0, 0, true);
188#endif
189
190#ifdef STARFORMATION
191 init_field("SSFR", "SubhaloSFR", MEM_MY_FLOAT, FILE_MY_IO_FLOAT, READ_IF_PRESENT, 1, A_S, &FoF->Subhalo[0].Sfr, NULL, SUBGROUPS, 0,
192 0, 0, 0, 0, 0, 0);
193 init_field("SSFG", "SubhaloGasMassSFR", MEM_MY_FLOAT, FILE_MY_IO_FLOAT, READ_IF_PRESENT, 1, A_S, &FoF->Subhalo[0].GasMassSfr, NULL,
194 SUBGROUPS, 0, 0, 0, 0, 0, 0, 0);
195#endif
196
197#endif
198}
199
200template <typename partset>
202{
203 return 0; /* not needed here */
204}
205
206template <typename partset>
207void fof_io<partset>::set_type_of_element(int index, int type)
208{
209 /* empty */
210}
211
212/* prepare list of ids with assigned group numbers */
213
214template <typename partset>
215void fof_io<partset>::fof_subfind_save_groups(int num, const char *basename, const char *grpcat_dirbasename)
216{
217#ifdef DO_NOT_PRODUCE_BIG_OUTPUT
218 mpi_printf("FOF/SUBFIND: We skip saving group catalogues.\n");
219 return;
220#endif
221
222 char buf[MAXLEN_PATH_EXTRA];
223
224 double t0 = Logs.second();
225 reset_io_byte_count();
226
228 {
229 if(ThisTask == 0)
230 {
231 sprintf(buf, "%s/%s_%03d", All.OutputDir, grpcat_dirbasename, num);
232 mkdir(buf, 02755);
233 }
234 MPI_Barrier(Communicator);
235 }
236
238 sprintf(buf, "%s/%s_%03d/%s_%03d", All.OutputDir, grpcat_dirbasename, num, basename, num);
239 else
240 sprintf(buf, "%s%s_%03d", All.OutputDir, basename, num);
241
242 write_multiple_files(buf, All.NumFilesPerSnapshot);
243
244 long long byte_count = get_io_byte_count(), byte_count_all;
245 sumup_longs(1, &byte_count, &byte_count_all, Communicator);
246
247 double t1 = Logs.second();
248
249 mpi_printf("FOF/SUBFIND: Group catalogues saved. took = %g sec, total size %g MB, corresponds to effective I/O rate of %g MB/sec\n",
250 Logs.timediff(t0, t1), byte_count_all / (1024.0 * 1024.0), byte_count_all / (1024.0 * 1024.0) / Logs.timediff(t0, t1));
251}
252
253template <typename partset>
255{
256 FoF->TotNgroups = 0;
257
258 char fname[MAXLEN_PATH_EXTRA], fname_multiple[MAXLEN_PATH_EXTRA];
259
260 sprintf(fname_multiple, "%s/groups_%03d/%s_%03d", All.OutputDir, num, "fof_subhalo_tab", num);
261 sprintf(fname, "%s%s_%03d", All.OutputDir, "fof_subhalo_tab", num);
262
263 int num_files = find_files(fname, fname_multiple);
264
265 if(num_files > 1)
266 strcpy(fname, fname_multiple);
267
268 /* we repeat reading the headers of the files two times. In the first iteration, only the
269 * particle numbers ending up on each processor are assembled, followed by memory allocation.
270 * In the second iteration, the data is actually read in.
271 */
272
273 for(int rep = 0; rep < 2; rep++)
274 {
275 FoF->Ngroups = 0;
276 FoF->Nsubhalos = 0;
277
278 read_files_driver(fname, rep, num_files);
279
280 /* now do the memory allocation */
281 if(rep == 0)
282 {
283 FoF->Group = (typename fof<partset>::group_properties *)Mem.mymalloc_movable(
284 &FoF->Group, "Group", FoF->Ngroups * sizeof(typename fof<partset>::group_properties));
285#ifdef SUBFIND
286 FoF->Subhalo = (typename fof<partset>::subhalo_properties *)Mem.mymalloc_movable(
287 &FoF->Subhalo, "Subhalo", FoF->Nsubhalos * sizeof(typename fof<partset>::subhalo_properties));
288#endif
289 }
290 }
291
292 MPI_Barrier(Communicator);
293
294 mpi_printf("\nFOF/SUBFIND: reading done.\n");
295 mpi_printf("FOF/SUBFIND: Total number of groups read: %lld\n", (long long int)FoF->TotNgroups);
296 mpi_printf("FOF/SUBFIND: Total number of subhalos read: %lld\n\n", (long long int)FoF->TotNsubhalos);
297
298 MPI_Allreduce(MPI_IN_PLACE, &LegacyFormat, 1, MPI_INT, MPI_MAX, Communicator);
299}
300
301template <typename partset>
302void fof_io<partset>::fill_file_header(int writeTask, int lastTask, long long *n_type, long long *ntot_type)
303{
304 /* determine group/id numbers of each type in file */
305
306 n_type[0] = FoF->Ngroups;
307 n_type[1] = FoF->Nsubhalos;
308 n_type[2] = FoF->Nids;
309
310 if(ThisTask == writeTask)
311 {
312 for(int n = 0; n < 3; n++)
313 ntot_type[n] = n_type[n];
314
315 for(int task = writeTask + 1; task <= lastTask; task++)
316 {
317 long long nn[3];
318 MPI_Recv(&nn[0], 3, MPI_LONG_LONG, task, TAG_LOCALN, Communicator, MPI_STATUS_IGNORE);
319 for(int n = 0; n < 3; n++)
320 ntot_type[n] += nn[n];
321 }
322
323 for(int task = writeTask + 1; task <= lastTask; task++)
324 MPI_Send(&ntot_type[0], 3, MPI_LONG_LONG, task, TAG_N, Communicator);
325 }
326 else
327 {
328 MPI_Send(&n_type[0], 3, MPI_LONG_LONG, writeTask, TAG_LOCALN, Communicator);
329 MPI_Recv(&ntot_type[0], 3, MPI_LONG_LONG, writeTask, TAG_N, Communicator, MPI_STATUS_IGNORE);
330 }
331
332 /* fill file header */
333
334 catalogue_header.Ngroups = ntot_type[0];
335 catalogue_header.Nsubhalos = ntot_type[1];
336 catalogue_header.Nids = ntot_type[2];
337
338 catalogue_header.TotNgroups = FoF->TotNgroups;
339 catalogue_header.TotNsubhalos = FoF->TotNsubhalos;
340 catalogue_header.TotNids = FoF->TotNids;
341
342 catalogue_header.num_files = All.NumFilesPerSnapshot;
343
344 catalogue_header.time = All.Time;
346 catalogue_header.redshift = 1.0 / All.Time - 1;
347 else
348 catalogue_header.redshift = 0;
349 catalogue_header.BoxSize = All.BoxSize;
350}
351
352template <typename partset>
353void fof_io<partset>::read_file_header(const char *fname, int filenr, int readTask, int lastTask, long long *n_type,
354 long long *ntot_type, int *nstart)
355{
356 if(ThisTask == readTask)
357 {
358 if(filenr == 0 && nstart == NULL)
359 {
360 mpi_printf(
361 "\nREAD-FOF: filenr=%d, '%s' contains:\n"
362 "READ-FOF: Type 0 (groups): %8lld\n"
363 "READ-FOF: Type 1 (subhalos): %8lld\n"
364 "READ-FOF: Type 2 (ids): %8lld\n",
365 filenr, fname, catalogue_header.Ngroups, catalogue_header.Nsubhalos, catalogue_header.Nids);
366 }
367 }
368
369 if(FoF->TotNgroups == 0)
370 {
371 FoF->TotNgroups = catalogue_header.TotNgroups;
372 FoF->TotNsubhalos = catalogue_header.TotNsubhalos;
373 FoF->TotNids = catalogue_header.TotNids;
374 }
375
376 FoF->Redshift = catalogue_header.redshift;
377 FoF->Time = catalogue_header.time;
378
379 for(int k = 0; k < 3; k++)
380 n_type[k] = ntot_type[k] = 0;
381
382 /* to collect the gas particles all at the beginning (in case several
383 snapshot files are read on the current CPU) we move the collisionless
384 particles such that a gap of the right size is created */
385
386 {
387 ntot_type[0] = catalogue_header.Ngroups;
388
389 long long n_in_file = catalogue_header.Ngroups;
390 int ntask = lastTask - readTask + 1;
391 int n_for_this_task = n_in_file / ntask;
392 if((ThisTask - readTask) < (n_in_file % ntask))
393 n_for_this_task++;
394
395 n_type[0] = n_for_this_task;
396
397 if(nstart)
398 memmove(&FoF->Group[n_for_this_task], &FoF->Group[0], FoF->Ngroups * sizeof(typename fof<partset>::group_properties));
399 }
400
401#ifdef SUBFIND
402 {
403 ntot_type[1] = catalogue_header.Nsubhalos;
404
405 long long n_in_file = catalogue_header.Nsubhalos;
406 int ntask = lastTask - readTask + 1;
407 int n_for_this_task = n_in_file / ntask;
408 if((ThisTask - readTask) < (n_in_file % ntask))
409 n_for_this_task++;
410
411 n_type[1] = n_for_this_task;
412
413 if(nstart)
414 memmove(&FoF->Subhalo[n_for_this_task], &FoF->Subhalo[0], FoF->Nsubhalos * sizeof(typename fof<partset>::subhalo_properties));
415 }
416#endif
417
418 if(nstart)
419 *nstart = 0;
420}
421
422template <typename partset>
424{
425 write_scalar_attribute(handle, "Ngroups_ThisFile", &catalogue_header.Ngroups, H5T_NATIVE_UINT64);
426 write_scalar_attribute(handle, "Nsubhalos_ThisFile", &catalogue_header.Nsubhalos, H5T_NATIVE_UINT64);
427 write_scalar_attribute(handle, "Nids_ThisFile", &catalogue_header.Nids, H5T_NATIVE_UINT64);
428
429 write_scalar_attribute(handle, "Ngroups_Total", &catalogue_header.TotNgroups, H5T_NATIVE_UINT64);
430 write_scalar_attribute(handle, "Nsubhalos_Total", &catalogue_header.TotNsubhalos, H5T_NATIVE_UINT64);
431 write_scalar_attribute(handle, "Nids_Total", &catalogue_header.TotNids, H5T_NATIVE_UINT64);
432
433 write_scalar_attribute(handle, "NumFiles", &catalogue_header.num_files, H5T_NATIVE_INT);
434
435 write_scalar_attribute(handle, "Time", &catalogue_header.time, H5T_NATIVE_DOUBLE);
436 write_scalar_attribute(handle, "Redshift", &catalogue_header.redshift, H5T_NATIVE_DOUBLE);
437 write_scalar_attribute(handle, "BoxSize", &catalogue_header.BoxSize, H5T_NATIVE_DOUBLE);
438
439 write_string_attribute(handle, "Git_commit", GIT_COMMIT);
440 write_string_attribute(handle, "Git_date", GIT_DATE);
441}
442
447template <typename partset>
448void fof_io<partset>::read_header_fields(const char *fname)
449{
450 memset(&catalogue_header, 0, sizeof(fof_subfind_header));
451
452 hid_t hdf5_file = my_H5Fopen(fname, H5F_ACC_RDONLY, H5P_DEFAULT);
453 hid_t handle = my_H5Gopen(hdf5_file, "/Header");
454
455 /* now read the header fields */
456 read_scalar_attribute(handle, "Ngroups_ThisFile", &catalogue_header.Ngroups, H5T_NATIVE_UINT64);
457
458 if(read_scalar_attribute(handle, "Nsubhalos_ThisFile", "Nsubgroups_ThisFile", &catalogue_header.Nsubhalos, H5T_NATIVE_UINT64))
459 LegacyFormat = 1;
460
461 read_scalar_attribute(handle, "Nids_ThisFile", &catalogue_header.Nids, H5T_NATIVE_UINT64);
462
463 read_scalar_attribute(handle, "Ngroups_Total", &catalogue_header.TotNgroups, H5T_NATIVE_UINT64);
464 read_scalar_attribute(handle, "Nsubhalos_Total", "Nsubgroups_Total", &catalogue_header.TotNsubhalos, H5T_NATIVE_UINT64);
465 read_scalar_attribute(handle, "Nids_Total", &catalogue_header.TotNids, H5T_NATIVE_UINT64);
466
467 read_scalar_attribute(handle, "NumFiles", &catalogue_header.num_files, H5T_NATIVE_INT);
468
469 read_scalar_attribute(handle, "Time", &catalogue_header.time, H5T_NATIVE_DOUBLE);
470 read_scalar_attribute(handle, "Redshift", &catalogue_header.redshift, H5T_NATIVE_DOUBLE);
471
472 read_scalar_attribute(handle, "BoxSize", &catalogue_header.BoxSize, H5T_NATIVE_DOUBLE);
473
474 my_H5Gclose(handle, "/Header");
475 my_H5Fclose(hdf5_file, fname);
476}
477
478template <typename partset>
480{
481 return catalogue_header.num_files;
482}
483
484template <typename partset>
486{
487 catalogue_header.num_files = numfiles;
488}
489
490template <typename partset>
491void fof_io<partset>::read_increase_numbers(int type, int n_for_this_task)
492{
493 switch(type)
494 {
495 case 0:
496 FoF->Ngroups += n_for_this_task;
497 break;
498 case 1:
499 FoF->Nsubhalos += n_for_this_task;
500 break;
501 case 2:
502 FoF->Nids += n_for_this_task;
503 break;
504 default:
505 Terminate("wrong group");
506 break;
507 }
508}
509
510template <typename partset>
511void fof_io<partset>::get_datagroup_name(int type, char *buf)
512{
513 switch(type)
514 {
515 case 0:
516 sprintf(buf, "/Group");
517 break;
518 case 1:
519 sprintf(buf, "/Subhalo");
520 break;
521 case 2:
522 sprintf(buf, "/IDs");
523 break;
524 default:
525 Terminate("wrong group: type=%d", type);
526 break;
527 }
528}
529
530template <typename partset>
532{
533 switch(array)
534 {
535 case A_G:
536 return (void *)(FoF->Group + index);
537
538#ifdef SUBFIND
539 case A_S:
540 return (void *)(FoF->Subhalo + index);
541#endif
542
543 default:
544 Terminate("we don't expect to get here");
545 }
546
547 return NULL;
548}
549
550template <typename partset>
552{
553 double t0 = Logs.second();
554
555 ID_list = (id_list *)Mem.mymalloc("ID_list", sizeof(id_list) * FoF->Nids);
556
557 long long nids = 0;
558 for(int i = 0; i < FoF->Tp->NumPart; i++)
559 {
560 if(FoF->Tp->PS[i].GroupNr.get() < HALONR_MAX)
561 {
562 if(nids >= FoF->Nids)
563 Terminate("nids >= Nids");
564
565 ID_list[nids].GroupNr = FoF->Tp->PS[i].GroupNr;
566 ID_list[nids].Type = FoF->Tp->P[i].getType();
567 ID_list[nids].ID = FoF->Tp->P[i].ID.get();
568#ifdef SUBFIND
569 ID_list[nids].SubRankInGr = FoF->Tp->PS[i].SubRankInGr;
570 ID_list[nids].BindingEgy = FoF->Tp->PS[i].v.DM_BindingEnergy;
571#endif
572 nids++;
573 }
574 }
575
576 long long totNids;
577 sumup_longs(1, &nids, &totNids, Communicator);
578 if(totNids != FoF->TotNids)
579 Terminate("Task=%d Nids=%lld totNids=%lld TotNids=%lld\n", ThisTask, FoF->Nids, totNids, FoF->TotNids);
580
581 /* sort the particle IDs according to group-number, and optionally subhalo number and binding energy */
582 mycxxsort_parallel(ID_list, ID_list + FoF->Nids, fof_subfind_compare_ID_list, Communicator);
583
584 double t1 = Logs.second();
585 mpi_printf("FOF/SUBFIND: Particle/cell IDs in groups globally sorted. took = %g sec\n", Logs.timediff(t0, t1));
586}
587
588/* now make sure that the following classes are really instantiated, otherwise we may get a linking problem */
589#include "../data/simparticles.h"
590template class fof_io<simparticles>;
591
592#if defined(LIGHTCONE) && defined(LIGHTCONE_PARTICLES_GROUPS)
593#include "../data/lcparticles.h"
594template class fof_io<lcparticles>;
595#endif
596
597#endif
global_data_all_processes All
Definition: main.cc:40
Definition: io.h:129
Definition: fof_io.h:20
int get_filenr_from_header(void)
void read_header_fields(const char *fname)
void * get_base_address_of_structure(enum arrays array, int index)
void fof_subfind_load_groups(int num)
void read_increase_numbers(int type, int n_for_this_task)
void fill_file_header(int writeTask, int lastTask, long long *nloc_part, long long *npart)
void fof_subfind_save_groups(int num, const char *basename, const char *grpcat_dirbasename)
void read_file_header(const char *fname, int filenr, int readTask, int lastTask, long long *nloc_part, long long *npart, int *nstart)
void set_type_of_element(int index, int type)
void write_header_fields(hid_t)
fof_io(fof< partset > *FoF_ptr, MPI_Comm comm, int format)
void get_datagroup_name(int grnr, char *gname)
void set_filenr_in_header(int)
int get_type_of_element(int index)
double timediff(double t0, double t1)
Definition: logs.cc:488
double second(void)
Definition: logs.cc:471
#define MAXLEN_PATH_EXTRA
Definition: constants.h:301
#define NTYPES
Definition: constants.h:308
hid_t my_H5Gopen(hid_t loc_id, const char *groupname)
Definition: hdf5_util.cc:311
hid_t my_H5Fopen(const char *fname, unsigned int flags, hid_t fapl_id)
Definition: hdf5_util.cc:297
herr_t my_H5Gclose(hid_t group_id, const char *groupname)
Definition: hdf5_util.cc:457
herr_t my_H5Fclose(hid_t file_id, const char *fname)
Definition: hdf5_util.cc:482
void write_scalar_attribute(hid_t handle, const char *attr_name, const void *buf, hid_t mem_type_id)
Definition: hdf5_util.cc:621
void read_scalar_attribute(hid_t handle, const char *attr_name, void *buf, hid_t mem_type_id)
Definition: hdf5_util.cc:589
void write_string_attribute(hid_t handle, const char *attr_name, const char *buf)
Definition: hdf5_util.cc:654
#define HALONR_MAX
Definition: idstorage.h:20
const types_in_memory mem_len_type
Definition: io.h:94
@ READ_IF_PRESENT
Definition: io.h:124
@ SKIP_ON_READ
Definition: io.h:125
@ GROUPS
Definition: io.h:107
@ SUBGROUPS
Definition: io.h:108
arrays
Definition: io.h:30
@ A_G
Definition: io.h:35
@ A_S
Definition: io.h:36
@ FILE_NONE
Definition: io.h:66
@ FILE_MY_IO_FLOAT
Definition: io.h:69
@ FILE_INT64
Definition: io.h:68
@ FILE_MY_ID_TYPE
Definition: io.h:71
@ FILE_INT
Definition: io.h:67
const types_in_file file_len_type
Definition: io.h:95
@ MEM_INT
Definition: io.h:79
@ MEM_MY_FLOAT
Definition: io.h:85
@ MEM_MY_FILEOFFSET
Definition: io.h:87
@ MEM_MY_ID_TYPE
Definition: io.h:81
@ MEM_MY_DOUBLE
Definition: io.h:86
@ MEM_INT64
Definition: io.h:80
@ FILE_IS_GROUPCAT
Definition: io.h:54
logs Logs
Definition: main.cc:43
#define Terminate(...)
Definition: macros.h:15
#define TAG_LOCALN
Definition: mpi_utils.h:52
#define TAG_N
Definition: mpi_utils.h:25
void sumup_longs(int n, long long *src, long long *res, MPI_Comm comm)
memory Mem
Definition: main.cc:44
double mycxxsort_parallel(T *begin, T *end, Comp comp, MPI_Comm comm)
char OutputDir[MAXLEN_PATH]
Definition: allvars.h:272
const char * GIT_COMMIT
const char * GIT_DATE