12#include "gadgetconfig.h"
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"
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"
45template <
typename partset>
50 this->N_IO_Fields = 0;
51 this->N_DataGroups = 3;
52 this->header_size =
sizeof(catalogue_header);
53 this->header_buf = &catalogue_header;
55 sprintf(this->info,
"FOF/SUBFIND: writing group catalogue");
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,
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,
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,
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,
70 GROUPS, 0, 0, 0, 0, 0, 0, 0,
true);
73 GROUPS, 0, 0, 0, 0, 0, 0, 0,
true);
76 GROUPS, 0, 0, 0, 0, 0, 0, 0);
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,
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);
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,
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);
100 GROUPS, 0, 0, 0, 0, 0, 0, 0);
102 GROUPS, 0, 0, 0, 0, 0, 0, 0);
105 GROUPS, 0, 0, 0, 0, 0, 0, 0);
107 GROUPS, 0, 0, 0, 0, 0, 0, 0);
110 GROUPS, 0, 0, 0, 0, 0, 0, 0);
112 GROUPS, 0, 0, 0, 0, 0, 0, 0);
115 GROUPS, 0, 0, 0, 0, 0, 0, 0);
117 GROUPS, 0, 0, 0, 0, 0, 0, 0);
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,
122 init_field(
"FFSH",
"GroupFirstSub",
MEM_INT64,
FILE_INT64,
READ_IF_PRESENT, 1,
A_G, &FoF->Group[0].FirstSub, NULL,
GROUPS, 0, 0, 0,
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);
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);
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);
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,
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);
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);
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,
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);
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);
172 &FoF->Subhalo[0].SubHalfMassRadType[0], NULL,
SUBGROUPS, 0, 0, 0, 0, 0, 0, 0);
175 NULL,
SUBGROUPS, 0, 0, 0, 0, 0, 0, 0,
true);
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);
185#if defined(SUBFIND_ORPHAN_TREATMENT)
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,
200template <
typename partset>
206template <
typename partset>
214template <
typename partset>
217#ifdef DO_NOT_PRODUCE_BIG_OUTPUT
218 mpi_printf(
"FOF/SUBFIND: We skip saving group catalogues.\n");
225 reset_io_byte_count();
231 sprintf(buf,
"%s/%s_%03d",
All.
OutputDir, grpcat_dirbasename, num);
234 MPI_Barrier(Communicator);
238 sprintf(buf,
"%s/%s_%03d/%s_%03d",
All.
OutputDir, grpcat_dirbasename, num, basename, num);
240 sprintf(buf,
"%s%s_%03d",
All.
OutputDir, basename, num);
244 long long byte_count = get_io_byte_count(), byte_count_all;
245 sumup_longs(1, &byte_count, &byte_count_all, Communicator);
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));
253template <
typename partset>
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);
263 int num_files = find_files(fname, fname_multiple);
266 strcpy(fname, fname_multiple);
273 for(
int rep = 0; rep < 2; rep++)
278 read_files_driver(fname, rep, num_files);
283 FoF->Group = (
typename fof<partset>::group_properties *)
Mem.mymalloc_movable(
284 &FoF->Group,
"Group", FoF->Ngroups *
sizeof(
typename fof<partset>::group_properties));
286 FoF->Subhalo = (
typename fof<partset>::subhalo_properties *)
Mem.mymalloc_movable(
287 &FoF->Subhalo,
"Subhalo", FoF->Nsubhalos *
sizeof(
typename fof<partset>::subhalo_properties));
292 MPI_Barrier(Communicator);
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);
298 MPI_Allreduce(MPI_IN_PLACE, &LegacyFormat, 1, MPI_INT, MPI_MAX, Communicator);
301template <
typename partset>
306 n_type[0] = FoF->Ngroups;
307 n_type[1] = FoF->Nsubhalos;
308 n_type[2] = FoF->Nids;
310 if(ThisTask == writeTask)
312 for(
int n = 0; n < 3; n++)
313 ntot_type[n] = n_type[n];
315 for(
int task = writeTask + 1; task <= lastTask; task++)
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];
323 for(
int task = writeTask + 1; task <= lastTask; task++)
324 MPI_Send(&ntot_type[0], 3, MPI_LONG_LONG, task,
TAG_N, Communicator);
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);
334 catalogue_header.Ngroups = ntot_type[0];
335 catalogue_header.Nsubhalos = ntot_type[1];
336 catalogue_header.Nids = ntot_type[2];
338 catalogue_header.TotNgroups = FoF->TotNgroups;
339 catalogue_header.TotNsubhalos = FoF->TotNsubhalos;
340 catalogue_header.TotNids = FoF->TotNids;
344 catalogue_header.time =
All.
Time;
346 catalogue_header.redshift = 1.0 /
All.
Time - 1;
348 catalogue_header.redshift = 0;
352template <
typename partset>
354 long long *ntot_type,
int *nstart)
356 if(ThisTask == readTask)
358 if(filenr == 0 && nstart == NULL)
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);
369 if(FoF->TotNgroups == 0)
371 FoF->TotNgroups = catalogue_header.TotNgroups;
372 FoF->TotNsubhalos = catalogue_header.TotNsubhalos;
373 FoF->TotNids = catalogue_header.TotNids;
376 FoF->Redshift = catalogue_header.redshift;
377 FoF->Time = catalogue_header.time;
379 for(
int k = 0; k < 3; k++)
380 n_type[k] = ntot_type[k] = 0;
387 ntot_type[0] = catalogue_header.Ngroups;
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))
395 n_type[0] = n_for_this_task;
398 memmove(&FoF->Group[n_for_this_task], &FoF->Group[0], FoF->Ngroups *
sizeof(
typename fof<partset>::group_properties));
403 ntot_type[1] = catalogue_header.Nsubhalos;
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))
411 n_type[1] = n_for_this_task;
414 memmove(&FoF->Subhalo[n_for_this_task], &FoF->Subhalo[0], FoF->Nsubhalos *
sizeof(
typename fof<partset>::subhalo_properties));
422template <
typename partset>
447template <
typename partset>
450 memset(&catalogue_header, 0,
sizeof(fof_subfind_header));
452 hid_t hdf5_file =
my_H5Fopen(fname, H5F_ACC_RDONLY, H5P_DEFAULT);
453 hid_t handle =
my_H5Gopen(hdf5_file,
"/Header");
458 if(
read_scalar_attribute(handle,
"Nsubhalos_ThisFile",
"Nsubgroups_ThisFile", &catalogue_header.Nsubhalos, H5T_NATIVE_UINT64))
464 read_scalar_attribute(handle,
"Nsubhalos_Total",
"Nsubgroups_Total", &catalogue_header.TotNsubhalos, H5T_NATIVE_UINT64);
478template <
typename partset>
481 return catalogue_header.num_files;
484template <
typename partset>
487 catalogue_header.num_files = numfiles;
490template <
typename partset>
496 FoF->Ngroups += n_for_this_task;
499 FoF->Nsubhalos += n_for_this_task;
502 FoF->Nids += n_for_this_task;
510template <
typename partset>
516 sprintf(buf,
"/Group");
519 sprintf(buf,
"/Subhalo");
522 sprintf(buf,
"/IDs");
530template <
typename partset>
536 return (
void *)(FoF->Group + index);
540 return (
void *)(FoF->Subhalo + index);
544 Terminate(
"we don't expect to get here");
550template <
typename partset>
555 ID_list = (id_list *)
Mem.mymalloc(
"ID_list",
sizeof(id_list) * FoF->Nids);
558 for(
int i = 0; i < FoF->Tp->NumPart; i++)
562 if(nids >= FoF->Nids)
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();
569 ID_list[nids].SubRankInGr = FoF->Tp->PS[i].SubRankInGr;
570 ID_list[nids].BindingEgy = FoF->Tp->PS[i].v.DM_BindingEnergy;
578 if(totNids != FoF->TotNids)
579 Terminate(
"Task=%d Nids=%lld totNids=%lld TotNids=%lld\n", ThisTask, FoF->Nids, totNids, FoF->TotNids);
582 mycxxsort_parallel(ID_list, ID_list + FoF->Nids, fof_subfind_compare_ID_list, Communicator);
585 mpi_printf(
"FOF/SUBFIND: Particle/cell IDs in groups globally sorted. took = %g sec\n",
Logs.
timediff(t0, t1));
589#include "../data/simparticles.h"
592#if defined(LIGHTCONE) && defined(LIGHTCONE_PARTICLES_GROUPS)
593#include "../data/lcparticles.h"
global_data_all_processes All
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)
#define MAXLEN_PATH_EXTRA
hid_t my_H5Gopen(hid_t loc_id, const char *groupname)
hid_t my_H5Fopen(const char *fname, unsigned int flags, hid_t fapl_id)
herr_t my_H5Gclose(hid_t group_id, const char *groupname)
herr_t my_H5Fclose(hid_t file_id, const char *fname)
void write_scalar_attribute(hid_t handle, const char *attr_name, const void *buf, hid_t mem_type_id)
void read_scalar_attribute(hid_t handle, const char *attr_name, void *buf, hid_t mem_type_id)
void write_string_attribute(hid_t handle, const char *attr_name, const char *buf)
const types_in_memory mem_len_type
const types_in_file file_len_type
void sumup_longs(int n, long long *src, long long *res, MPI_Comm comm)
double mycxxsort_parallel(T *begin, T *end, Comp comp, MPI_Comm comm)
int ComovingIntegrationOn
char OutputDir[MAXLEN_PATH]