GADGET-4
io_readsnap.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 <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 "../system/system.h"
39 
40 readsnap_io::readsnap_io(mergertree *MergerTree_ptr, MPI_Comm comm, int format) : IO_Def(comm, format)
41 {
42  MergerTree = MergerTree_ptr;
43 
44  this->N_IO_Fields = 0;
45  this->N_DataGroups = NTYPES;
46  this->header_size = sizeof(header);
47  this->header_buf = &header;
48  this->type_of_file = FILE_IS_SNAPSHOT;
49  sprintf(this->info, "MERGERTREE: reading snapshot IDs");
50 
51  init_field("ID ", "ParticleIDs", MEM_MY_ID_TYPE, FILE_MY_ID_TYPE, READ_IF_PRESENT, 1, A_MTRP, &MergerTree->MtrP[0].ID, NULL,
52  ALL_TYPES, 0, 0, 0, 0, 0, 0, 0);
53 
54  init_field("FLOF", "FileOffset", MEM_MY_FILEOFFSET, FILE_NONE, SKIP_ON_READ, 1, A_MTRP, &MergerTree->MtrP[0].FileOffset, NULL,
55  ALL_TYPES, 0, 0, 0, 0, 0, 0, 0);
56 }
57 
79 void readsnap_io::mergertree_read_snap_ids(int num)
80 {
81  if(All.ICFormat < 1 || All.ICFormat > 4)
82  Terminate("ICFormat=%d not supported.\n", All.ICFormat);
83 
84  char fname[MAXLEN_PATH_EXTRA], fname_multiple[MAXLEN_PATH_EXTRA];
85  sprintf(fname_multiple, "%s/snapdir_%03d/%s_%03d", All.OutputDir, num, All.SnapshotFileBase, num);
86  sprintf(fname, "%s%s_%03d", All.OutputDir, All.SnapshotFileBase, num);
87 
88  TIMER_START(CPU_SNAPSHOT);
89 
90  int num_files = find_files(fname, fname_multiple);
91 
92  if(num_files > 1)
93  strcpy(fname, fname_multiple);
94 
95  /* we repeat reading the headers of the files two times. In the first iteration, only the
96  * particle numbers ending up on each processor are assembled, followed by memory allocation.
97  * In the second iteration, the data is actually read in.
98  */
99  for(int rep = 0; rep < 2; rep++)
100  {
101  MergerTree->MtrP_NumPart = 0;
102 
103  read_files_driver(fname, rep, num_files);
104 
105  /* now do the memory allocation */
106  if(rep == 0)
107  {
108  MergerTree->MtrP = (mergertree::mergertree_particle_data *)Mem.mymalloc_movable_clear(
109  &MergerTree->MtrP, "MtrP", (MergerTree->MtrP_NumPart + 1) * sizeof(mergertree::mergertree_particle_data));
110  }
111  }
112 
113  MPI_Barrier(Communicator);
114 
115  mpi_printf("READSNAPID: reading done.\n");
116 
117  TIMER_STOP(CPU_SNAPSHOT);
118 }
119 
120 void readsnap_io::fill_file_header(int writeTask, int lastTask, long long *n_type, long long *ntot_type)
121 { /* empty */
122 }
123 
124 void readsnap_io::read_file_header(const char *fname, int filenr, int readTask, int lastTask, long long *n_type, long long *ntot_type,
125  int *nstart)
126 {
127  if(ThisTask == readTask)
128  {
129  if(filenr == 0 && nstart == NULL)
130  {
131  mpi_printf(
132  "\nREADSNAPID: filenr=%d, '%s' contains:\n"
133  "READSNAPID: Type 0 (gas): %8lld (tot=%15lld) masstab= %g\n",
134  filenr, fname, (long long)header.npart[0], (long long)header.npartTotal[0], All.MassTable[0]);
135 
136  for(int type = 1; type < NTYPES; type++)
137  {
138  mpi_printf("READSNAPID: Type %d: %8lld (tot=%15lld) masstab= %g\n", type, (long long)header.npart[type],
139  (long long)header.npartTotal[type], All.MassTable[type]);
140  }
141  mpi_printf("\n");
142  }
143  }
144 
145  /* to collect the gas particles all at the beginning (in case several
146  snapshot files are read on the current CPU) we move the collisionless
147  particles such that a gap of the right size is created */
148 
149  long long nall = 0;
150  for(int type = 0; type < NTYPES; type++)
151  {
152  ntot_type[type] = header.npart[type];
153 
154  long long n_in_file = header.npart[type];
155  int ntask = lastTask - readTask + 1;
156  int n_for_this_task = n_in_file / ntask;
157  if((ThisTask - readTask) < (n_in_file % ntask))
158  n_for_this_task++;
159 
160  n_type[type] = n_for_this_task;
161 
162  nall += n_for_this_task;
163  }
164 
165  if(nstart)
166  {
167  memmove(&MergerTree->MtrP[nall], &MergerTree->MtrP[0], MergerTree->MtrP_NumPart * sizeof(mergertree::mergertree_particle_data));
168  *nstart = 0;
169  }
170 }
171 
172 void readsnap_io::write_header_fields(hid_t handle)
173 { /* empty */
174 }
175 
180 void readsnap_io::read_header_fields(const char *fname)
181 {
182  for(int i = 0; i < NTYPES; i++)
183  {
184  header.npart[i] = 0;
185  header.npartTotal[i] = 0;
186  header.mass[i] = 0;
187  }
188 
189  hsize_t ntypes = NTYPES;
190 
191  hid_t hdf5_file = my_H5Fopen(fname, H5F_ACC_RDONLY, H5P_DEFAULT);
192  hid_t handle = my_H5Gopen(hdf5_file, "/Header");
193 
194  /* check if the file in question actually has this number of types */
195  hid_t hdf5_attribute = my_H5Aopen_name(handle, "NumPart_ThisFile");
196  hid_t space = H5Aget_space(hdf5_attribute);
197  hsize_t dims, len;
198  H5Sget_simple_extent_dims(space, &dims, &len);
199  H5Sclose(space);
200  if(len != ntypes)
201  Terminate("Length of NumPart_ThisFile attribute (%d) does not match NTYPES(ICS) (%d)", (int)len, (int)ntypes);
202  my_H5Aclose(hdf5_attribute, "NumPart_ThisFile");
203 
204  /* now read the header fields */
205 
206 #ifdef GADGET2_HEADER
207  read_vector_attribute(handle, "NumPart_ThisFile", header.npart, H5T_NATIVE_UINT, ntypes);
208 #else
209  read_vector_attribute(handle, "NumPart_ThisFile", header.npart, H5T_NATIVE_UINT64, ntypes);
210 #endif
211 
212  read_vector_attribute(handle, "NumPart_Total", header.npartTotal, H5T_NATIVE_UINT64, ntypes);
213 
214  read_scalar_attribute(handle, "BoxSize", &header.BoxSize, H5T_NATIVE_DOUBLE);
215  read_vector_attribute(handle, "MassTable", header.mass, H5T_NATIVE_DOUBLE, ntypes);
216  read_scalar_attribute(handle, "Time", &header.time, H5T_NATIVE_DOUBLE);
217  read_scalar_attribute(handle, "Redshift", &header.redshift, H5T_NATIVE_DOUBLE);
218  read_scalar_attribute(handle, "NumFilesPerSnapshot", &header.num_files, H5T_NATIVE_INT);
219 
220  my_H5Gclose(handle, "/Header");
221  my_H5Fclose(hdf5_file, fname);
222 }
223 
224 int readsnap_io::get_filenr_from_header(void) { return header.num_files; }
225 
226 void readsnap_io::set_filenr_in_header(int numfiles) { header.num_files = numfiles; }
227 
228 void readsnap_io::read_increase_numbers(int type, int n_for_this_task) { MergerTree->MtrP_NumPart += n_for_this_task; }
229 
230 void readsnap_io::get_datagroup_name(int type, char *buf) { sprintf(buf, "/PartType%d", type); }
231 
232 int readsnap_io::get_type_of_element(int index) { return MergerTree->MtrP[index].Type; }
233 
234 void readsnap_io::set_type_of_element(int index, int type) { MergerTree->MtrP[index].Type = type; }
235 
236 void *readsnap_io::get_base_address_of_structure(enum arrays array, int index)
237 {
238  switch(array)
239  {
240  case A_MTRP:
241  return (void *)(MergerTree->MtrP + index);
242  default:
243  Terminate("we don't expect to get here");
244  }
245 
246  return NULL;
247 }
248 #endif
global_data_all_processes All
Definition: main.cc:40
Definition: io.h:129
#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
void read_vector_attribute(hid_t handle, const char *attr_name, void *buf, hid_t mem_type_id, int length)
Definition: hdf5_util.cc:614
herr_t my_H5Fclose(hid_t file_id, const char *fname)
Definition: hdf5_util.cc:482
hid_t my_H5Aopen_name(hid_t loc_id, const char *attr_name)
Definition: hdf5_util.cc:371
herr_t my_H5Aclose(hid_t attr_id, const char *attr_name)
Definition: hdf5_util.cc:429
void read_scalar_attribute(hid_t handle, const char *attr_name, void *buf, hid_t mem_type_id)
Definition: hdf5_util.cc:589
@ READ_IF_PRESENT
Definition: io.h:124
@ SKIP_ON_READ
Definition: io.h:125
@ ALL_TYPES
Definition: io.h:103
arrays
Definition: io.h:30
@ A_MTRP
Definition: io.h:40
@ FILE_NONE
Definition: io.h:66
@ FILE_MY_ID_TYPE
Definition: io.h:71
@ MEM_MY_FILEOFFSET
Definition: io.h:87
@ MEM_MY_ID_TYPE
Definition: io.h:81
@ FILE_IS_SNAPSHOT
Definition: io.h:53
#define TIMER_START(counter)
Starts the timer counter.
Definition: logs.h:197
#define TIMER_STOP(counter)
Stops the timer counter.
Definition: logs.h:220
#define Terminate(...)
Definition: macros.h:19
memory Mem
Definition: main.cc:44
double MassTable[NTYPES]
Definition: allvars.h:269
char OutputDir[MAXLEN_PATH]
Definition: allvars.h:272
char SnapshotFileBase[MAXLEN_PATH]
Definition: allvars.h:272