12#include "gadgetconfig.h"
16#include <gsl/gsl_rng.h>
27#include "../data/allvars.h"
28#include "../data/dtypes.h"
29#include "../data/mymalloc.h"
30#include "../fof/fof.h"
31#include "../io/hdf5_util.h"
33#include "../logs/timer.h"
34#include "../main/simulation.h"
35#include "../mergertree/io_descendant.h"
36#include "../mergertree/mergertree.h"
37#include "../mpi_utils/mpi_utils.h"
38#include "../sort/parallel_sort.h"
39#include "../subfind/subfind.h"
40#include "../system/system.h"
42descendant_io::descendant_io(mergertree *MergerTree_ptr, MPI_Comm comm,
int format) :
IO_Def(comm, format)
44 MergerTree = MergerTree_ptr;
46 this->N_IO_Fields = 0;
47 this->N_DataGroups = 1;
48 this->header_size =
sizeof(header);
49 this->header_buf = &header;
51 sprintf(this->info,
"MERGERTREE: writing descendant information");
53 init_field(
"DSNR",
"DescSubhaloNr",
MEM_INT64,
FILE_INT64,
READ_IF_PRESENT, 1,
A_DESC, NULL, io_func_descsubhalonr,
PREVSUBS, 0, 0,
57 PREVSUBS, 0, 0, 0, 0, 0, 0, 0,
true);
59 init_field(
"NSNR",
"NextProgSubhaloNr",
MEM_INT64,
FILE_INT64,
READ_IF_PRESENT, 1,
A_DESC, NULL, io_func_nextsubhalonr,
PREVSUBS, 0,
60 0, 0, 0, 0, 0, 0,
true);
63 PREVSUBS, 0, 0, 0, 0, 0, 0, 0,
true);
66 NULL,
PREVSUBS, 0, 0, 0, 0, 0, 0, 0,
true);
69void descendant_io::mergertree_save_descendants(
int num)
82 MPI_Barrier(Communicator);
86 sprintf(buf,
"%s/groups_%03d/%s_%03d",
All.
OutputDir, num - 1,
"subhalo_desc", num - 1);
88 sprintf(buf,
"%s%s_%03d",
All.
OutputDir,
"subhalo_desc", num - 1);
93void descendant_io::mergertree_read_descendants(
int num)
97 sprintf(fname_multiple,
"%s/groups_%03d/%s_%03d",
All.
OutputDir, num,
"subhalo_desc", num);
98 sprintf(fname,
"%s%s_%03d",
All.
OutputDir,
"subhalo_desc", num);
102 int num_files = find_files(fname, fname_multiple);
105 strcpy(fname, fname_multiple);
111 for(
int rep = 0; rep < 2; rep++)
115 read_files_driver(fname, rep, num_files);
120 MergerTree->Descendants = (mergertree::desc_list *)
Mem.mymalloc_movable(&MergerTree->Descendants,
"Descendants",
121 Nsubhalos *
sizeof(mergertree::desc_list));
125 MPI_Barrier(Communicator);
128void descendant_io::fill_file_header(
int writeTask,
int lastTask,
long long *n_type,
long long *ntot_type)
132 n_type[0] = MergerTree->PrevNsubhalos;
134 if(ThisTask == writeTask)
136 for(
int n = 0; n < 1; n++)
137 ntot_type[n] = n_type[n];
139 for(
int task = writeTask + 1; task <= lastTask; task++)
142 MPI_Recv(&nn[0], 1, MPI_LONG_LONG, task,
TAG_LOCALN, Communicator, MPI_STATUS_IGNORE);
143 for(
int n = 0; n < 1; n++)
144 ntot_type[n] += nn[n];
147 for(
int task = writeTask + 1; task <= lastTask; task++)
148 MPI_Send(&ntot_type[0], 1, MPI_LONG_LONG, task,
TAG_N, Communicator);
152 MPI_Send(&n_type[0], 1, MPI_LONG_LONG, writeTask,
TAG_LOCALN, Communicator);
153 MPI_Recv(&ntot_type[0], 1, MPI_LONG_LONG, writeTask,
TAG_N, Communicator, MPI_STATUS_IGNORE);
158 header.Nsubhalos = ntot_type[0];
159 header.TotNsubhalos = MergerTree->PrevTotNsubhalos;
163void descendant_io::read_file_header(
const char *fname,
int filenr,
int readTask,
int lastTask,
long long *n_type,
164 long long *ntot_type,
int *nstart)
166 if(ThisTask == readTask)
168 if(filenr == 0 && nstart == NULL)
170 mpi_printf(
"\nREAD-DESCENDANTS: filenr=%d, '%s' contains (subhalos): %8d\n", filenr, fname, header.Nsubhalos);
174 if(TotNsubhalos == 0)
175 TotNsubhalos = header.TotNsubhalos;
177 for(
int k = 0; k < 1; k++)
178 n_type[k] = ntot_type[k] = 0;
184 ntot_type[0] = header.Nsubhalos;
186 long long n_in_file = header.Nsubhalos;
187 int ntask = lastTask - readTask + 1;
188 int n_for_this_task = n_in_file / ntask;
189 if((ThisTask - readTask) < (n_in_file % ntask))
192 n_type[0] = n_for_this_task;
196 memmove(&MergerTree->Descendants[n_for_this_task], &MergerTree->Descendants[0], Nsubhalos *
sizeof(mergertree::desc_list));
201void descendant_io::write_header_fields(hid_t handle)
214void descendant_io::read_header_fields(
const char *fname)
216 memset(&header, 0,
sizeof(io_header));
218 hid_t hdf5_file =
my_H5Fopen(fname, H5F_ACC_RDONLY, H5P_DEFAULT);
219 hid_t handle =
my_H5Gopen(hdf5_file,
"/Header");
222 read_scalar_attribute(handle,
"Nsubhalos_ThisFile",
"Nsubgroups_ThisFile", &header.Nsubhalos, H5T_NATIVE_UINT64);
224 read_scalar_attribute(handle,
"Nsubhalos_Total",
"Nsubgroups_Total", &header.TotNsubhalos, H5T_NATIVE_UINT64);
232int descendant_io::get_filenr_from_header(
void) {
return header.num_files; }
234void descendant_io::set_filenr_in_header(
int numfiles) { header.num_files = numfiles; }
236void descendant_io::read_increase_numbers(
int type,
int n_for_this_task)
241 Nsubhalos += n_for_this_task;
249void descendant_io::get_datagroup_name(
int type,
char *buf)
254 sprintf(buf,
"/Subhalo");
262int descendant_io::get_type_of_element(
int index)
268void descendant_io::set_type_of_element(
int index,
int type)
272void *descendant_io::get_base_address_of_structure(
enum arrays array,
int index)
277 return (
void *)(MergerTree->Descendants + index);
279 Terminate(
"strange, we don't expect to get here");
global_data_all_processes All
#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)
char OutputDir[MAXLEN_PATH]