GADGET-4
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 #include <errno.h>
15 #include <hdf5.h>
16 #include <math.h>
17 #include <mpi.h>
18 #include <stdio.h>
19 #include <stdlib.h>
20 #include <string.h>
21 #include <sys/stat.h>
22 #include <algorithm>
23 
24 #include "../cooling_sfr/cooling.h"
25 #include "../data/allvars.h"
26 #include "../data/dtypes.h"
27 #include "../data/mymalloc.h"
28 #include "../fof/fof.h"
29 #include "../io/hdf5_util.h"
30 #include "../io/io.h"
31 #include "../io/parameters.h"
32 #include "../lightcone/lightcone.h"
33 #include "../logs/timer.h"
34 #include "../main/main.h"
35 #include "../main/simulation.h"
36 #include "../mergertree/mergertree.h"
37 #include "../mpi_utils/mpi_utils.h"
38 #include "../subfind/subfind.h"
39 #include "../system/system.h"
40 
41 #define HALF_ROUND_STYLE 1
42 #include "../half/half.hpp"
43 using half_float::half;
44 
45 /* local functions */
46 
48 {
49  if(N_IO_Fields > 0)
50  Mem.myfree_movable(IO_Fields);
51 }
52 
66 void IO_Def::init_field(const char *label, const char *datasetname, enum types_in_memory type_in_memory,
67  enum types_in_file type_in_file_output, enum read_flags read_flag, int values_per_block, enum arrays array,
68  void *pointer_to_field, void (*io_func)(IO_Def *, int, int, void *, int), int typelist_bitmask, int flagunit,
69  double a, double h, double L, double M, double V, double c, bool compression_on)
70 {
71  const int alloc_step = 5;
72 
73  if(values_per_block < 1) // if we have no values, we don't register the field
74  return;
75 
76  if(Max_IO_Fields == 0)
77  {
78  IO_Fields = (IO_Field *)Mem.mymalloc_movable(&IO_Fields, "IO_Fields", alloc_step * sizeof(IO_Field));
79  Max_IO_Fields = alloc_step;
80  }
81  else if(Max_IO_Fields == N_IO_Fields)
82  {
83  Max_IO_Fields = ((Max_IO_Fields / alloc_step) + 1) * alloc_step;
84  IO_Fields = (IO_Field *)Mem.myrealloc_movable(IO_Fields, Max_IO_Fields * sizeof(IO_Field));
85  }
86 
87  int n = N_IO_Fields++;
88  IO_Field *field = &IO_Fields[n];
89 
90  strncpy(field->label, label, LABEL_LEN);
91  field->label[LABEL_LEN] = 0;
92 
93  strncpy(field->datasetname, datasetname, DATASETNAME_LEN);
94  field->datasetname[DATASETNAME_LEN] = 0;
95 
96  field->type_in_memory = type_in_memory;
97  field->type_in_file_output = type_in_file_output;
98  field->read_flag = read_flag;
99  field->values_per_block = values_per_block;
100  field->typelist = typelist_bitmask;
101 #ifdef ALLOW_HDF5_COMPRESSION
102  field->compression_on = compression_on;
103 #else
104  field->compression_on = false;
105 #endif
106 
107  field->array = array;
108  field->io_func = io_func;
109 
110  if(array == A_NONE)
111  {
112  field->offset = 0;
113  }
114  else
115  {
116  field->offset = (size_t)pointer_to_field - (size_t)get_base_address_of_structure(array, 0);
117  }
118 
119  field->hasunit = flagunit;
120  field->a = a;
121  field->h = h;
122  field->L = L;
123  field->M = M;
124  field->V = V;
125  field->c = c;
126 }
127 
132 int IO_Def::find_files(const char *fname, const char *fname_multiple)
133 {
134  FILE *fd;
135  char buf[200], buf1[200];
136  int dummy, files_found = 0;
137 
138  if(file_format == FILEFORMAT_HDF5)
139  {
140  sprintf(buf, "%s.%d.hdf5", fname_multiple, 0);
141  sprintf(buf1, "%s.hdf5", fname);
142  }
143  else
144  {
145  sprintf(buf, "%s.%d", fname_multiple, 0);
146  sprintf(buf1, "%s", fname);
147  }
148 
149  memset(header_buf, 0, header_size);
150 
151  if(ThisTask == 0)
152  {
153  if((fd = fopen(buf, "r")))
154  {
155  if(file_format == FILEFORMAT_LEGACY1 || file_format == FILEFORMAT_LEGACY2)
156  {
157  if(file_format == FILEFORMAT_LEGACY2)
158  {
159  my_fread(&dummy, sizeof(dummy), 1, fd);
160  my_fread(&dummy, sizeof(dummy), 1, fd);
161  my_fread(&dummy, sizeof(dummy), 1, fd);
162  my_fread(&dummy, sizeof(dummy), 1, fd);
163  }
164 
165  my_fread(&dummy, sizeof(dummy), 1, fd);
167  my_fread(&dummy, sizeof(dummy), 1, fd);
168  }
169  fclose(fd);
170 
171  if(file_format == FILEFORMAT_HDF5)
172  read_header_fields(buf);
173 
174  files_found = 1;
175  }
176  }
177 
178  MPI_Bcast(header_buf, header_size, MPI_BYTE, 0, Communicator);
179  MPI_Bcast(&files_found, 1, MPI_INT, 0, Communicator);
180 
181  if(get_filenr_from_header() > 0)
182  return get_filenr_from_header();
183 
184  if(ThisTask == 0)
185  {
186  if((fd = fopen(buf1, "r")))
187  {
188  if(file_format == FILEFORMAT_LEGACY1 || file_format == FILEFORMAT_LEGACY2)
189  {
190  if(file_format == FILEFORMAT_LEGACY2)
191  {
192  my_fread(&dummy, sizeof(dummy), 1, fd);
193  my_fread(&dummy, sizeof(dummy), 1, fd);
194  my_fread(&dummy, sizeof(dummy), 1, fd);
195  my_fread(&dummy, sizeof(dummy), 1, fd);
196  }
197 
198  my_fread(&dummy, sizeof(dummy), 1, fd);
200  my_fread(&dummy, sizeof(dummy), 1, fd);
201  }
202  fclose(fd);
203 
204  if(file_format == FILEFORMAT_HDF5)
205  read_header_fields(buf1);
206 
208 
209  files_found = 1;
210  }
211  }
212 
213  MPI_Bcast(header_buf, header_size, MPI_BYTE, 0, Communicator);
214  MPI_Bcast(&files_found, 1, MPI_INT, 0, Communicator);
215 
216  if(get_filenr_from_header() > 0)
217  return get_filenr_from_header();
218 
219  if(files_found != 0)
220  Terminate("Have found IC files, but number of files in header seems to be zero\n");
221 
222  Terminate("\nCan't find files, neither as '%s'\nnor as '%s'\n", buf, buf1);
223 
224  return 0;
225 }
226 
227 void IO_Def::read_files_driver(const char *fname, int rep, int num_files)
228 {
229  if(rep == 0)
230  {
232  (long long *)Mem.mymalloc_movable(&ntype_in_files, "ntype_in_files", num_files * N_DataGroups * sizeof(long long));
233  memset(ntype_in_files, 0, num_files * N_DataGroups * sizeof(long long));
234  }
235 
236  void *CommBuffer = Mem.mymalloc("CommBuffer", COMMBUFFERSIZE);
237 
238  int rest_files = num_files;
239 
240  while(rest_files > NTask)
241  {
242  char buf[MAXLEN_PATH];
243 
244  sprintf(buf, "%s.%d", fname, ThisTask + (rest_files - NTask));
245  if(file_format == FILEFORMAT_HDF5)
246  sprintf(buf, "%s.%d.hdf5", fname, ThisTask + (rest_files - NTask));
247 
248  int ngroups = NTask / All.MaxFilesWithConcurrentIO;
250  ngroups++;
251  int groupMaster = (ThisTask / ngroups) * ngroups;
252 
253  for(int gr = 0; gr < ngroups; gr++)
254  {
255  if(ThisTask == (groupMaster + gr)) /* ok, it's this processor's turn */
256  {
257  if(rep == 0)
258  share_particle_number_in_file(buf, ThisTask + (rest_files - NTask), ThisTask, ThisTask);
259  else
260  read_file(buf, ThisTask + (rest_files - NTask), ThisTask, ThisTask, CommBuffer);
261  }
262  MPI_Barrier(Communicator);
263  }
264 
265  rest_files -= NTask;
266  }
267 
268  if(rest_files > 0)
269  {
270  int masterTask, filenr, lastTask;
271 
272  distribute_file(rest_files, &filenr, &masterTask, &lastTask);
273 
274  char buf[MAXLEN_PATH];
275 
276  if(num_files > 1)
277  {
278  sprintf(buf, "%s.%d", fname, filenr);
279  if(file_format == FILEFORMAT_HDF5)
280  sprintf(buf, "%s.%d.hdf5", fname, filenr);
281  }
282  else
283  {
284  sprintf(buf, "%s", fname);
285  if(file_format == FILEFORMAT_HDF5)
286  sprintf(buf, "%s.hdf5", fname);
287  }
288 
289  int ngroups = rest_files / All.MaxFilesWithConcurrentIO;
290  if((rest_files % All.MaxFilesWithConcurrentIO))
291  ngroups++;
292 
293  for(int gr = 0; gr < ngroups; gr++)
294  {
295  if((filenr / All.MaxFilesWithConcurrentIO) == gr) /* ok, it's this processor's turn */
296  {
297  if(rep == 0)
298  share_particle_number_in_file(buf, filenr, masterTask, lastTask);
299  else
300  read_file(buf, filenr, masterTask, lastTask, CommBuffer);
301  }
302  MPI_Barrier(Communicator);
303  }
304  }
305 
306  Mem.myfree(CommBuffer);
307 
308  if(rep == 0)
309  MPI_Allreduce(MPI_IN_PLACE, ntype_in_files, num_files * N_DataGroups, MPI_LONG_LONG, MPI_MAX, Communicator);
310  else
311  {
312  /* we are done */
313  Mem.myfree_movable(ntype_in_files);
314  ntype_in_files = NULL;
315  }
316 }
317 
330 void IO_Def::share_particle_number_in_file(const char *fname, int filenr, int readTask, int lastTask)
331 {
332  long long n_type[N_DataGroups], npart[N_DataGroups];
333  unsigned int blksize1, blksize2;
334 
335  if(ThisTask == readTask)
336  {
337  if(file_format == FILEFORMAT_HDF5)
338  {
339  read_header_fields(fname);
340  }
341  else if(file_format == FILEFORMAT_LEGACY1 || file_format == FILEFORMAT_LEGACY2)
342  {
343  FILE *fd = 0;
344 
345  if(!(fd = fopen(fname, "r")))
346  Terminate("can't open file `%s' for reading initial conditions.\n", fname);
347 
348  if(file_format == FILEFORMAT_LEGACY2)
349  {
350  char label[LABEL_LEN];
351  int nextblock;
352  my_fread(&blksize1, sizeof(int), 1, fd);
353  my_fread(&label, sizeof(char), LABEL_LEN, fd);
354  my_fread(&nextblock, sizeof(int), 1, fd);
355  printf("%s Reading header => '%c%c%c%c' (%d byte)\n", info, label[0], label[1], label[2], label[3], nextblock);
356  my_fread(&blksize2, sizeof(int), 1, fd);
357  }
358 
359  my_fread(&blksize1, sizeof(int), 1, fd);
361  my_fread(&blksize2, sizeof(int), 1, fd);
362 
363 #ifdef GADGET2_HEADER
364  if(blksize1 != 256 || blksize2 != 256)
365  Terminate("incorrect GADGET2 header format, blksize1=%d blksize2=%d header_size=%d\n", blksize1, blksize2,
366  (int)header_size);
367 #else
368  if(blksize1 != blksize2)
369  Terminate("incorrect header format, blksize1=%d blksize2=%d header_size=%d \n%s \n", blksize1, blksize2, (int)header_size,
370  blksize1 == 256 ? "You may need to set GADGET2_HEADER" : "");
371 #endif
372  fclose(fd);
373  }
374  else
375  Terminate("illegal format");
376 
377  for(int task = readTask + 1; task <= lastTask; task++)
378  MPI_Ssend(header_buf, header_size, MPI_BYTE, task, TAG_HEADER, Communicator);
379  }
380  else
381  MPI_Recv(header_buf, header_size, MPI_BYTE, readTask, TAG_HEADER, Communicator, MPI_STATUS_IGNORE);
382 
383  read_file_header(fname, filenr, readTask, lastTask, n_type, npart, NULL);
384 
385  if(ThisTask == readTask)
386  {
387  mpi_printf("READIC: Reading file `%s' on task=%d and distribute it to %d to %d.\n", fname, ThisTask, readTask, lastTask);
388  myflush(stdout);
389  }
390 
391  for(int type = 0; type < N_DataGroups; type++)
392  {
393  ntype_in_files[filenr * N_DataGroups + type] = npart[type];
394 
395  long long n_in_file = npart[type];
396  int ntask = lastTask - readTask + 1;
397  int n_for_this_task = n_in_file / ntask;
398  if((ThisTask - readTask) < (n_in_file % ntask))
399  n_for_this_task++;
400 
401  read_increase_numbers(type, n_for_this_task);
402  }
403 }
404 
414 void IO_Def::fill_write_buffer(int blocknr, int *startindex, int pc, int type, void *CommBuffer)
415 {
416  if(blocknr < 0 || blocknr >= N_IO_Fields)
417  Terminate("something is wrong here: blocknr=%d N_IO_Fields=%d", blocknr, N_IO_Fields);
418 
419  IO_Field *field = &IO_Fields[blocknr];
420 
421  int *intp = (int *)CommBuffer;
422  long long *longp = (long long *)CommBuffer;
423  float *floatp = (float *)CommBuffer;
424  double *doublep = (double *)CommBuffer;
425 
426  MyIDType *ip = (MyIDType *)CommBuffer;
427  MyFloat *fp = (MyFloat *)CommBuffer;
428  MyDouble *dp = (MyDouble *)CommBuffer;
429  MyIntPosType *intposp = (MyIntPosType *)CommBuffer;
430 
431  int pindex = *startindex;
432 
433  for(int n = 0; n < pc; pindex++)
434  {
435  if(type_of_file == FILE_IS_SNAPSHOT && type < NTYPES && get_type_of_element(pindex) != type)
436  continue;
437 #ifdef LIGHTCONE
438  if(type_of_file == FILE_IS_LIGHTCONE && type < NTYPES && get_type_of_element(pindex) != type)
439  continue;
440 #endif
441  if(field->io_func)
442  {
443  switch(field->type_in_memory)
444  {
445  case MEM_INT:
446  field->io_func(this, pindex, field->values_per_block, intp, 0);
447  intp += field->values_per_block;
448  break;
449  case MEM_INT64:
450  field->io_func(this, pindex, field->values_per_block, longp, 0);
451  longp += field->values_per_block;
452  break;
453  case MEM_MY_ID_TYPE:
454  field->io_func(this, pindex, field->values_per_block, ip, 0);
455  ip += field->values_per_block;
456  break;
457  case MEM_MY_INTPOS_TYPE:
458  field->io_func(this, pindex, field->values_per_block, intposp, 0);
459  intposp += field->values_per_block;
460  break;
461  case MEM_FLOAT:
462  field->io_func(this, pindex, field->values_per_block, floatp, 0);
463  floatp += field->values_per_block;
464  break;
465  case MEM_DOUBLE:
466  field->io_func(this, pindex, field->values_per_block, doublep, 0);
467  doublep += field->values_per_block;
468  break;
469  case MEM_MY_FLOAT:
470  field->io_func(this, pindex, field->values_per_block, fp, 0);
471  fp += field->values_per_block;
472  break;
473  case MEM_MY_DOUBLE:
474  field->io_func(this, pindex, field->values_per_block, dp, 0);
475  dp += field->values_per_block;
476  break;
477  default:
478  Terminate("ERROR in fill_write_buffer: Type not found!\n");
479 
480  break;
481  }
482  }
483  else
484  {
485  void *array_pos;
486 
487  switch(field->array)
488  {
489  case A_NONE:
490  array_pos = NULL;
491  break;
492 
493  default:
494  array_pos = get_base_address_of_structure(field->array, pindex);
495  break;
496  }
497 
498  for(int k = 0; k < field->values_per_block; k++)
499  {
500  switch(field->type_in_memory)
501  {
502  case MEM_INT:
503  *intp++ = *((int *)((size_t)array_pos + field->offset + k * sizeof(int)));
504  break;
505 
506  case MEM_INT64:
507  *longp++ = *((long long *)((size_t)array_pos + field->offset + k * sizeof(long long)));
508  break;
509 
510  case MEM_MY_ID_TYPE:
511  *ip++ = *((MyIDType *)((size_t)array_pos + field->offset + k * sizeof(MyIDType)));
512  break;
513 
514  case MEM_MY_INTPOS_TYPE:
515  *intposp++ = *((MyIntPosType *)((size_t)array_pos + field->offset + k * sizeof(MyIntPosType)));
516  break;
517 
518  case MEM_FLOAT:
519  *floatp++ = *((float *)((size_t)array_pos + field->offset + k * sizeof(float)));
520  break;
521 
522  case MEM_DOUBLE:
523  *doublep++ = *((double *)((size_t)array_pos + field->offset + k * sizeof(double)));
524  break;
525 
526  case MEM_MY_FLOAT:
527  *fp++ = *((MyFloat *)((size_t)array_pos + field->offset + k * sizeof(MyFloat)));
528  break;
529 
530  case MEM_MY_DOUBLE:
531  *dp++ = *((MyDouble *)((size_t)array_pos + field->offset + k * sizeof(MyDouble)));
532  break;
533 
534  default:
535  Terminate("ERROR in fill_write_buffer: Type not found!\n");
536  break;
537  }
538  }
539  }
540 
541  n++;
542  }
543 
544  *startindex = pindex;
545 }
546 
556 void IO_Def::empty_read_buffer(int blocknr, int offset, int pc, int type, long long nprevious, void *CommBuffer)
557 {
558  IO_Field *field = &IO_Fields[blocknr];
559 
560  int *intp = (int *)CommBuffer;
561  long long *longp = (long long *)CommBuffer;
562  float *floatp = (float *)CommBuffer;
563  double *doublep = (double *)CommBuffer;
564 
565  MyIDType *ip = (MyIDType *)CommBuffer;
566  MyFloat *fp = (MyFloat *)CommBuffer;
567  MyDouble *dp = (MyDouble *)CommBuffer;
568  MyIntPosType *intposp = (MyIntPosType *)CommBuffer;
569 
570  if(field->read_flag != SKIP_ON_READ || field->type_in_memory == MEM_MY_FILEOFFSET)
571  {
572  for(int n = 0; n < pc; n++)
573  {
574  if(field->io_func)
575  {
576  switch(field->type_in_memory)
577  {
578  case MEM_INT:
579  field->io_func(this, offset + n, field->values_per_block, intp, 1);
580  intp += field->values_per_block;
581  break;
582  case MEM_INT64:
583  field->io_func(this, offset + n, field->values_per_block, longp, 1);
584  longp += field->values_per_block;
585  break;
586  case MEM_MY_ID_TYPE:
587  field->io_func(this, offset + n, field->values_per_block, ip, 1);
588  ip += field->values_per_block;
589  break;
590  case MEM_MY_INTPOS_TYPE:
591  field->io_func(this, offset + n, field->values_per_block, intposp, 1);
592  intposp += field->values_per_block;
593  break;
594  case MEM_FLOAT:
595  field->io_func(this, offset + n, field->values_per_block, floatp, 1);
596  floatp += field->values_per_block;
597  break;
598  case MEM_DOUBLE:
599  field->io_func(this, offset + n, field->values_per_block, doublep, 1);
600  doublep += field->values_per_block;
601  break;
602  case MEM_MY_FLOAT:
603  field->io_func(this, offset + n, field->values_per_block, fp, 1);
604  fp += field->values_per_block;
605  break;
606  case MEM_MY_DOUBLE:
607  field->io_func(this, offset + n, field->values_per_block, dp, 1);
608  dp += field->values_per_block;
609  break;
610  case MEM_MY_FILEOFFSET:
611  Terminate("undefined");
612  break;
613  }
614  }
615  else
616  {
617  void *array_pos;
618  switch(field->array)
619  {
620  case A_NONE:
621  array_pos = 0;
622  break;
623 
624  default:
625  array_pos = get_base_address_of_structure(field->array, offset + n);
626  break;
627  }
628 
629  for(int k = 0; k < field->values_per_block; k++)
630  {
631  switch(field->type_in_memory)
632  {
633  case MEM_INT:
634  *((int *)((size_t)array_pos + field->offset + k * sizeof(int))) = *intp++;
635  break;
636  case MEM_INT64:
637  *((long long *)((size_t)array_pos + field->offset + k * sizeof(long long))) = *longp++;
638  break;
639  case MEM_MY_ID_TYPE:
640  *((MyIDType *)((size_t)array_pos + field->offset + k * sizeof(MyIDType))) = *ip++;
641  break;
642  case MEM_MY_INTPOS_TYPE:
643  *((MyIntPosType *)((size_t)array_pos + field->offset + k * sizeof(MyIntPosType))) = *intposp++;
644  break;
645  case MEM_FLOAT:
646  *((float *)((size_t)array_pos + field->offset + k * sizeof(float))) = *floatp++;
647  break;
648  case MEM_DOUBLE:
649  *((double *)((size_t)array_pos + field->offset + k * sizeof(double))) = *doublep++;
650  break;
651  case MEM_MY_FLOAT:
652  *((MyFloat *)((size_t)array_pos + field->offset + k * sizeof(MyFloat))) = *fp++;
653  break;
654  case MEM_MY_DOUBLE:
655  *((MyDouble *)((size_t)array_pos + field->offset + k * sizeof(MyDouble))) = *dp++;
656  break;
657  case MEM_MY_FILEOFFSET:
658  *((long long *)((size_t)array_pos + field->offset + k * sizeof(long long))) = nprevious++;
659  break;
660  default:
661  Terminate("ERROR: Type not found!\n");
662  break;
663  }
664  }
665  }
666  }
667  }
668 }
669 
670 void IO_Def::polling(int numfilesperdump)
671 {
672  if(ThisTask == 0)
673  if(files_completed < numfilesperdump)
674  {
675  MPI_Status status;
676  int flag;
677 
678  /* now check for a completion message */
679  MPI_Iprobe(MPI_ANY_SOURCE, TAG_KEY, Communicator, &flag, &status);
680 
681  if(flag)
682  {
683  int source = status.MPI_SOURCE;
684 
685  int dummy;
686  MPI_Recv(&dummy, 1, MPI_INT, source, TAG_KEY, Communicator, MPI_STATUS_IGNORE);
687  files_completed++;
688 
689  if(files_started < numfilesperdump)
690  {
691  /* send start signal */
692  MPI_Ssend(&ThisTask, 1, MPI_INT, seq[files_started++].thistask, TAG_N, Communicator);
693  }
694  }
695  }
696 }
697 
698 /* driver routine for outputting multiple files, scheduled in optimal order under the constraint not to write more than a certain
699  * number of files simultanuously */
700 void IO_Def::write_multiple_files(char *fname, int numfilesperdump, int append_flag, int chunksize)
701 {
702  if(ThisTask == 0)
703  if(!(seq = (seq_data *)Mem.mymalloc("seq", NTask * sizeof(seq_data))))
704  Terminate("can't allocate seq_data");
705 
706  void *CommBuffer = Mem.mymalloc("CommBuffer", COMMBUFFERSIZE);
707 
708  /* assign processors to output files */
709  int filenr, masterTask, lastTask;
710  distribute_file(numfilesperdump, &filenr, &masterTask, &lastTask);
711 
712  char buf[MAXLEN_PATH];
713  if(numfilesperdump > 1)
714  sprintf(buf, "%s.%d", fname, filenr);
715  else
716  sprintf(buf, "%s", fname);
717 
718  seq_data seq_loc;
719  seq_loc.thistask = ThisTask;
720  seq_loc.rankinnode = RankInThisNode;
721  seq_loc.thisnode = ThisNode;
722 
723  if(masterTask != ThisTask)
724  seq_loc.thistask = -1;
725 
726  MPI_Gather(&seq_loc, sizeof(seq_data), MPI_BYTE, seq, sizeof(seq_data), MPI_BYTE, 0, Communicator);
727 
728  if(ThisTask == 0)
729  {
730  int count = NTask;
731  for(int i = 0; i < count; i++)
732  {
733  if(seq[i].thistask < 0)
734  {
735  count--;
736  seq[i] = seq[count];
737  i--;
738  }
739  }
740  if(count != numfilesperdump)
741  Terminate("count=%d != numfilesperdump=%d", count, numfilesperdump);
742 
743  std::sort(seq, seq + numfilesperdump);
744 
745  files_started = 0;
746  files_completed = 0;
747 
748  for(int i = 1; i < std::min<int>(All.MaxFilesWithConcurrentIO, numfilesperdump); i++)
749  {
750  files_started++;
751  MPI_Ssend(&ThisTask, 1, MPI_INT, seq[i].thistask, TAG_N, Communicator);
752  }
753 
754  files_started++;
755  if(append_flag)
756  append_file(buf, masterTask, lastTask, CommBuffer, numfilesperdump, chunksize);
757  else
758  write_file(buf, masterTask, lastTask, CommBuffer, numfilesperdump, chunksize);
759  files_completed++;
760 
761  if(files_started < numfilesperdump)
762  {
763  /* send start signal */
764  MPI_Ssend(&ThisTask, 1, MPI_INT, seq[files_started++].thistask, TAG_N, Communicator);
765  }
766 
767  while(files_completed < numfilesperdump)
768  polling(numfilesperdump);
769  }
770  else if(masterTask == ThisTask)
771  {
772  /* wait for start signal */
773  int dummy;
774  MPI_Recv(&dummy, 1, MPI_INT, 0, TAG_N, Communicator, MPI_STATUS_IGNORE); /* wait until we are told to start */
775 
776  if(append_flag)
777  append_file(buf, masterTask, lastTask, CommBuffer, numfilesperdump, chunksize);
778  else
779  write_file(buf, masterTask, lastTask, CommBuffer, numfilesperdump, chunksize);
780 
781  /* send back completion notice */
782  MPI_Ssend(&ThisTask, 1, MPI_INT, 0, TAG_KEY, Communicator);
783  }
784  else
785  {
786  if(append_flag)
787  append_file(buf, masterTask, lastTask, CommBuffer, numfilesperdump, chunksize);
788  else
789  write_file(buf, masterTask, lastTask, CommBuffer, numfilesperdump, chunksize);
790  }
791 
792  Mem.myfree(CommBuffer);
793 
794  if(ThisTask == 0)
795  Mem.myfree(seq);
796 }
797 
817 void IO_Def::write_file(char *fname, int writeTask, int lastTask, void *CommBuffer, int numfilesperdump, int chunksize)
818 {
819  int typelist[N_DataGroups];
820  long long n_type[N_DataGroups], npart[N_DataGroups];
821  char label[LABEL_LEN + 1];
822  unsigned int blksize, bytes_per_blockelement_in_file = 0;
823  FILE *fd = 0;
824  hid_t hdf5_file = 0, hdf5_grp[N_DataGroups], hdf5_headergrp = 0, hdf5_dataspace_memory;
825  hid_t hdf5_dataspace_in_file = 0, hdf5_dataset = 0, hdf5_prop = 0;
826  hsize_t dims[2], count[2], start[2];
827  int rank = 0, pcsum = 0;
828  hid_t hdf5_paramsgrp = 0;
829  hid_t hdf5_configgrp = 0;
830 
831 #define SKIP \
832  { \
833  my_fwrite(&blksize, sizeof(int), 1, fd); \
834  }
835 
836  fill_file_header(writeTask, lastTask, n_type, npart);
837 
838  /* open file and write header */
839  if(ThisTask == writeTask)
840  {
841  if(file_format == FILEFORMAT_HDF5)
842  {
843  char buf[MAXLEN_PATH];
844  sprintf(buf, "%s.hdf5", fname);
845  mpi_printf("%s file: '%s' (file 1 of %d)\n", info, fname, numfilesperdump);
846 
847  rename_file_to_bak_if_it_exists(buf);
848 
849  hdf5_file = my_H5Fcreate(buf, H5F_ACC_TRUNC, H5P_DEFAULT, H5P_DEFAULT);
850 
851  hdf5_headergrp = my_H5Gcreate(hdf5_file, "/Header", 0);
852  write_header_fields(hdf5_headergrp);
853 
854  hdf5_paramsgrp = my_H5Gcreate(hdf5_file, "/Parameters", 0);
855  write_parameters_attributes_in_hdf5(hdf5_paramsgrp);
856 
857  hdf5_configgrp = my_H5Gcreate(hdf5_file, "/Config", 0);
858  write_compile_time_options_in_hdf5(hdf5_configgrp);
859 
860  for(int type = 0; type < N_DataGroups; type++)
861  {
862  if(npart[type] > 0)
863  {
864  get_datagroup_name(type, buf);
865  hdf5_grp[type] = my_H5Gcreate(hdf5_file, buf, 0);
866  }
867  }
868  }
869  else
870  {
871  rename_file_to_bak_if_it_exists(fname);
872 
873  if(!(fd = fopen(fname, "w")))
874  Terminate("can't open file `%s' for writing.\n", fname);
875 
876  mpi_printf("%s file: '%s' (file 1 of %d)\n", info, fname, numfilesperdump);
877 
878  if(file_format == FILEFORMAT_LEGACY2)
879  {
880  blksize = sizeof(int) + 4 * sizeof(char);
881  SKIP;
882  my_fwrite((const void *)"HEAD", sizeof(char), 4, fd);
883  int nextblock = header_size + 2 * sizeof(int);
884  my_fwrite(&nextblock, sizeof(int), 1, fd);
885  SKIP;
886  }
887 
888  blksize = header_size;
889  SKIP;
891  SKIP;
892  }
893  }
894 
895  for(int blocknr = 0; blocknr < N_IO_Fields; blocknr++)
896  {
897  polling(numfilesperdump);
898 
899  if(IO_Fields[blocknr].type_in_file_output != FILE_NONE)
900  {
901  unsigned int bytes_per_blockelement = get_bytes_per_memory_blockelement(blocknr, 0);
902  int blockmaxlen = (int)(COMMBUFFERSIZE / bytes_per_blockelement);
903  long long npart_in_block = get_particles_in_block(blocknr, npart, typelist);
904  hid_t hdf5_memory_datatype = get_hdf5_memorytype_of_block(blocknr);
905  char dname[MAXLEN_PATH];
906  get_dataset_name(blocknr, dname);
907 
908  if(npart_in_block > 0)
909  {
910  mpi_printf("%s block %d (%s)...\n", info, blocknr, dname);
911 
912  if(ThisTask == writeTask)
913  {
914  if(file_format == FILEFORMAT_LEGACY1 || file_format == FILEFORMAT_LEGACY2)
915  {
916  bytes_per_blockelement_in_file =
917  IO_Fields[blocknr].values_per_block * H5Tget_size(get_hdf5_outputtype_of_block(blocknr));
918 
919  if(file_format == FILEFORMAT_LEGACY2)
920  {
921  blksize = sizeof(int) + LABEL_LEN * sizeof(char);
922  SKIP;
923  get_Tab_IO_Label(blocknr, label);
924  my_fwrite(label, sizeof(char), LABEL_LEN, fd);
925  int nextblock = npart_in_block * bytes_per_blockelement_in_file + 2 * sizeof(int);
926  my_fwrite(&nextblock, sizeof(int), 1, fd);
927  SKIP;
928  }
929 
930  blksize = npart_in_block * bytes_per_blockelement_in_file;
931  SKIP;
932  }
933  }
934 
935  for(int type = 0; type < N_DataGroups; type++)
936  {
937  if(typelist[type])
938  {
939  if(ThisTask == writeTask && file_format == FILEFORMAT_HDF5 && npart[type] > 0)
940  {
941  hid_t hdf5_file_datatype = get_hdf5_outputtype_of_block(blocknr);
942 
943  dims[0] = npart[type];
944  dims[1] = get_values_per_blockelement(blocknr);
945  if(dims[1] == 1)
946  rank = 1;
947  else
948  rank = 2;
949 
950  if(chunksize > 0)
951  {
952  hsize_t maxdims[2] = {H5S_UNLIMITED, dims[1]};
953  hdf5_dataspace_in_file = my_H5Screate_simple(rank, dims, maxdims);
954 
955  /* Modify dataset creation properties, i.e. enable chunking */
956  hdf5_prop = H5Pcreate(H5P_DATASET_CREATE);
957  hsize_t chunk_dims[2] = {0, dims[1]};
958  chunk_dims[0] = chunksize;
959  H5Pset_chunk(hdf5_prop, rank, chunk_dims);
960 
961  hdf5_dataset =
962  my_H5Dcreate(hdf5_grp[type], dname, hdf5_file_datatype, hdf5_dataspace_in_file, hdf5_prop);
963  }
964  else if(IO_Fields[blocknr].compression_on)
965  {
966  hdf5_dataspace_in_file = my_H5Screate_simple(rank, dims, NULL);
967 
968  /* Modify dataset creation properties, i.e. enable compression */
969  hdf5_prop = H5Pcreate(H5P_DATASET_CREATE);
970  hsize_t chunk_dims[2] = {COMPRESSION_CHUNKSIZE, dims[1]};
971  if(chunk_dims[0] > dims[0])
972  chunk_dims[0] = dims[0];
973  H5Pset_chunk(hdf5_prop, rank, chunk_dims); /* set chunk size */
974  H5Pset_shuffle(hdf5_prop); /* reshuffle bytes to get better compression ratio */
975  H5Pset_deflate(hdf5_prop, 9); /* gzip compression level 9 */
976  if(H5Pall_filters_avail(hdf5_prop))
977  hdf5_dataset =
978  my_H5Dcreate(hdf5_grp[type], dname, hdf5_file_datatype, hdf5_dataspace_in_file, hdf5_prop);
979  else
980  Terminate("HDF5: Compression not available!\n");
981  }
982  else
983  {
984  hdf5_dataspace_in_file = my_H5Screate_simple(rank, dims, NULL);
985  hdf5_dataset =
986  my_H5Dcreate(hdf5_grp[type], dname, hdf5_file_datatype, hdf5_dataspace_in_file, H5P_DEFAULT);
987  }
988 
989  write_dataset_attributes(hdf5_dataset, blocknr);
990 
991  byte_count += dims[0] * dims[1] * my_H5Tget_size(hdf5_file_datatype); /* for I/O performance measurement */
992 
993  pcsum = 0;
994  }
995 
996  for(int task = writeTask, offset = 0; task <= lastTask; task++)
997  {
998  int n_for_this_task;
999 
1000  if(task == ThisTask)
1001  {
1002  n_for_this_task = n_type[type];
1003 
1004  for(int p = writeTask; p <= lastTask; p++)
1005  if(p != ThisTask)
1006  MPI_Send(&n_for_this_task, 1, MPI_INT, p, TAG_NFORTHISTASK, Communicator);
1007  }
1008  else
1009  MPI_Recv(&n_for_this_task, 1, MPI_INT, task, TAG_NFORTHISTASK, Communicator, MPI_STATUS_IGNORE);
1010 
1011  while(n_for_this_task > 0)
1012  {
1013  int pc = n_for_this_task;
1014 
1015  if(pc > blockmaxlen)
1016  pc = blockmaxlen;
1017 
1018  if(ThisTask == task)
1019  fill_write_buffer(blocknr, &offset, pc, type, CommBuffer);
1020 
1021  if(ThisTask == writeTask && task != writeTask)
1022  MPI_Recv(CommBuffer, bytes_per_blockelement * pc, MPI_BYTE, task, TAG_PDATA, Communicator,
1023  MPI_STATUS_IGNORE);
1024 
1025  if(ThisTask != writeTask && task == ThisTask)
1026  MPI_Ssend(CommBuffer, bytes_per_blockelement * pc, MPI_BYTE, writeTask, TAG_PDATA, Communicator);
1027 
1028  if(ThisTask == writeTask)
1029  {
1030  if(file_format == FILEFORMAT_HDF5)
1031  {
1032  start[0] = pcsum;
1033  start[1] = 0;
1034 
1035  count[0] = pc;
1036  count[1] = get_values_per_blockelement(blocknr);
1037  pcsum += pc;
1038 
1039  my_H5Sselect_hyperslab(hdf5_dataspace_in_file, H5S_SELECT_SET, start, NULL, count, NULL);
1040 
1041  dims[0] = pc;
1042  dims[1] = get_values_per_blockelement(blocknr);
1043  hdf5_dataspace_memory = my_H5Screate_simple(rank, dims, NULL);
1044 
1045  my_H5Dwrite(hdf5_dataset, hdf5_memory_datatype, hdf5_dataspace_memory, hdf5_dataspace_in_file,
1046  H5P_DEFAULT, CommBuffer, dname);
1047 
1048  my_H5Sclose(hdf5_dataspace_memory, H5S_SIMPLE);
1049  }
1050  else
1051  {
1052  if(bytes_per_blockelement_in_file != bytes_per_blockelement)
1053  {
1054  char *CommAuxBuffer =
1055  (char *)Mem.mymalloc("CommAuxBuffer", bytes_per_blockelement_in_file * pc);
1056  type_cast_data((char *)CommBuffer, bytes_per_blockelement, (char *)CommAuxBuffer,
1057  bytes_per_blockelement_in_file, pc, blocknr);
1058  my_fwrite(CommAuxBuffer, bytes_per_blockelement_in_file, pc, fd);
1059  Mem.myfree(CommAuxBuffer);
1060  }
1061  else
1062  my_fwrite(CommBuffer, bytes_per_blockelement, pc, fd);
1063  }
1064  }
1065 
1066  n_for_this_task -= pc;
1067  }
1068  }
1069 
1070  if(ThisTask == writeTask && file_format == FILEFORMAT_HDF5 && npart[type] > 0)
1071  {
1072  if(file_format == FILEFORMAT_HDF5)
1073  {
1074  my_H5Dclose(hdf5_dataset, dname);
1075  if(chunksize > 0 || IO_Fields[blocknr].compression_on)
1076  my_H5Pclose(hdf5_prop);
1077  my_H5Sclose(hdf5_dataspace_in_file, H5S_SIMPLE);
1078  }
1079  }
1080  }
1081  }
1082 
1083  if(ThisTask == writeTask)
1084  {
1085  if(file_format == FILEFORMAT_LEGACY1 || file_format == FILEFORMAT_LEGACY2)
1086  SKIP;
1087  }
1088  }
1089  }
1090  }
1091 
1092  if(ThisTask == writeTask)
1093  {
1094  if(file_format == FILEFORMAT_HDF5)
1095  {
1096  char buf[MAXLEN_PATH];
1097 
1098  for(int type = N_DataGroups - 1; type >= 0; type--)
1099  if(npart[type] > 0)
1100  {
1101  get_datagroup_name(type, buf);
1102  my_H5Gclose(hdf5_grp[type], buf);
1103  }
1104 
1105  my_H5Gclose(hdf5_headergrp, "/Header");
1106  my_H5Gclose(hdf5_paramsgrp, "/Parameters");
1107  my_H5Gclose(hdf5_configgrp, "/Config");
1108 
1109  sprintf(buf, "%s.hdf5", fname);
1110  my_H5Fclose(hdf5_file, buf);
1111  }
1112  else
1113  fclose(fd);
1114  }
1115 }
1116 
1117 void IO_Def::append_file(char *fname, int writeTask, int lastTask, void *CommBuffer, int numfilesperdump, int chunksize)
1118 {
1119  int typelist[N_DataGroups];
1120  long long n_type[N_DataGroups], npart[N_DataGroups], n_previous[N_DataGroups];
1121  hid_t hdf5_file = 0, hdf5_grp[N_DataGroups], hdf5_headergrp = 0, hdf5_dataspace_memory;
1122  hid_t hdf5_dataspace_in_file = 0, hdf5_dataset = 0;
1123  hsize_t dims[2], count[2], start[2];
1124  int rank = 0, pcsum = 0;
1125 
1126  if(file_format != FILEFORMAT_HDF5)
1127  Terminate("appending to files only works with HDF5 format\n");
1128 
1129  read_file_header(NULL, 0, writeTask, lastTask, n_type, npart, NULL);
1130 
1131  for(int n = 0; n < N_DataGroups; n++)
1132  n_previous[n] = n_type[n];
1133 
1134  fill_file_header(writeTask, lastTask, n_type, npart);
1135 
1136  /* open file and write header */
1137  if(ThisTask == writeTask)
1138  {
1139  char buf[MAXLEN_PATH];
1140  sprintf(buf, "%s.hdf5", fname);
1141 
1142  hdf5_file = my_H5Fopen(buf, H5F_ACC_RDWR, H5P_DEFAULT);
1143 
1144  hdf5_headergrp = my_H5Gopen(hdf5_file, "/Header");
1145  write_header_fields(hdf5_headergrp);
1146 
1147  for(int type = 0; type < N_DataGroups; type++)
1148  {
1149  if(npart[type] > 0)
1150  {
1151  get_datagroup_name(type, buf);
1152  if(n_previous[type] == 0)
1153  hdf5_grp[type] = my_H5Gcreate(hdf5_file, buf, 0);
1154  else
1155  hdf5_grp[type] = my_H5Gopen(hdf5_file, buf);
1156  }
1157  }
1158  }
1159 
1160  for(int blocknr = 0; blocknr < N_IO_Fields; blocknr++)
1161  {
1162  polling(numfilesperdump);
1163 
1164  if(IO_Fields[blocknr].type_in_file_output != FILE_NONE)
1165  {
1166  unsigned int bytes_per_blockelement = get_bytes_per_memory_blockelement(blocknr, 0);
1167  int blockmaxlen = (int)(COMMBUFFERSIZE / bytes_per_blockelement);
1168  long long npart_in_block = get_particles_in_block(blocknr, npart, typelist);
1169  hid_t hdf5_memory_datatype = get_hdf5_memorytype_of_block(blocknr);
1170  char dname[MAXLEN_PATH];
1171  get_dataset_name(blocknr, dname);
1172 
1173  if(npart_in_block > 0)
1174  {
1175  mpi_printf("%s block %d (%s)...\n", info, blocknr, dname);
1176 
1177  for(int type = 0; type < N_DataGroups; type++)
1178  {
1179  if(typelist[type])
1180  {
1181  if(ThisTask == writeTask && file_format == FILEFORMAT_HDF5 && npart[type] > 0)
1182  {
1183  hid_t hdf5_file_datatype = get_hdf5_outputtype_of_block(blocknr);
1184 
1185  dims[0] = npart[type] + n_previous[type];
1186  dims[1] = get_values_per_blockelement(blocknr);
1187  if(dims[1] == 1)
1188  rank = 1;
1189  else
1190  rank = 2;
1191 
1192  hdf5_dataspace_in_file = my_H5Screate_simple(rank, dims, NULL);
1193 
1194  if(n_previous[type] == 0)
1195  {
1196  if(chunksize > 0)
1197  {
1198  hsize_t maxdims[2] = {H5S_UNLIMITED, dims[1]};
1199  hdf5_dataspace_in_file = my_H5Screate_simple(rank, dims, maxdims);
1200 
1201  /* Modify dataset creation properties, i.e. enable chunking */
1202  hid_t prop = H5Pcreate(H5P_DATASET_CREATE);
1203  hsize_t chunk_dims[2] = {0, dims[1]};
1204  chunk_dims[0] = chunksize;
1205  H5Pset_chunk(prop, rank, chunk_dims);
1206 
1207  hdf5_dataset = my_H5Dcreate(hdf5_grp[type], dname, hdf5_file_datatype, hdf5_dataspace_in_file, prop);
1208  }
1209  else
1210  {
1211  hdf5_dataset =
1212  my_H5Dcreate(hdf5_grp[type], dname, hdf5_file_datatype, hdf5_dataspace_in_file, H5P_DEFAULT);
1213  write_dataset_attributes(hdf5_dataset, blocknr);
1214  }
1215  }
1216  else
1217  {
1218  hdf5_dataset = my_H5Dopen_if_existing(hdf5_grp[type], dname);
1219  my_H5Dset_extent(hdf5_dataset, dims);
1220  }
1221 
1222  byte_count += dims[0] * dims[1] * my_H5Tget_size(hdf5_file_datatype); /* for I/O performance measurement */
1223 
1224  pcsum = 0;
1225  }
1226 
1227  for(int task = writeTask, offset = 0; task <= lastTask; task++)
1228  {
1229  int n_for_this_task;
1230 
1231  if(task == ThisTask)
1232  {
1233  n_for_this_task = n_type[type];
1234 
1235  for(int p = writeTask; p <= lastTask; p++)
1236  if(p != ThisTask)
1237  MPI_Send(&n_for_this_task, 1, MPI_INT, p, TAG_NFORTHISTASK, Communicator);
1238  }
1239  else
1240  MPI_Recv(&n_for_this_task, 1, MPI_INT, task, TAG_NFORTHISTASK, Communicator, MPI_STATUS_IGNORE);
1241 
1242  while(n_for_this_task > 0)
1243  {
1244  int pc = n_for_this_task;
1245 
1246  if(pc > blockmaxlen)
1247  pc = blockmaxlen;
1248 
1249  if(ThisTask == task)
1250  fill_write_buffer(blocknr, &offset, pc, type, CommBuffer);
1251 
1252  if(ThisTask == writeTask && task != writeTask)
1253  MPI_Recv(CommBuffer, bytes_per_blockelement * pc, MPI_BYTE, task, TAG_PDATA, Communicator,
1254  MPI_STATUS_IGNORE);
1255 
1256  if(ThisTask != writeTask && task == ThisTask)
1257  MPI_Ssend(CommBuffer, bytes_per_blockelement * pc, MPI_BYTE, writeTask, TAG_PDATA, Communicator);
1258 
1259  if(ThisTask == writeTask)
1260  {
1261  start[0] = pcsum + n_previous[type];
1262  start[1] = 0;
1263 
1264  count[0] = pc;
1265  count[1] = get_values_per_blockelement(blocknr);
1266  pcsum += pc;
1267 
1268  my_H5Sselect_hyperslab(hdf5_dataspace_in_file, H5S_SELECT_SET, start, NULL, count, NULL);
1269 
1270  dims[0] = pc;
1271  dims[1] = get_values_per_blockelement(blocknr);
1272  hdf5_dataspace_memory = my_H5Screate_simple(rank, dims, NULL);
1273 
1274  my_H5Dwrite(hdf5_dataset, hdf5_memory_datatype, hdf5_dataspace_memory, hdf5_dataspace_in_file,
1275  H5P_DEFAULT, CommBuffer, dname);
1276 
1277  my_H5Sclose(hdf5_dataspace_memory, H5S_SIMPLE);
1278  }
1279 
1280  n_for_this_task -= pc;
1281  }
1282  }
1283 
1284  if(ThisTask == writeTask && npart[type] > 0)
1285  {
1286  my_H5Dclose(hdf5_dataset, dname);
1287  my_H5Sclose(hdf5_dataspace_in_file, H5S_SIMPLE);
1288  }
1289  }
1290  }
1291  }
1292  }
1293  }
1294 
1295  if(ThisTask == writeTask)
1296  {
1297  char buf[MAXLEN_PATH];
1298 
1299  for(int type = N_DataGroups - 1; type >= 0; type--)
1300  if(npart[type] > 0)
1301  {
1302  get_datagroup_name(type, buf);
1303  my_H5Gclose(hdf5_grp[type], buf);
1304  }
1305 
1306  my_H5Gclose(hdf5_headergrp, "/Header");
1307 
1308  sprintf(buf, "%s.hdf5", fname);
1309  my_H5Fclose(hdf5_file, buf);
1310  }
1311 }
1312 
1313 void IO_Def::type_cast_data(char *src, int src_bytes_per_element, char *target, int target_bytes_per_element, int len, int blocknr)
1314 {
1315  switch(IO_Fields[blocknr].type_in_memory)
1316  {
1317  case MEM_INT:
1318  case MEM_INT64:
1319  case MEM_MY_ID_TYPE:
1320  case MEM_MY_INTPOS_TYPE:
1321  case MEM_MY_FILEOFFSET:
1322  if(target_bytes_per_element > src_bytes_per_element)
1323  {
1324  if(target_bytes_per_element != 2 * src_bytes_per_element)
1325  Terminate("something is odd here: target_bytes_per_element=%d != 2 * src_bytes_per_element=%d", target_bytes_per_element,
1326  src_bytes_per_element);
1327 
1328  int fac = src_bytes_per_element / sizeof(int);
1329  int *sp = (int *)src;
1330  long long *tp = (long long *)target;
1331  for(int i = 0; i < len * fac; i++)
1332  *tp++ = *sp++;
1333  }
1334  else
1335  {
1336  if(src_bytes_per_element != 2 * target_bytes_per_element)
1337  Terminate("something is odd here: src_bytes_per_element=%d != 2 * target_bytes_per_element=%d", src_bytes_per_element,
1338  target_bytes_per_element);
1339  int fac = src_bytes_per_element / sizeof(long long);
1340  long long *sp = (long long *)src;
1341  int *tp = (int *)target;
1342  for(int i = 0; i < len * fac; i++)
1343  *tp++ = *sp++;
1344  }
1345  break;
1346 
1347  case MEM_FLOAT:
1348  case MEM_DOUBLE:
1349  case MEM_MY_FLOAT:
1350  case MEM_MY_DOUBLE:
1351  if((target_bytes_per_element % 8) == 0 &&
1352  (src_bytes_per_element % 4) == 0) // target_bytes_per_element multiply of 8, src_bytes_per_element mulitple of 4
1353  {
1354  int fac = src_bytes_per_element / sizeof(float);
1355  float *sp = (float *)src;
1356  double *tp = (double *)target;
1357  for(int i = 0; i < len * fac; i++)
1358  *tp++ = *sp++;
1359  }
1360  else if((target_bytes_per_element % 4) == 0 && (src_bytes_per_element % 8) == 0)
1361  {
1362  int fac = src_bytes_per_element / sizeof(double);
1363  double *sp = (double *)src;
1364  float *tp = (float *)target;
1365  for(int i = 0; i < len * fac; i++)
1366  *tp++ = *sp++;
1367  }
1368  else if((target_bytes_per_element & 8) == 0 && (src_bytes_per_element % 2) == 0)
1369  {
1370  int fac = src_bytes_per_element / sizeof(half);
1371  half *sp = (half *)src;
1372  double *tp = (double *)target;
1373  for(int i = 0; i < len * fac; i++)
1374  *tp++ = *sp++;
1375  }
1376  else if((target_bytes_per_element % 2) == 0 && (src_bytes_per_element % 8) == 0)
1377  {
1378  int fac = src_bytes_per_element / sizeof(double);
1379  double *sp = (double *)src;
1380  half *tp = (half *)target;
1381  for(int i = 0; i < len * fac; i++)
1382  *tp++ = *sp++;
1383  }
1384  else
1385  {
1386  Terminate("Strange conversion requested: target_bytes_per_element=%d src_bytes_per_element=%d\n",
1387  target_bytes_per_element, src_bytes_per_element);
1388  }
1389  break;
1390  }
1391 }
1392 
1402 void IO_Def::read_file(const char *fname, int filenr, int readTask, int lastTask, void *CommBuffer)
1403 {
1404  long long n_type[N_DataGroups], npart[N_DataGroups];
1405  int typelist[N_DataGroups];
1406  unsigned int blksize1, blksize2, bytes_per_blockelement_in_file = 0;
1407  hid_t hdf5_file = 0, hdf5_grp[N_DataGroups], hdf5_dataspace_in_file;
1408  hid_t hdf5_dataspace_in_memory, hdf5_dataset;
1409  FILE *fd = 0;
1410 
1411  /* open file and read header */
1412 
1413  if(ThisTask == readTask)
1414  {
1415  if(file_format == FILEFORMAT_HDF5)
1416  {
1417  read_header_fields(fname);
1418  }
1419  else if(file_format == FILEFORMAT_LEGACY1 || file_format == FILEFORMAT_LEGACY2)
1420  {
1421  if(!(fd = fopen(fname, "r")))
1422  Terminate("can't open file `%s' for reading initial conditions.\n", fname);
1423 
1424  if(file_format == FILEFORMAT_LEGACY2)
1425  {
1426  int nextblock;
1427  char label[LABEL_LEN];
1428  my_fread(&blksize1, sizeof(int), 1, fd);
1429  my_fread(&label, sizeof(char), LABEL_LEN, fd);
1430  my_fread(&nextblock, sizeof(int), 1, fd);
1431  my_fread(&blksize2, sizeof(int), 1, fd);
1432  }
1433 
1434  my_fread(&blksize1, sizeof(int), 1, fd);
1435  my_fread(header_buf, header_size, 1, fd);
1436  my_fread(&blksize2, sizeof(int), 1, fd);
1437  }
1438  else
1439  Terminate("unknown ICFormat");
1440 
1441  for(int task = readTask + 1; task <= lastTask; task++)
1442  MPI_Ssend(header_buf, header_size, MPI_BYTE, task, TAG_HEADER, Communicator);
1443  }
1444  else
1445  MPI_Recv(header_buf, header_size, MPI_BYTE, readTask, TAG_HEADER, Communicator, MPI_STATUS_IGNORE);
1446 
1447  int nstart;
1448  read_file_header(fname, filenr, readTask, lastTask, n_type, npart, &nstart);
1449 
1450  if(ThisTask == readTask)
1451  {
1452  if(file_format == FILEFORMAT_HDF5)
1453  {
1454  hdf5_file = my_H5Fopen(fname, H5F_ACC_RDONLY, H5P_DEFAULT);
1455 
1456  for(int type = 0; type < N_DataGroups; type++)
1457  {
1458  if(npart[type] > 0)
1459  {
1460  char buf[MAXLEN_PATH];
1461  get_datagroup_name(type, buf);
1462  hdf5_grp[type] = my_H5Gopen(hdf5_file, buf);
1463  }
1464  }
1465  }
1466  }
1467 
1468  for(int blocknr = 0; blocknr < N_IO_Fields; blocknr++)
1469  {
1470  if((IO_Fields[blocknr].read_flag != SKIP_ON_READ &&
1472  blocknr > 4) /* this second conditions allows short legacy ICs to be read in */
1473  ) ||
1474  IO_Fields[blocknr].type_in_memory == MEM_MY_FILEOFFSET)
1475  {
1476  unsigned int bytes_per_blockelement = get_bytes_per_memory_blockelement(blocknr, 1);
1477  int blockmaxlen = (int)(COMMBUFFERSIZE / bytes_per_blockelement);
1478  long long npart_in_block = get_particles_in_block(blocknr, npart, &typelist[0]);
1479  hid_t hdf5_memory_datatype = get_hdf5_memorytype_of_block(blocknr);
1480  char dname[MAXLEN_PATH];
1481  get_dataset_name(blocknr, dname);
1482 
1483  if(npart_in_block > 0)
1484  {
1485  if(filenr == 0)
1486  mpi_printf("READIC: reading block %d (%s)...\n", blocknr, dname);
1487 
1488  if(ThisTask == readTask && IO_Fields[blocknr].type_in_memory != MEM_MY_FILEOFFSET)
1489  {
1490  if(file_format == FILEFORMAT_LEGACY2)
1491  {
1492  char expected_label[LABEL_LEN + 1];
1493  get_Tab_IO_Label(blocknr, expected_label);
1494  find_block(expected_label, fd);
1495  }
1496 
1497  if(file_format == FILEFORMAT_LEGACY1 || file_format == FILEFORMAT_LEGACY2)
1498  {
1499  my_fread(&blksize1, sizeof(int), 1, fd);
1500  bytes_per_blockelement_in_file = blksize1 / npart_in_block;
1501  }
1502  }
1503 
1504  int offset = 0;
1505 
1506  for(int type = 0; type < N_DataGroups; type++)
1507  {
1509  offset = 0;
1510 
1511  int n_in_file = npart[type];
1512  int pcsum = 0;
1513 
1514  long long nprevious = 0;
1515  for(int t = 0; t < type; t++)
1516  for(int f = 0; f < get_filenr_from_header(); f++)
1517  nprevious += ntype_in_files[f * N_DataGroups + t];
1518 
1519  for(int nr = 0; nr < filenr; nr++)
1520  nprevious += ntype_in_files[nr * N_DataGroups + type];
1521 
1522  if(typelist[type] == 0)
1523  {
1524  int ntask = lastTask - readTask + 1;
1525  int n_for_this_task = n_in_file / ntask;
1526  if((ThisTask - readTask) < (n_in_file % ntask))
1527  n_for_this_task++;
1528 
1529  offset += n_for_this_task;
1530  }
1531  else
1532  {
1533  for(int task = readTask; task <= lastTask; task++)
1534  {
1535  int ntask = lastTask - readTask + 1;
1536  int n_for_this_task = n_in_file / ntask;
1537  if((task - readTask) < (n_in_file % ntask))
1538  n_for_this_task++;
1539 
1540  do
1541  {
1542  int pc = n_for_this_task;
1543 
1544  if(pc > blockmaxlen)
1545  pc = blockmaxlen;
1546 
1547  if(ThisTask == readTask && IO_Fields[blocknr].type_in_memory != MEM_MY_FILEOFFSET)
1548  {
1549  if(file_format == FILEFORMAT_LEGACY1 || file_format == FILEFORMAT_LEGACY2)
1550  {
1551  if(bytes_per_blockelement_in_file != bytes_per_blockelement)
1552  {
1553  char *CommAuxBuffer =
1554  (char *)Mem.mymalloc("CommAuxBuffer", bytes_per_blockelement_in_file * pc);
1555  my_fread(CommAuxBuffer, bytes_per_blockelement_in_file, pc, fd);
1556  type_cast_data((char *)CommAuxBuffer, bytes_per_blockelement_in_file, (char *)CommBuffer,
1557  bytes_per_blockelement, pc, blocknr);
1558  Mem.myfree(CommAuxBuffer);
1559  }
1560  else
1561  my_fread(CommBuffer, bytes_per_blockelement, pc, fd);
1562  }
1563 
1564  if(file_format == FILEFORMAT_HDF5 && pc > 0)
1565  {
1566  hdf5_dataset = my_H5Dopen_if_existing(hdf5_grp[type], dname);
1567  int rank;
1568  hsize_t dims[2], count[2], start[2];
1569 
1570  dims[0] = npart[type];
1571  dims[1] = get_values_per_blockelement(blocknr);
1572  if(dims[1] == 1)
1573  rank = 1;
1574  else
1575  rank = 2;
1576 
1577  hdf5_dataspace_in_file = my_H5Screate_simple(rank, dims, NULL);
1578 
1579  dims[0] = pc;
1580  hdf5_dataspace_in_memory = my_H5Screate_simple(rank, dims, NULL);
1581 
1582  start[0] = pcsum;
1583  start[1] = 0;
1584 
1585  count[0] = pc;
1586  count[1] = get_values_per_blockelement(blocknr);
1587  pcsum += pc;
1588 
1589  my_H5Sselect_hyperslab(hdf5_dataspace_in_file, H5S_SELECT_SET, start, NULL, count, NULL);
1590 
1591  // Test if dataset was present
1592  if(hdf5_dataset < 0)
1593  {
1594  // no, pad with zeros
1595  if((ThisTask == readTask) && (task == ThisTask))
1596  printf("\tDataset %s not present for particle type %d, using zero.\n", dname, type);
1597  memset(CommBuffer, 0, dims[0] * dims[1] * my_H5Tget_size(hdf5_memory_datatype));
1598  }
1599  else
1600  {
1601  hid_t hdf5_file_datatype = H5Dget_type(hdf5_dataset);
1602  byte_count += dims[0] * dims[1] *
1603  my_H5Tget_size(hdf5_file_datatype); /* for I/O performance measurement */
1604  H5Tclose(hdf5_file_datatype);
1605 
1606  my_H5Dread(hdf5_dataset, hdf5_memory_datatype, hdf5_dataspace_in_memory,
1607  hdf5_dataspace_in_file, H5P_DEFAULT, CommBuffer, dname);
1608  my_H5Dclose(hdf5_dataset, dname);
1609  }
1610  my_H5Sclose(hdf5_dataspace_in_memory, H5S_SIMPLE);
1611  my_H5Sclose(hdf5_dataspace_in_file, H5S_SIMPLE);
1612  }
1613  }
1614 
1615  if(ThisTask == readTask && task != readTask && pc > 0)
1616  MPI_Ssend(CommBuffer, bytes_per_blockelement * pc, MPI_BYTE, task, TAG_PDATA, Communicator);
1617 
1618  if(ThisTask != readTask && task == ThisTask && pc > 0)
1619  MPI_Recv(CommBuffer, bytes_per_blockelement * pc, MPI_BYTE, readTask, TAG_PDATA, Communicator,
1620  MPI_STATUS_IGNORE);
1621 
1622  if(ThisTask == task)
1623  {
1624  if(blocknr == 0 && IO_Fields[blocknr].array == A_P)
1625  {
1626  for(int n = 0; n < pc; n++)
1627  set_type_of_element(nstart + offset + n, type); /* initialize type */
1628  }
1629 #ifdef MERGERTREE
1630  if(blocknr == 0 && IO_Fields[blocknr].array == A_MTRP)
1631  {
1632  for(int n = 0; n < pc; n++)
1633  {
1634  set_type_of_element(nstart + offset + n, type); /* initialize type */
1635  // MtrP[nstart + offset + n].Type = type; /* initialize type */
1636  }
1637  }
1638 #endif
1639  empty_read_buffer(blocknr, nstart + offset, pc, type, nprevious, CommBuffer);
1640 
1641  offset += pc;
1642  }
1643 
1644  n_for_this_task -= pc;
1645  nprevious += pc;
1646  }
1647  while(n_for_this_task > 0);
1648  }
1649  }
1650  }
1651  if(ThisTask == readTask && IO_Fields[blocknr].type_in_memory != MEM_MY_FILEOFFSET)
1652  {
1653  if(file_format == FILEFORMAT_LEGACY1 || file_format == FILEFORMAT_LEGACY2)
1654  {
1655  my_fread(&blksize2, sizeof(int), 1, fd);
1656  ;
1657  if(blksize1 != blksize2)
1658  {
1659  char buf[MAXLEN_PATH];
1660  sprintf(buf, "incorrect block-sizes detected!\n Task=%d blocknr=%d blksize1=%d blksize2=%d\n", ThisTask,
1661  blocknr, blksize1, blksize2);
1662  if(blocknr == 2) /* block number 2 is always IDs */
1663  strcat(buf, "Possible mismatch of 32bit and 64bit ID's in IC file and GADGET compilation !\n");
1664  Terminate(buf);
1665  }
1666  }
1667  }
1668  }
1669  }
1670  }
1671 
1672  for(int type = 0; type < N_DataGroups; type++)
1673  {
1674  long long n_in_file = npart[type];
1675  int ntask = lastTask - readTask + 1;
1676  int n_for_this_task = n_in_file / ntask;
1677  if((ThisTask - readTask) < (n_in_file % ntask))
1678  n_for_this_task++;
1679 
1680  read_increase_numbers(type, n_for_this_task);
1681  }
1682 
1683  if(ThisTask == readTask)
1684  {
1685  if(file_format == FILEFORMAT_LEGACY1 || file_format == FILEFORMAT_LEGACY2)
1686  fclose(fd);
1687 
1688  if(file_format == FILEFORMAT_HDF5)
1689  {
1690  for(int type = N_DataGroups - 1; type >= 0; type--)
1691  if(npart[type] > 0)
1692  {
1693  char buf[MAXLEN_PATH];
1694  get_datagroup_name(type, buf);
1695  my_H5Gclose(hdf5_grp[type], buf);
1696  }
1697  my_H5Fclose(hdf5_file, fname);
1698  }
1699  }
1700 }
1701 
1713 void IO_Def::distribute_file(int nfiles, int *filenr, int *master, int *last)
1714 {
1715  int tasks_per_file = NTask / nfiles;
1716  int tasks_left = NTask % nfiles;
1717 
1718  if(tasks_left == 0)
1719  {
1720  int group = ThisTask / tasks_per_file;
1721  *master = group * tasks_per_file;
1722  *last = (group + 1) * tasks_per_file - 1;
1723  *filenr = group;
1724  return;
1725  }
1726 
1727  double tpf = ((double)NTask) / nfiles;
1728 
1729  *last = -1;
1730  for(int i = 0; i < nfiles; i++)
1731  {
1732  *master = *last + 1;
1733  *last = (i + 1) * tpf;
1734  if(*last >= NTask)
1735  *last = *last - 1;
1736  if(*last < *master)
1737  Terminate("last < master");
1738  *filenr = i;
1739 
1740  if(i == nfiles - 1)
1741  *last = NTask - 1;
1742 
1743  if(ThisTask >= *master && ThisTask <= *last)
1744  return;
1745  }
1746 }
1747 
1757 int IO_Def::get_bytes_per_memory_blockelement(int blocknr, int mode)
1758 {
1759  if(blocknr < 0 || blocknr >= N_IO_Fields)
1760  Terminate("something is wrong here: blocknr=%d N_IO_Fields=%d", blocknr, N_IO_Fields);
1761 
1762  int bytes_per_blockelement = 0;
1763 
1764  IO_Field *field = &IO_Fields[blocknr];
1765 
1766  switch(field->type_in_memory)
1767  {
1768  case MEM_INT:
1769  bytes_per_blockelement = field->values_per_block * sizeof(int);
1770  break;
1771  case MEM_INT64:
1772  bytes_per_blockelement = field->values_per_block * sizeof(long long);
1773  break;
1774  case MEM_MY_ID_TYPE:
1775  bytes_per_blockelement = field->values_per_block * sizeof(MyIDType);
1776  break;
1777  case MEM_MY_INTPOS_TYPE:
1778  bytes_per_blockelement = field->values_per_block * sizeof(MyIntPosType);
1779  break;
1780  case MEM_FLOAT:
1781  bytes_per_blockelement = field->values_per_block * sizeof(float);
1782  break;
1783  case MEM_DOUBLE:
1784  bytes_per_blockelement = field->values_per_block * sizeof(double);
1785  break;
1786  case MEM_MY_FLOAT:
1787  bytes_per_blockelement = field->values_per_block * sizeof(MyFloat);
1788  break;
1789  case MEM_MY_DOUBLE:
1790  bytes_per_blockelement = field->values_per_block * sizeof(MyDouble);
1791  break;
1792  case MEM_MY_FILEOFFSET:
1793  bytes_per_blockelement = field->values_per_block * sizeof(long long);
1794  break;
1795  }
1796 
1797  return bytes_per_blockelement;
1798 }
1799 
1800 hid_t IO_Def::get_hdf5_outputtype_of_block(int blocknr)
1801 {
1802  hid_t hdf5_datatype = 0;
1803 
1804  switch(IO_Fields[blocknr].type_in_file_output)
1805  {
1806  case FILE_INT:
1807  hdf5_datatype = H5T_NATIVE_INT;
1808  break;
1809  case FILE_INT64:
1810  hdf5_datatype = H5T_NATIVE_INT64;
1811  break;
1812  case FILE_MY_IO_FLOAT:
1813 #ifdef OUTPUT_IN_DOUBLEPRECISION
1814  hdf5_datatype = H5T_NATIVE_DOUBLE;
1815 #else
1816  hdf5_datatype = H5T_NATIVE_FLOAT;
1817 #endif
1818  break;
1819  case FILE_MY_ID_TYPE:
1820 #if defined(IDS_32BIT)
1821  hdf5_datatype = H5T_NATIVE_UINT32;
1822 #elif defined(IDS_48BIT)
1823  hdf5_datatype = Int48_memtype;
1824 #else
1825  hdf5_datatype = H5T_NATIVE_UINT64;
1826 #endif
1827  break;
1828  case FILE_MY_INTPOS_TYPE:
1829 #if defined(POSITIONS_IN_32BIT)
1830  hdf5_datatype = H5T_NATIVE_UINT32;
1831 #elif defined(POSITIONS_IN_64BIT)
1832  hdf5_datatype = H5T_NATIVE_UINT64;
1833 #else
1834  hdf5_datatype = Int128_memtype;
1835 #endif
1836  break;
1837  case FILE_DOUBLE:
1838  hdf5_datatype = H5T_NATIVE_DOUBLE;
1839  break;
1840  case FILE_FLOAT:
1841  hdf5_datatype = H5T_NATIVE_FLOAT;
1842  break;
1843  case FILE_HALF:
1844  hdf5_datatype = Halfprec_memtype;
1845  break;
1846  case FILE_NONE:
1847  Terminate("undefined type");
1848  break;
1849  }
1850 
1851  return hdf5_datatype;
1852 }
1853 
1854 hid_t IO_Def::get_hdf5_memorytype_of_block(int blocknr)
1855 {
1856  hid_t hdf5_datatype = 0;
1857 
1858  switch(IO_Fields[blocknr].type_in_memory)
1859  {
1860  case MEM_INT:
1861  hdf5_datatype = H5T_NATIVE_INT;
1862  break;
1863  case MEM_INT64:
1864  hdf5_datatype = H5T_NATIVE_INT64;
1865  break;
1866  case MEM_MY_ID_TYPE:
1867 #ifdef IDS_32BIT
1868  hdf5_datatype = H5T_NATIVE_UINT32;
1869 #else
1870  hdf5_datatype = H5T_NATIVE_UINT64;
1871 #endif
1872  break;
1873  case MEM_MY_INTPOS_TYPE:
1874 #ifdef POSITIONS_IN_32BIT
1875  hdf5_datatype = H5T_NATIVE_UINT32;
1876 #elif defined(POSITIONS_IN_64BIT)
1877  hdf5_datatype = H5T_NATIVE_UINT64;
1878 #else
1879  hdf5_datatype = Int128_memtype;
1880 #endif
1881  break;
1882  case MEM_FLOAT:
1883  hdf5_datatype = H5T_NATIVE_FLOAT;
1884  break;
1885  case MEM_DOUBLE:
1886  hdf5_datatype = H5T_NATIVE_DOUBLE;
1887  break;
1888  case MEM_MY_FLOAT:
1889  hdf5_datatype = H5T_NATIVE_MYFLOAT;
1890  break;
1891  case MEM_MY_DOUBLE:
1892  hdf5_datatype = H5T_NATIVE_MYDOUBLE;
1893  break;
1894  case MEM_MY_FILEOFFSET:
1895  hdf5_datatype = H5T_NATIVE_INT64;
1896  break;
1897  }
1898 
1899  return hdf5_datatype;
1900 }
1901 
1910 int IO_Def::get_values_per_blockelement(int blocknr)
1911 {
1912  if(blocknr < 0 || blocknr >= N_IO_Fields)
1913  Terminate("something is wrong here: blocknr=%d N_IO_Fields=%d", blocknr, N_IO_Fields);
1914 
1915  return IO_Fields[blocknr].values_per_block;
1916 }
1917 
1928 long long IO_Def::get_particles_in_block(int blocknr, long long *npart_file, int *typelist)
1929 {
1930  long long npart = 0;
1931 
1932  for(int i = 0; i < N_DataGroups; i++)
1933  typelist[i] = 0;
1934 
1935  switch(IO_Fields[blocknr].typelist)
1936  {
1937  case MASS_BLOCK:
1938  for(int i = 0; i < NTYPES; i++)
1939  {
1940  typelist[i] = (All.MassTable[i] == 0 && npart_file[i] > 0);
1941  npart += npart_file[i] * typelist[i];
1942  }
1943  return npart;
1944  break;
1945 
1946  case AGE_BLOCK:
1947  for(int i = 0; i < NTYPES; i++)
1948  {
1949  typelist[i] = (npart_file[i] > 0);
1950 
1951  if((file_format == FILEFORMAT_HDF5 && (i == 0 || i == 1 || i == 5)) || (file_format != FILEFORMAT_HDF5 && i != 4))
1952  typelist[i] = 0;
1953 
1954  npart += npart_file[i] * typelist[i];
1955  }
1956  return npart;
1957  break;
1958 
1959  case Z_BLOCK:
1960  for(int i = 0; i < NTYPES; i++)
1961  {
1962  typelist[i] = (npart_file[i] > 0);
1963  if((file_format == FILEFORMAT_HDF5 && (i == 1 || i == 5)) || (file_format != FILEFORMAT_HDF5 && (i != 0 && i != 4)))
1964  typelist[i] = 0;
1965 
1966  npart += npart_file[i] * typelist[i];
1967  }
1968  return npart;
1969  break;
1970 
1971  case GROUPS:
1972  npart = npart_file[0];
1973  typelist[0] = 1;
1974  return npart;
1975  break;
1976 
1977  case SUBGROUPS:
1978  npart = npart_file[1];
1979  typelist[1] = 1;
1980  return npart;
1981  break;
1982 
1983  case ID_BLOCK:
1984  npart = npart_file[2];
1985  typelist[2] = 1;
1986  return npart;
1987  break;
1988 
1989  case CURRSUBS:
1990  case PREVSUBS:
1991  case TREELINK:
1992  case TREELENGTH:
1993  case MASSMAPS:
1994  case GALSNAPS:
1995  npart = npart_file[0];
1996  typelist[0] = 1;
1997  return npart;
1998  break;
1999 
2000  case TREEHALOS:
2001  npart = npart_file[1];
2002  typelist[1] = 1;
2003  return npart;
2004  break;
2005 
2006  case TREETIMES:
2007  npart = npart_file[2];
2008  typelist[2] = 1;
2009  return npart;
2010  break;
2011 
2012  case TREETABLE:
2013  npart = npart_file[NTYPES];
2014  typelist[NTYPES] = 1;
2015  return npart;
2016  break;
2017 
2018  case HEALPIXTAB:
2019  npart = npart_file[NTYPES + 1];
2020  typelist[NTYPES + 1] = 1;
2021  return npart;
2022  break;
2023 
2024  default:
2025  for(int i = 0; i < N_DataGroups; i++)
2026  {
2027  if((IO_Fields[blocknr].typelist & (1 << i)) && npart_file[i] > 0)
2028  {
2029  typelist[i] = 1;
2030  npart += npart_file[i];
2031  }
2032  else
2033  typelist[i] = 0;
2034  }
2035  return npart;
2036  break;
2037  }
2038 
2039  return 0;
2040 }
2041 
2049 void IO_Def::get_Tab_IO_Label(int blocknr, char *label)
2050 {
2051  if(blocknr < 0 || blocknr >= N_IO_Fields)
2052  Terminate("something is wrong here: blocknr=%d N_IO_Fields=%d", blocknr, N_IO_Fields);
2053 
2054  strcpy(label, IO_Fields[blocknr].label);
2055 }
2056 
2064 void IO_Def::get_dataset_name(int blocknr, char *buf)
2065 {
2066  if(blocknr < 0 || blocknr >= N_IO_Fields)
2067  Terminate("something is wrong here: blocknr=%d N_IO_Fields=%d", blocknr, N_IO_Fields);
2068 
2069  strcpy(buf, IO_Fields[blocknr].datasetname);
2070 }
2071 
2072 void IO_Def::write_dataset_attributes(hid_t hdf5_dataset, int blocknr)
2073 {
2074  if(blocknr < 0 || blocknr >= N_IO_Fields)
2075  Terminate("something is wrong here: blocknr=%d N_IO_Fields=%d", blocknr, N_IO_Fields);
2076 
2077  if(IO_Fields[blocknr].hasunit == 0)
2078  return;
2079 
2081  {
2082  write_scalar_attribute(hdf5_dataset, "a_scaling", &IO_Fields[blocknr].a, H5T_NATIVE_DOUBLE);
2083  write_scalar_attribute(hdf5_dataset, "h_scaling", &IO_Fields[blocknr].h, H5T_NATIVE_DOUBLE);
2084  }
2085  else
2086  {
2087  double zero = 0;
2088  write_scalar_attribute(hdf5_dataset, "a_scaling", &zero, H5T_NATIVE_DOUBLE);
2089  write_scalar_attribute(hdf5_dataset, "h_scaling", &zero, H5T_NATIVE_DOUBLE);
2090  }
2091 
2092  write_scalar_attribute(hdf5_dataset, "length_scaling", &IO_Fields[blocknr].L, H5T_NATIVE_DOUBLE);
2093  write_scalar_attribute(hdf5_dataset, "mass_scaling", &IO_Fields[blocknr].M, H5T_NATIVE_DOUBLE);
2094  write_scalar_attribute(hdf5_dataset, "velocity_scaling", &IO_Fields[blocknr].V, H5T_NATIVE_DOUBLE);
2095 
2096  write_scalar_attribute(hdf5_dataset, "to_cgs", &IO_Fields[blocknr].c, H5T_NATIVE_DOUBLE);
2097 }
2098 
2106 void IO_Def::write_parameters_attributes_in_hdf5(hid_t handle)
2107 {
2108  for(int i = 0; i < All.NParameters; i++)
2109  {
2110  switch(All.ParametersType[i])
2111  {
2112  case PARAM_DOUBLE:
2113  write_scalar_attribute(handle, All.ParametersTag[i], All.ParametersValue[i], H5T_NATIVE_DOUBLE);
2114  break;
2115  case PARAM_STRING:
2116  write_string_attribute(handle, All.ParametersTag[i], (const char *)All.ParametersValue[i]);
2117  break;
2118  case PARAM_INT:
2119  write_scalar_attribute(handle, All.ParametersTag[i], All.ParametersValue[i], H5T_NATIVE_INT);
2120  break;
2121  }
2122  }
2123 }
2124 
2125 /*---------------------- Routine find a block in a snapfile -------------------*/
2126 void IO_Def::find_block(char *label, FILE *fd)
2127 {
2128  unsigned int blocksize = 0, blksize;
2129  char blocklabel[5] = {" "};
2130 
2131 #define FBSKIP \
2132  { \
2133  my_fread(&blksize, sizeof(int), 1, fd); \
2134  }
2135 
2136  rewind(fd);
2137 
2138  while(!feof(fd) && blocksize == 0)
2139  {
2140  FBSKIP;
2141 
2142  if(blksize != 8)
2143  {
2144  Terminate("Incorrect Format (blksize=%u)!\n", blksize);
2145  }
2146  else
2147  {
2148  my_fread(blocklabel, LABEL_LEN * sizeof(char), 1, fd);
2149  my_fread(&blocksize, sizeof(int), 1, fd);
2150 
2151  FBSKIP;
2152 
2153  if(strncmp(label, blocklabel, LABEL_LEN) != 0)
2154  {
2155  fseek(fd, blocksize, 1);
2156  blocksize = 0;
2157  }
2158  }
2159  }
2160  if(feof(fd))
2161  Terminate("Block '%c%c%c%c' not found !\n", label[0], label[1], label[2], label[3]);
2162 }
2163 
2164 void IO_Def::read_segment(const char *fname, int type, long long offset, long long count, int num_files)
2165 {
2166  long long nleft = count;
2167  long long offleft = offset;
2168  long long nskip = 0;
2169 
2170  for(int filenr = 0; filenr < num_files && nleft > 0; filenr++)
2171  {
2172  long long nloc = ntype_in_files[filenr * N_DataGroups + type];
2173 
2174  if(nloc > offleft) // we may have something in this file
2175  {
2176  long long nread;
2177 
2178  if(nloc - offleft > nleft) // there are more particles in the file then we need
2179  nread = nleft;
2180  else
2181  nread = nloc - offleft;
2182 
2183  /* now read partial list */
2184  read_single_file_segment(fname, filenr, type, offleft, nread, nskip, num_files);
2185 
2186  nleft -= nread;
2187  nskip += nread;
2188  offleft += nread;
2189  }
2190 
2191  offleft -= nloc;
2192  }
2193 
2194  if(nleft > 0)
2195  {
2196  for(int filenr = 0; filenr < num_files; filenr++)
2197  {
2198  long long nloc = ntype_in_files[filenr * N_DataGroups + type];
2199  printf("filenr=%d: nloc=%lld\n", filenr, nloc);
2200  }
2201  Terminate("Not all desired entries read: nleft=%lld type=%d\n", nleft, type);
2202  }
2203 }
2204 
2205 void IO_Def::read_single_file_segment(const char *basename, int filenr, int type, long long offset, unsigned long long count,
2206  long long storage_offset, int num_files)
2207 {
2208  int bytes_per_blockelement_in_file = 0;
2209  hid_t hdf5_file = 0, hdf5_grp = 0, hdf5_dataspace_in_file;
2210  hid_t hdf5_dataspace_in_memory, hdf5_dataset;
2211  FILE *fd = 0;
2212  char fname[MAXLEN_PATH];
2213 
2214  if(num_files > 1)
2215  {
2216  if(file_format == FILEFORMAT_HDF5)
2217  sprintf(fname, "%s.%d.hdf5", basename, filenr);
2218  else
2219  sprintf(fname, "%s.%d", basename, filenr);
2220  }
2221  else
2222  {
2223  if(file_format == FILEFORMAT_HDF5)
2224  sprintf(fname, "%s.hdf5", basename);
2225  else
2226  sprintf(fname, "%s", basename);
2227  }
2228 
2229  /* open file */
2230  if(file_format == FILEFORMAT_HDF5)
2231  {
2232  hdf5_file = my_H5Fopen(fname, H5F_ACC_RDONLY, H5P_DEFAULT);
2233  char buf[MAXLEN_PATH];
2234  get_datagroup_name(type, buf);
2235  hdf5_grp = my_H5Gopen(hdf5_file, buf);
2236  }
2237  else if(file_format == FILEFORMAT_LEGACY1 || file_format == FILEFORMAT_LEGACY2)
2238  {
2239  if(!(fd = fopen(fname, "r")))
2240  Terminate("can't open file `%s' for reading initial conditions.\n", fname);
2241 
2242  unsigned int blksize1, blksize2;
2243 
2244  if(file_format == FILEFORMAT_LEGACY2)
2245  {
2246  int nextblock;
2247  char label[LABEL_LEN];
2248  my_fread(&blksize1, sizeof(int), 1, fd);
2249  my_fread(&label, sizeof(char), LABEL_LEN, fd);
2250  my_fread(&nextblock, sizeof(int), 1, fd);
2251  my_fread(&blksize2, sizeof(int), 1, fd);
2252  }
2253 
2254  my_fread(&blksize1, sizeof(int), 1, fd);
2255  my_fread(header_buf, header_size, 1, fd);
2256  my_fread(&blksize2, sizeof(int), 1, fd);
2257  }
2258  else
2259  Terminate("unknown ICFormat");
2260 
2261  long long npart[N_DataGroups];
2262  for(int i = 0; i < N_DataGroups; i++)
2263  npart[i] = ntype_in_files[filenr * N_DataGroups + i];
2264 
2265  for(int blocknr = 0; blocknr < N_IO_Fields; blocknr++)
2266  {
2267  if(IO_Fields[blocknr].type_in_memory != MEM_MY_FILEOFFSET)
2268  {
2269  unsigned int blksize1, blksize2;
2270  int typelist[N_DataGroups];
2271  int bytes_per_blockelement = get_bytes_per_memory_blockelement(blocknr, 1);
2272  long long npart_in_block = get_particles_in_block(blocknr, npart, &typelist[0]);
2273  hid_t hdf5_memory_datatype = get_hdf5_memorytype_of_block(blocknr);
2274  char dname[MAXLEN_PATH];
2275  get_dataset_name(blocknr, dname);
2276 
2277  if(npart_in_block > 0 && typelist[type] > 0 && IO_Fields[blocknr].read_flag != SKIP_ON_READ)
2278  {
2279  if(file_format == FILEFORMAT_LEGACY2)
2280  {
2281  char expected_label[LABEL_LEN + 1];
2282  get_Tab_IO_Label(blocknr, expected_label);
2283  find_block(expected_label, fd);
2284  }
2285 
2286  if(file_format == FILEFORMAT_LEGACY1 || file_format == FILEFORMAT_LEGACY2)
2287  {
2288  my_fread(&blksize1, sizeof(int), 1, fd);
2289  bytes_per_blockelement_in_file = blksize1 / npart_in_block;
2290  }
2291 
2292  void *CommBuffer = (char *)Mem.mymalloc("CommBuffer", bytes_per_blockelement * count);
2293 
2294  if(file_format == FILEFORMAT_LEGACY1 || file_format == FILEFORMAT_LEGACY2)
2295  {
2296  fseek(fd, bytes_per_blockelement_in_file * offset, SEEK_CUR);
2297 
2298  if(bytes_per_blockelement_in_file != bytes_per_blockelement)
2299  {
2300  char *CommAuxBuffer = (char *)Mem.mymalloc("CommAuxBuffer", bytes_per_blockelement_in_file * count);
2301  my_fread(CommAuxBuffer, bytes_per_blockelement_in_file, count, fd);
2302  type_cast_data((char *)CommAuxBuffer, bytes_per_blockelement_in_file, (char *)CommBuffer, bytes_per_blockelement,
2303  count, blocknr);
2304  Mem.myfree(CommAuxBuffer);
2305  }
2306  else
2307  my_fread(CommBuffer, bytes_per_blockelement, count, fd);
2308 
2309  fseek(fd, bytes_per_blockelement_in_file * (npart_in_block - offset - count), SEEK_CUR);
2310 
2311  my_fread(&blksize2, sizeof(int), 1, fd);
2312  if(blksize1 != blksize2)
2313  {
2314  char buf[MAXLEN_PATH];
2315  sprintf(buf, "incorrect block-sizes detected!\n Task=%d blocknr=%d blksize1=%d blksize2=%d\n", ThisTask,
2316  blocknr, blksize1, blksize2);
2317  if(blocknr == 2) /* block number 2 is always IDs */
2318  strcat(buf, "Possible mismatch of 32bit and 64bit ID's in IC file and GADGET compilation !\n");
2319  Terminate(buf);
2320  }
2321  }
2322 
2323  if(file_format == FILEFORMAT_HDF5)
2324  {
2325  hdf5_dataset = my_H5Dopen_if_existing(hdf5_grp, dname);
2326  int rank;
2327  hsize_t dims[2], nelem[2], start[2];
2328 
2329  dims[0] = npart[type];
2330  dims[1] = get_values_per_blockelement(blocknr);
2331  if(dims[1] == 1)
2332  rank = 1;
2333  else
2334  rank = 2;
2335 
2336  hdf5_dataspace_in_file = my_H5Screate_simple(rank, dims, NULL);
2337 
2338  dims[0] = count;
2339  hdf5_dataspace_in_memory = my_H5Screate_simple(rank, dims, NULL);
2340 
2341  start[0] = offset;
2342  start[1] = 0;
2343 
2344  nelem[0] = count;
2345  nelem[1] = get_values_per_blockelement(blocknr);
2346 
2347  my_H5Sselect_hyperslab(hdf5_dataspace_in_file, H5S_SELECT_SET, start, NULL, nelem, NULL);
2348 
2349  // Test if dataset was present
2350  if(hdf5_dataset < 0)
2351  {
2352  // no, pad with zeros
2353  printf("\tDataset %s not present for particle type %d, using zero.\n", dname, type);
2354  memset(CommBuffer, 0, dims[0] * dims[1] * my_H5Tget_size(hdf5_memory_datatype));
2355  }
2356  else
2357  {
2358  hid_t hdf5_file_datatype = H5Dget_type(hdf5_dataset);
2359  byte_count += dims[0] * dims[1] * my_H5Tget_size(hdf5_file_datatype); /* for I/O performance measurement */
2360  H5Tclose(hdf5_file_datatype);
2361 
2362  my_H5Dread(hdf5_dataset, hdf5_memory_datatype, hdf5_dataspace_in_memory, hdf5_dataspace_in_file, H5P_DEFAULT,
2363  CommBuffer, dname);
2364  my_H5Dclose(hdf5_dataset, dname);
2365  }
2366  my_H5Sclose(hdf5_dataspace_in_memory, H5S_SIMPLE);
2367  my_H5Sclose(hdf5_dataspace_in_file, H5S_SIMPLE);
2368  }
2369 
2370  empty_read_buffer(blocknr, storage_offset, count, type, 0, CommBuffer);
2371 
2372  Mem.myfree(CommBuffer);
2373  }
2374  else
2375  {
2376  if(file_format == FILEFORMAT_LEGACY1 && npart_in_block > 0)
2377  {
2378  my_fread(&blksize1, sizeof(int), 1, fd);
2379  bytes_per_blockelement_in_file = blksize1 / npart_in_block;
2380  fseek(fd, bytes_per_blockelement_in_file * npart_in_block, SEEK_CUR);
2381  my_fread(&blksize2, sizeof(int), 1, fd);
2382  if(blksize1 != blksize2)
2383  {
2384  char buf[MAXLEN_PATH];
2385  sprintf(buf, "incorrect block-sizes detected!\n Task=%d blocknr=%d blksize1=%d blksize2=%d\n", ThisTask,
2386  blocknr, blksize1, blksize2);
2387  if(blocknr == 2) /* block number 2 is always IDs */
2388  strcat(buf, "Possible mismatch of 32bit and 64bit ID's in IC file and GADGET compilation !\n");
2389  Terminate(buf);
2390  }
2391  }
2392  }
2393  }
2394  }
2395 
2396  if(file_format == FILEFORMAT_LEGACY1 || file_format == FILEFORMAT_LEGACY2)
2397  fclose(fd);
2398 
2399  if(file_format == FILEFORMAT_HDF5)
2400  {
2401  char buf[MAXLEN_PATH];
2402  get_datagroup_name(type, buf);
2403  my_H5Gclose(hdf5_grp, buf);
2404  my_H5Fclose(hdf5_file, fname);
2405  }
2406 
2407  read_increase_numbers(type, count);
2408 }
2409 
2410 void IO_Def::rename_file_to_bak_if_it_exists(char *fname)
2411 {
2412  char fin[MAXLEN_PATH], buf[2 * MAXLEN_PATH];
2413 
2414  strcpy(fin, fname);
2415 
2416  char *p = strrchr(fin, '/');
2417  if(p)
2418  {
2419  *p = 0;
2420  sprintf(buf, "%s/bak-%s", fin, p + 1);
2421  }
2422  else
2423  sprintf(buf, "bak-%s", fname);
2424 
2425  if(FILE *fcheck = fopen(fname, "r")) // check if file already exists, if yes, try to rename the existing file
2426  {
2427  fclose(fcheck);
2428 
2429  if(!fopen(buf, "r")) // only do the rename of the old file if the back-up file doesn't exist yet
2430  {
2431  mpi_printf("%s rename '%s' to '%s'\n", info, fname, buf);
2432  rename(fname, buf);
2433  }
2434  }
2435 }
2436 
2437 void IO_Def::alloc_and_read_ntype_in_files(const char *fname, int num_files)
2438 {
2439  ntype_in_files = (long long *)Mem.mymalloc_movable(&ntype_in_files, "ntype_in_files", num_files * N_DataGroups * sizeof(long long));
2440 
2441  for(int filenr = 0; filenr < num_files; filenr++)
2442  {
2443  char buf[3 * MAXLEN_PATH];
2444 
2445  if(num_files > 1)
2446  {
2447  sprintf(buf, "%s.%d", fname, filenr);
2448  if(file_format == 3)
2449  sprintf(buf, "%s.%d.hdf5", fname, filenr);
2450  }
2451  else
2452  {
2453  sprintf(buf, "%s", fname);
2454  if(file_format == 3)
2455  sprintf(buf, "%s.hdf5", fname);
2456  }
2457 
2458  if(file_format == 3)
2459  {
2460  read_header_fields(buf);
2461  }
2462  else if(file_format == 1 || file_format == 2)
2463  {
2464  FILE *fd = 0;
2465 
2466  if(!(fd = fopen(buf, "r")))
2467  Terminate("can't open file `%s' for reading initial conditions.\n", fname);
2468 
2469  int blksize1, blksize2;
2470 
2471  if(file_format == 2)
2472  {
2473  char label[4];
2474  int nextblock;
2475  my_fread(&blksize1, sizeof(int), 1, fd);
2476  my_fread(&label, sizeof(char), 4, fd);
2477  my_fread(&nextblock, sizeof(int), 1, fd);
2478  printf("READIC: Reading header => '%c%c%c%c' (%d byte)\n", label[0], label[1], label[2], label[3], nextblock);
2479  my_fread(&blksize2, sizeof(int), 1, fd);
2480  }
2481 
2482  my_fread(&blksize1, sizeof(int), 1, fd);
2483  my_fread(header_buf, header_size, 1, fd);
2484  my_fread(&blksize2, sizeof(int), 1, fd);
2485 
2486  if(blksize1 != blksize2)
2487  Terminate("incorrect header format, blksize1=%d blksize2=%d header_size=%d\n", blksize1, blksize2, (int)header_size);
2488 
2489  fclose(fd);
2490  }
2491 
2492  long long n_type[N_DataGroups], npart[N_DataGroups];
2493 
2494  read_file_header(fname, filenr, 0, 0, n_type, npart, NULL);
2495 
2496  for(int type = 0; type < N_DataGroups; type++)
2497  ntype_in_files[filenr * N_DataGroups + type] = npart[type];
2498  }
2499 }
global_data_all_processes All
Definition: main.cc:40
Definition: io.h:129
virtual int get_filenr_from_header(void)=0
int N_DataGroups
Definition: io.h:139
void alloc_and_read_ntype_in_files(const char *fname, int num_files)
Definition: io.cc:2437
virtual void * get_base_address_of_structure(enum arrays array, int index)=0
void write_multiple_files(char *fname, int numfilesperdump, int append_flag=0, int chunk_size=0)
Definition: io.cc:700
long long * ntype_in_files
Definition: io.h:176
int Max_IO_Fields
Definition: io.h:141
virtual ~IO_Def()
Definition: io.cc:47
void * header_buf
Definition: io.h:174
void read_segment(const char *fname, int type, long long offset, long long count, int numfiles)
Definition: io.cc:2164
virtual void set_type_of_element(int index, int type)=0
virtual void read_increase_numbers(int type, int n_for_this_task)=0
virtual void fill_file_header(int writeTask, int lastTask, long long *nloc_part, long long *npart)=0
size_t header_size
Definition: io.h:173
void init_field(const char *label, const char *datasetname, enum types_in_memory type_in_memory, enum types_in_file type_in_file_output, enum read_flags read_flag, int values_per_block, enum arrays array, void *pointer_to_field, void(*io_func)(IO_Def *, int, int, void *, int), int typelist_bitmask, int hasunits, double a, double h, double L, double M, double V, double c, bool compression_on=false)
Definition: io.cc:66
int N_IO_Fields
Definition: io.h:140
virtual void read_header_fields(const char *fname)=0
virtual void read_file_header(const char *fname, int filenr, int readTask, int lastTask, long long *nloc_part, long long *npart, int *nstart)=0
char info[100]
Definition: io.h:177
virtual void set_filenr_in_header(int)=0
enum file_contents type_of_file
Definition: io.h:187
virtual void write_header_fields(hid_t)=0
void read_files_driver(const char *fname, int rep, int numfiles)
Definition: io.cc:227
virtual void get_datagroup_name(int grnr, char *gname)=0
void read_single_file_segment(const char *fname, int filenr, int type, long long offset, unsigned long long count, long long storage_offset, int numfiles)
Definition: io.cc:2205
virtual int get_type_of_element(int index)=0
int find_files(const char *fname, const char *fname_multiple)
This function determines on how many files a given snapshot or group/desc catalogue is distributed.
Definition: io.cc:132
void write_compile_time_options_in_hdf5(hid_t handle)
long long byte_count
size_t my_fread(void *ptr, size_t size, size_t nmemb, FILE *stream)
A wrapper for the fread() function.
size_t my_fwrite(const void *ptr, size_t size, size_t nmemb, FILE *stream)
A wrapper for the fwrite() function.
char ParametersType[MAX_PARAMETERS]
Definition: parameters.h:46
char ParametersTag[MAX_PARAMETERS][MAXLEN_PARAM_TAG]
Definition: parameters.h:44
int NParameters
Definition: parameters.h:42
void * ParametersValue[MAX_PARAMETERS]
Definition: parameters.h:45
int ThisNode
Definition: setcomm.h:36
void mpi_printf(const char *fmt,...)
Definition: setcomm.h:55
int ThisTask
Definition: setcomm.h:33
int NTask
Definition: setcomm.h:32
int RankInThisNode
Definition: setcomm.h:39
MPI_Comm Communicator
Definition: setcomm.h:31
#define FILEFORMAT_LEGACY2
Definition: constants.h:18
#define FILEFORMAT_LEGACY1
Definition: constants.h:17
#define COMMBUFFERSIZE
Definition: constants.h:47
#define NTYPES
Definition: constants.h:308
#define FILEFORMAT_HDF5
Definition: constants.h:19
#define MAXLEN_PATH
Definition: constants.h:300
#define H5T_NATIVE_MYDOUBLE
Definition: dtypes.h:92
float MyDouble
Definition: dtypes.h:87
float MyFloat
Definition: dtypes.h:86
uint32_t MyIntPosType
Definition: dtypes.h:35
#define H5T_NATIVE_MYFLOAT
Definition: dtypes.h:91
@ RST_BEGIN
Definition: dtypes.h:313
unsigned int MyIDType
Definition: dtypes.h:68
hid_t my_H5Gopen(hid_t loc_id, const char *groupname)
Definition: hdf5_util.cc:311
herr_t my_H5Sselect_hyperslab(hid_t space_id, H5S_seloper_t op, const hsize_t *start, const hsize_t *stride, const hsize_t *count, const hsize_t *block)
Definition: hdf5_util.cc:549
herr_t my_H5Dread(hid_t dataset_id, hid_t mem_type_id, hid_t mem_space_id, hid_t file_space_id, hid_t xfer_plist_id, void *buf, const char *datasetname)
Definition: hdf5_util.cc:385
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
hid_t Halfprec_memtype
Definition: hdf5_util.cc:29
size_t my_H5Tget_size(hid_t datatype_id)
Definition: hdf5_util.cc:564
herr_t my_H5Dset_extent(hid_t dset_id, const hsize_t size[])
Definition: hdf5_util.cc:336
hid_t Int48_memtype
Definition: hdf5_util.cc:30
hid_t my_H5Gcreate(hid_t loc_id, const char *groupname, size_t size_hint)
Definition: hdf5_util.cc:174
herr_t my_H5Pclose(hid_t plist)
Definition: hdf5_util.cc:468
hid_t my_H5Screate_simple(int rank, const hsize_t *current_dims, const hsize_t *maximum_dims)
Definition: hdf5_util.cc:283
herr_t my_H5Fclose(hid_t file_id, const char *fname)
Definition: hdf5_util.cc:482
herr_t my_H5Dclose(hid_t dataset_id, const char *datasetname)
Definition: hdf5_util.cc:443
hid_t my_H5Fcreate(const char *fname, unsigned int flags, hid_t fcpl_id, hid_t fapl_id)
Definition: hdf5_util.cc:158
hid_t my_H5Dopen_if_existing(hid_t file_id, const char *datasetname)
Definition: hdf5_util.cc:352
hid_t my_H5Dcreate(hid_t loc_id, const char *datasetname, hid_t type_id, hid_t space_id, hid_t dcpl_id)
Definition: hdf5_util.cc:190
herr_t my_H5Sclose(hid_t dataspace_id, H5S_class_t type)
Definition: hdf5_util.cc:496
hid_t Int128_memtype
Definition: hdf5_util.cc:31
void write_scalar_attribute(hid_t handle, const char *attr_name, const void *buf, hid_t mem_type_id)
Definition: hdf5_util.cc:621
herr_t my_H5Dwrite(hid_t dataset_id, hid_t mem_type_id, hid_t mem_space_id, hid_t file_space_id, hid_t xfer_plist_id, const void *buf, const char *datasetname)
Definition: hdf5_util.cc:206
void write_string_attribute(hid_t handle, const char *attr_name, const char *buf)
Definition: hdf5_util.cc:654
#define COMPRESSION_CHUNKSIZE
Definition: hdf5_util.h:17
#define SKIP
#define FBSKIP
read_flags
Definition: io.h:123
@ SKIP_ON_READ
Definition: io.h:125
@ GALSNAPS
Definition: io.h:118
@ MASSMAPS
Definition: io.h:115
@ TREELENGTH
Definition: io.h:112
@ PREVSUBS
Definition: io.h:110
@ AGE_BLOCK
Definition: io.h:105
@ TREETIMES
Definition: io.h:117
@ TREELINK
Definition: io.h:114
@ Z_BLOCK
Definition: io.h:106
@ GROUPS
Definition: io.h:107
@ MASS_BLOCK
Definition: io.h:104
@ CURRSUBS
Definition: io.h:111
@ SUBGROUPS
Definition: io.h:108
@ ID_BLOCK
Definition: io.h:109
@ TREEHALOS
Definition: io.h:113
@ TREETABLE
Definition: io.h:116
@ HEALPIXTAB
Definition: io.h:119
arrays
Definition: io.h:30
@ A_MTRP
Definition: io.h:40
@ A_NONE
Definition: io.h:31
@ A_P
Definition: io.h:33
types_in_file
Definition: io.h:65
@ FILE_FLOAT
Definition: io.h:74
@ FILE_MY_INTPOS_TYPE
Definition: io.h:72
@ FILE_DOUBLE
Definition: io.h:73
@ FILE_HALF
Definition: io.h:70
@ 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
#define LABEL_LEN
Definition: io.h:26
#define DATASETNAME_LEN
Definition: io.h:27
types_in_memory
Definition: io.h:78
@ MEM_INT
Definition: io.h:79
@ MEM_DOUBLE
Definition: io.h:84
@ MEM_MY_FLOAT
Definition: io.h:85
@ MEM_MY_FILEOFFSET
Definition: io.h:87
@ MEM_MY_ID_TYPE
Definition: io.h:81
@ MEM_MY_INTPOS_TYPE
Definition: io.h:82
@ MEM_FLOAT
Definition: io.h:83
@ MEM_MY_DOUBLE
Definition: io.h:86
@ MEM_INT64
Definition: io.h:80
@ FILE_IS_SNAPSHOT
Definition: io.h:53
@ FILE_IS_LIGHTCONE
Definition: io.h:59
#define Terminate(...)
Definition: macros.h:19
#define TAG_HEADER
Definition: mpi_utils.h:26
#define TAG_NFORTHISTASK
Definition: mpi_utils.h:39
#define TAG_PDATA
Definition: mpi_utils.h:27
#define TAG_N
Definition: mpi_utils.h:25
#define TAG_KEY
Definition: mpi_utils.h:29
memory Mem
Definition: main.cc:44
#define PARAM_INT
Definition: parameters.h:20
#define PARAM_STRING
Definition: parameters.h:19
#define PARAM_DOUBLE
Definition: parameters.h:18
enum restart_options RestartFlag
Definition: allvars.h:68
double MassTable[NTYPES]
Definition: allvars.h:269
void myflush(FILE *fstream)
Definition: system.cc:125