GADGET-4
snap_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 
23 #include "../cooling_sfr/cooling.h"
24 #include "../data/allvars.h"
25 #include "../data/dtypes.h"
26 #include "../data/intposconvert.h"
27 #include "../data/mymalloc.h"
28 #include "../fof/fof.h"
29 #include "../gitversion/version.h"
30 #include "../gravtree/gravtree.h"
31 #include "../io/hdf5_util.h"
32 #include "../io/io.h"
33 #include "../io/snap_io.h"
34 #include "../lightcone/lightcone.h"
35 #include "../logs/logs.h"
36 #include "../logs/timer.h"
37 #include "../main/main.h"
38 #include "../main/simulation.h"
39 #include "../mpi_utils/mpi_utils.h"
40 #include "../sort/peano.h"
41 #include "../src/pm/pm.h"
42 #include "../system/system.h"
43 
51 {
52  Sp = Sp_ptr;
53 
54  this->N_IO_Fields = 0;
55  this->N_DataGroups = NTYPES + 1; // the last data group is a tree table used only for storing/reading tree-reordered particle data
56  this->header_size = sizeof(header);
57  this->header_buf = &header;
59  sprintf(this->info, "SNAPSHOT: writing snapshot");
60 
61 #ifdef OUTPUT_COORDINATES_AS_INTEGERS
62  init_field("IPOS", "IntCoordinates", MEM_MY_INTPOS_TYPE, FILE_MY_INTPOS_TYPE, READ_IF_PRESENT, 3, A_P, NULL, io_func_intpos,
63  ALL_TYPES, 1, 1., -1., 1., 0., 0., All.UnitLength_in_cm * Sp->FacIntToCoord, true);
64 #else
65  init_field("POS ", "Coordinates", MEM_MY_DOUBLE, FILE_MY_IO_FLOAT, READ_IF_PRESENT, 3, A_P, NULL, io_func_pos, ALL_TYPES, 1, 1., -1.,
66  1., 0., 0., All.UnitLength_in_cm, true);
67 #endif
68 
69 #ifdef OUTPUT_VELOCITIES_IN_HALF_PRECISION
70  init_field("VEL ", "Velocities", MEM_MY_FLOAT, FILE_HALF, READ_IF_PRESENT, 3, A_NONE, NULL, io_func_vel, ALL_TYPES, 1, 0.5, 0., 0.,
72 #else
73  init_field("VEL ", "Velocities", MEM_MY_FLOAT, FILE_MY_IO_FLOAT, READ_IF_PRESENT, 3, A_NONE, NULL, io_func_vel, ALL_TYPES, 1, 0.5,
74  0., 0., 0., 1., All.UnitVelocity_in_cm_per_s);
75 #endif
76 
77 #ifdef OUTPUT_ACCELERATION
78 #ifdef OUTPUT_ACCELERATIONS_IN_HALF_PRECISION
79  All.accel_normalize_fac = 10.0 * All.Hubble * (100.0 * 1.0e5 / All.UnitVelocity_in_cm_per_s);
80 
81  init_field("ACCE", "Acceleration", MEM_MY_FLOAT, FILE_HALF, SKIP_ON_READ, 3, A_NONE, 0, io_func_accel, ALL_TYPES, 1, -2.0, 1, -1, 0,
83 #else
85 
86  init_field("ACCE", "Acceleration", MEM_MY_FLOAT, FILE_MY_IO_FLOAT, SKIP_ON_READ, 3, A_NONE, 0, io_func_accel, ALL_TYPES, 1, -2.0, 1,
88 #endif
89 
90  /* hydro acceleration */
91  init_field("HACC", "HydroAcceleration", MEM_MY_FLOAT, FILE_MY_IO_FLOAT, READ_IF_PRESENT, 3, A_SPHP, &Sp->SphP[0].HydroAccel, 0,
92  GAS_ONLY, 0, 0, 0, 0, 0, 0, 0);
93 #endif
94 
95  init_field("ID ", "ParticleIDs", MEM_MY_ID_TYPE, FILE_MY_ID_TYPE, READ_IF_PRESENT, 1, A_P, NULL, io_func_id, ALL_TYPES, 0, 0, 0, 0,
96  0, 0, 0, true);
97 
98 #ifndef LEAN
99  init_field("MASS", "Masses", MEM_MY_DOUBLE, FILE_MY_IO_FLOAT, READ_IF_PRESENT, 1, A_P, NULL, io_func_mass, MASS_BLOCK, 1, 0., -1.,
100  0., 1., 0., All.UnitMass_in_g);
101 #endif
102 
103 #ifdef SECOND_ORDER_LPT_ICS
104  init_field("LPTM", "SecondOrderICMasses", MEM_FLOAT, FILE_NONE, READ_IF_PRESENT, 1, A_P, &Sp->P[0].OldAcc, NULL, ALL_TYPES, 0, 0, 0,
105  0, 0, 0, 0);
106 #endif
107 
108  init_field("U ", "InternalEnergy", MEM_MY_FLOAT, FILE_MY_IO_FLOAT, READ_IF_PRESENT, 1, A_NONE, NULL, io_func_u, GAS_ONLY, 1, 0.,
110 
111  init_field("RHO ", "Density", MEM_MY_FLOAT, FILE_MY_IO_FLOAT, READ_IF_PRESENT, 1, A_SPHP, &Sp->SphP[0].Density, NULL, GAS_ONLY, 1,
112  -3., 2., -3., 1., 0., All.UnitDensity_in_cgs);
113 
114  init_field("HSML", "SmoothingLength", MEM_MY_FLOAT, FILE_MY_IO_FLOAT, READ_IF_PRESENT, 1, A_SPHP, &Sp->SphP[0].Hsml, NULL, GAS_ONLY,
115  1, 1., -1., 1., 0., 0., All.UnitLength_in_cm);
116 
117 #ifdef STARFORMATION
118 
120  1, A_NONE, 0, io_func_sfr, GAS_ONLY, 1, 0., 0., -1., 1., 1., SOLAR_MASS / SEC_PER_YEAR);
121 
122  init_field("AGE ", "StellarFormationTime", MEM_MY_FLOAT, FILE_MY_IO_FLOAT, READ_IF_PRESENT, 1, A_P, &Sp->P[0].StellarAge, NULL,
123  AGE_BLOCK, /* stellar formation time */
124  0, 0, 0, 0, 0, 0, 0);
125 
126  init_field("Z ", "Metallicity", MEM_MY_FLOAT, FILE_MY_IO_FLOAT, READ_IF_PRESENT, 1, A_NONE, 0, io_func_metallicity,
127  Z_BLOCK, /* gas and star metallicity */
128  0, 0, 0, 0, 0, 0, 0);
129 #endif
130 
131 #if defined(PRESSURE_ENTROPY_SPH) && defined(OUTPUT_PRESSURE_SPH_DENSITY)
132  init_field("PRHO", "PressureSphDensity", MEM_MY_FLOAT, FILE_MY_IO_FLOAT, READ_IF_PRESENT, 1, A_SPHP, &Sp->SphP[0].PressureSphDensity,
133  NULL, GAS_ONLY, /* Pressure density */
134  1, -3., 2., -3., 1., 0., All.UnitDensity_in_cgs);
135 #endif
136 
137 #ifdef OUTPUT_PRESSURE
138  init_field("PRES", "Pressure", MEM_MY_FLOAT, FILE_MY_IO_FLOAT, SKIP_ON_READ, 1, A_NONE, 0, io_func_pressure,
139  GAS_ONLY, /* particle pressure */
141 #endif
142 
143 #if defined(TIMEDEP_ART_VISC) && defined(OUTPUT_VISCOSITY_PARAMETER)
144  init_field("ALP ", "ArtificialViscosityParameter", MEM_MY_FLOAT, FILE_MY_IO_FLOAT, READ_IF_PRESENT, 1, A_SPHP, &Sp->SphP[0].Alpha,
145  NULL, GAS_ONLY, 1, -3., 2., -3., 1., 0., 1);
146 #endif
147 
148 #ifdef OUTPUT_ENTROPY
149  init_field("ENTR", "Entropy", MEM_MY_FLOAT, FILE_MY_IO_FLOAT, READ_IF_PRESENT, 1, A_SPHP, &Sp->SphP[0].Entropy, 0,
150  GAS_ONLY, /* particle entropy */
151  0, 0, 0, 0, 0, 0, 0);
152 #endif
153 
154 #ifdef COOLING
155  init_field("NE ", "ElectronAbundance", MEM_MY_FLOAT, FILE_MY_IO_FLOAT, READ_IF_PRESENT, 1, A_SPHP, &Sp->SphP[0].Ne, 0,
156  GAS_ONLY, /* electron abundance */
157  0, 0, 0, 0, 0, 0, 0);
158 #endif
159 
160 #ifdef OUTPUT_POTENTIAL
161  init_field("POT ", "Potential", MEM_MY_FLOAT, FILE_MY_IO_FLOAT, READ_IF_PRESENT, 1, A_P, &Sp->P[0].Potential, 0,
162  ALL_TYPES, /* potential */
163  1, -1., 0., 0., 0., 2., All.UnitVelocity_in_cm_per_s * All.UnitVelocity_in_cm_per_s);
164 #endif
165 
166 #ifdef OUTPUT_CHANGEOFENTROPY
167  init_field("ENDT", "RateOfChangeOfEntropy", MEM_MY_FLOAT, FILE_MY_IO_FLOAT, READ_IF_PRESENT, 1, A_SPHP, &Sp->SphP[0].DtEntropy, 0,
168  GAS_ONLY, /* particle entropy change */
169  0, 0, 0, 0, 0, 0, 0);
170 #endif
171 
172 #ifdef OUTPUT_TIMESTEP
173  init_field("TSTP", "TimeStep", MEM_MY_FLOAT, FILE_MY_IO_FLOAT, SKIP_ON_READ, 1, A_NONE, 0, io_func_timestep,
174  ALL_TYPES, /* time step */
175  0, 0, 0, 0, 0, 0, 0);
176 
177  init_field("TSTH", "TimeStepHydro", MEM_MY_FLOAT, FILE_MY_IO_FLOAT, SKIP_ON_READ, 1, A_NONE, 0, io_func_timestephydro,
178  GAS_ONLY, /* hydro time step */
179  0, 0, 0, 0, 0, 0, 0);
180 #endif
181 
182 #ifdef OUTPUT_DIVVEL
183  init_field("DIVV", "VelocityDivergence", MEM_MY_FLOAT, FILE_MY_IO_FLOAT, READ_IF_PRESENT, 1, A_SPHP, &Sp->SphP[0].DivVel, 0,
184  GAS_ONLY, /* hydro velocity divergence */
185  0, 0, 0, 0, 0, 0, 0);
186 #endif
187 
188 #ifdef OUTPUT_CURLVEL
189  init_field("ROTV", "VelocityCurl", MEM_MY_FLOAT, FILE_MY_IO_FLOAT, READ_IF_PRESENT, 1, A_SPHP, &Sp->SphP[0].CurlVel, 0,
190  GAS_ONLY, /* absolute value of rot v */
191  0, 0, 0, 0, 0, 0, 0);
192 #endif
193 
194 #ifdef OUTPUT_VELOCITY_GRADIENT
195 #ifndef IMPROVED_VELOCITY_GRADIENTS
196 #error "The option OUTPUT_VELOCITY_GRADIENT requires IMPROVED_VELOCITY_GRADIENTS"
197 #endif
198  init_field("GRAV", "VelocityGradient", MEM_MY_FLOAT, FILE_MY_IO_FLOAT, READ_IF_PRESENT, 9, A_SPHP, &Sp->SphP[0].dvel[0][0], 0,
199  GAS_ONLY, 0, 0, 0, 0, 0, 0, 0);
200 #endif
201 
202 #ifdef OUTPUT_COOLHEAT
203  init_field("COHE", "CoolingHeatingEnergy", MEM_MY_FLOAT, FILE_MY_IO_FLOAT, READ_IF_PRESENT, 1, A_SPHP, &Sp->SphP[0].CoolHeat, 0,
204  GAS_ONLY, 0, 0, 0, 0, 0, 0, 0);
205 #endif
206 
207 #if defined(SUBFIND) && defined(SUBFIND_STORE_LOCAL_DENSITY)
208 
210  A_PS, &Sp->PS[0].SubfindDensity, 0, ALL_TYPES, /* subfind density */
211  1, -3., 2., -3., 1., 0., All.UnitDensity_in_cgs);
212 
214  A_PS, &Sp->PS[0].SubfindHsml, 0, ALL_TYPES, /* subfind hsml */
215  1, 1., -1., 1., 0., 0., All.UnitLength_in_cm);
216 
218  A_PS, &Sp->PS[0].SubfindVelDisp, 0, ALL_TYPES, /* subfind velocity dispersion */
219  1, 0., 0., 0., 0., 1., All.UnitVelocity_in_cm_per_s);
220 #endif
221 
222 #if defined(GADGET2_HEADER) && defined(REARRANGE_OPTION) && defined(MERGERTREE)
223  Terminate("GADGET2_HEADER does not work together with REARRANGE_OPTION\n");
224 #endif
225 }
226 
227 #if defined(REARRANGE_OPTION) && defined(MERGERTREE)
228 void snap_io::init_extra(simparticles *Sp_ptr, mergertree *MergerTree_ptr)
229 {
230  Sp = Sp_ptr;
231  MergerTree = MergerTree_ptr;
232 
233  init_field("MTRI", "TreeID", MEM_INT64, FILE_INT64, READ_IF_PRESENT, 1, A_P, &Sp->P[0].TreeID, NULL, ALL_TYPES, 0, 0, 0, 0, 0, 0, 0);
234 
235  init_field("MTRL", "ParticleCount", MEM_INT, FILE_INT, READ_IF_PRESENT, 1, A_TT, &MergerTree->TreeTable[0].HaloCount, NULL,
236  TREETABLE, 0, 0, 0, 0, 0, 0, 0);
237  init_field("MTRS", "ParticleFirst", MEM_INT64, FILE_INT64, READ_IF_PRESENT, 1, A_TT, &MergerTree->TreeTable[0].FirstHalo, NULL,
238  TREETABLE, 0, 0, 0, 0, 0, 0, 0);
239  init_field("MTRI", "TreeID", MEM_INT64, FILE_INT64, READ_IF_PRESENT, 1, A_TT, &MergerTree->TreeTable[0].TreeID, NULL, TREETABLE, 0,
240  0, 0, 0, 0, 0, 0);
241 }
242 #endif
243 
244 void snap_io::read_snapshot(int num, mysnaptype loc_snap_type)
245 {
246  snap_type = loc_snap_type;
247 
248  char buf[MAXLEN_PATH_EXTRA];
249 
250  if(snap_type == MOST_BOUND_PARTICLE_SNAPHOT)
251  {
252  if(All.NumFilesPerSnapshot > 1)
253  sprintf(buf, "%s/snapdir_%03d/%s-prevmostboundonly_%03d", All.OutputDir, num, All.SnapshotFileBase, num);
254  else
255  sprintf(buf, "%s%s-prevmostboundonly_%03d", All.OutputDir, All.SnapshotFileBase, num);
256  }
257  else
258  {
259  if(All.NumFilesPerSnapshot > 1)
260  sprintf(buf, "%s/snapdir_%03d/%s_%03d", All.OutputDir, num, All.SnapshotFileBase, num);
261  else
262  sprintf(buf, "%s%s_%03d", All.OutputDir, All.SnapshotFileBase, num);
263  }
264 
265  read_ic(buf);
266 }
267 
289 void snap_io::read_ic(const char *fname)
290 {
291  if(All.ICFormat < 1 || All.ICFormat > 4)
292  Terminate("ICFormat = %d not supported.\n", All.ICFormat);
293 
294  TIMER_START(CPU_SNAPSHOT);
295 
296  double t0 = Logs.second();
297 
298  Sp->TotNumPart = 0;
299 
300  int num_files = find_files(fname, fname);
301 
303 
304  /* we repeat reading the headers of the files two times. In the first iteration, only the
305  * particle numbers ending up on each processor are assembled, followed by memory allocation.
306  * In the second iteration, the data is actually read in.
307  */
308  for(int rep = 0; rep < 2; rep++)
309  {
310  Sp->NumPart = 0;
311  Sp->NumGas = 0;
312 
313  read_files_driver(fname, rep, num_files);
314 
315  /* now do the memory allocation */
316  if(rep == 0)
317  {
318  int max_load, max_sphload;
319  MPI_Allreduce(&Sp->NumPart, &max_load, 1, MPI_INT, MPI_MAX, Communicator);
320  MPI_Allreduce(&Sp->NumGas, &max_sphload, 1, MPI_INT, MPI_MAX, Communicator);
321 
322 #ifdef TILING
323  max_load *= TILING * TILING * TILING;
324  max_sphload *= TILING * TILING * TILING;
325 #endif
326 
327  Sp->MaxPart = max_load / (1.0 - 2 * ALLOC_TOLERANCE);
328  Sp->MaxPartSph = max_sphload / (1.0 - 2 * ALLOC_TOLERANCE);
329 
330  Sp->allocate_memory(); /* allocate particle storage */
331 
332 #ifndef OUTPUT_COORDINATES_AS_INTEGERS
333  Ptmp = (ptmp_data *)Mem.mymalloc_movable(&Ptmp, "Ptmp", Sp->NumPart * sizeof(ptmp_data));
334 #endif
335  }
336  }
337 
338  MPI_Barrier(Communicator);
339 
340  long long byte_count = get_io_byte_count(), byte_count_all;
341  sumup_longs(1, &byte_count, &byte_count_all, Communicator);
342 
343  double t1 = Logs.second();
344 
345  mpi_printf("READIC: reading done. Took %g sec, total size %g MB, corresponds to effective I/O rate of %g MB/sec\n",
346  Logs.timediff(t0, t1), byte_count_all / (1024.0 * 1024.0), byte_count_all / (1024.0 * 1024.0) / Logs.timediff(t0, t1));
347 
348  mpi_printf("\nREADIC: Total number of particles : %lld\n\n", Sp->TotNumPart);
349 
350  snap_init_domain_mapping();
351 
352 #ifndef OUTPUT_COORDINATES_AS_INTEGERS
353  Mem.myfree(Ptmp);
354 #endif
355 
356 #ifdef TILING
357 
358  MyIntPosType halfboxlen = ((MyIntPosType)1) << ((BITS_FOR_POSITIONS - 1));
359 
360  MyIntPosType len = (halfboxlen / TILING) + (halfboxlen / TILING);
361 
362  Sp->TotNumGas *= TILING * TILING * TILING;
363  Sp->TotNumPart *= TILING * TILING * TILING;
364 
365  int add_numgas = Sp->NumGas * (TILING * TILING * TILING - 1);
366  /* create a gap behind the existing gas particles where we will insert the new gas particles */
367  memmove(static_cast<void *>(Sp->P + Sp->NumGas + add_numgas), static_cast<void *>(Sp->P + Sp->NumGas),
368  (Sp->NumPart - Sp->NumGas) * sizeof(simparticles::pdata));
369 
370  int off = 0;
371 
372  for(int i = 0; i < TILING; i++)
373  for(int j = 0; j < TILING; j++)
374  for(int k = 0; k < TILING; k++)
375  if(i != 0 || j != 0 || k != 0)
376  for(int n = 0; n < Sp->NumGas; n++)
377  {
378  Sp->P[Sp->NumGas + off] = Sp->P[n];
379 
380  Sp->P[Sp->NumGas + off].IntPos[0] += i * len;
381  Sp->P[Sp->NumGas + off].IntPos[1] += j * len;
382  Sp->P[Sp->NumGas + off].IntPos[2] += k * len;
383 
384  off++;
385  }
386 
387  if(off != add_numgas)
388  Terminate("this should not happen");
389 
390  Sp->NumGas += add_numgas;
391  Sp->NumPart += add_numgas;
392 
393  int add_numpart = (Sp->NumPart - Sp->NumGas) * (TILING * TILING * TILING - 1);
394 
395  off = 0;
396 
397  for(int i = 0; i < TILING; i++)
398  for(int j = 0; j < TILING; j++)
399  for(int k = 0; k < TILING; k++)
400  if(i != 0 || j != 0 || k != 0)
401  for(int n = Sp->NumGas; n < Sp->NumPart; n++)
402  {
403  Sp->P[Sp->NumPart + off] = Sp->P[n];
404 
405  Sp->P[Sp->NumPart + off].IntPos[0] += i * len;
406  Sp->P[Sp->NumPart + off].IntPos[1] += j * len;
407  Sp->P[Sp->NumPart + off].IntPos[2] += k * len;
408 
409  off++;
410  }
411 
412  if(off != add_numpart)
413  Terminate("this should not happen");
414 
415  Sp->NumPart += add_numpart;
416 #endif
417 
418 #ifdef GADGET2_HEADER
419 #ifndef INITIAL_CONDITIONS_CONTAIN_ENTROPY
420  if(header.flag_entropy_instead_u)
421  Terminate("Initial condition file contains entropy, but INITIAL_CONDITIONS_CONTAIN_ENTROPY is not set\n");
422 #else
423  if(!header.flag_entropy_instead_u)
424  Terminate("Initial condition file contains uthermal, but INITIAL_CONDITIONS_CONTAIN_ENTROPY is set\n");
425 #endif
426 #endif
427 
428  TIMER_STOP(CPU_SNAPSHOT);
429 }
430 
434 void snap_io::snap_init_domain_mapping(void)
435 {
436 #ifdef PERIODIC
437 
438  Sp->RegionLen = All.BoxSize;
439  Sp->FacCoordToInt = pow(2.0, BITS_FOR_POSITIONS) / Sp->RegionLen;
440  Sp->FacIntToCoord = Sp->RegionLen / pow(2.0, BITS_FOR_POSITIONS);
441 
442 #else
443 
444  double posmin[3], posmax[3];
445  for(int k = 0; k < 3; k++)
446  {
447  posmin[k] = MAX_REAL_NUMBER;
448  posmax[k] = -MAX_REAL_NUMBER;
449  }
450 
451  for(int i = 0; i < Sp->NumPart; i++)
452  for(int k = 0; k < 3; k++)
453  {
454  if(Ptmp[i].Pos[k] < posmin[k])
455  posmin[k] = Ptmp[i].Pos[k];
456 
457  if(Ptmp[i].Pos[k] > posmax[k])
458  posmax[k] = Ptmp[i].Pos[k];
459  }
460 
461  double xyz[6] = {posmin[0], posmin[1], posmin[2], -posmax[0], -posmax[1], -posmax[2]};
462  double xyz_glob[6];
463 
464  MPI_Allreduce(xyz, xyz_glob, 6, MPI_DOUBLE, MPI_MIN, Communicator);
465 
466  mpi_printf("READIC: Region covered with particles: (%g %g %g) -> (%g %g %g)\n", xyz_glob[0], xyz_glob[1], xyz_glob[2], -xyz_glob[3],
467  -xyz_glob[4], -xyz_glob[5]);
468 
469  Sp->RegionLen = 0;
470  for(int j = 0; j < 3; j++)
471  if(-xyz_glob[j + 3] - xyz_glob[j] > Sp->RegionLen)
472  Sp->RegionLen = -xyz_glob[j + 3] - xyz_glob[j];
473 
474  Sp->RegionLen *= 4.0;
475 
476  mpi_printf("READIC: Initial root node size: %g\n", Sp->RegionLen);
477 
478  for(int j = 0; j < 3; j++)
479  {
480  Sp->RegionCenter[j] = 0.5 * (xyz_glob[j] - xyz_glob[j + 3]);
481  Sp->RegionCorner[j] = Sp->RegionCenter[j] - 0.5 * Sp->RegionLen;
482  }
483 
484  Sp->FacCoordToInt = pow(2.0, BITS_FOR_POSITIONS) / Sp->RegionLen;
485  Sp->FacIntToCoord = Sp->RegionLen / pow(2.0, BITS_FOR_POSITIONS);
486 
487 #endif
488 
489 #ifndef OUTPUT_COORDINATES_AS_INTEGERS
490  for(int i = 0; i < Sp->NumPart; i++)
491  Sp->pos_to_intpos(Ptmp[i].Pos, Sp->P[i].IntPos); /* converts floating point representation to integers */
492 #endif
493 }
494 
496 {
497  if(snap_type == NORMAL_SNAPSHOT)
498  {
499  return Sp->P[index].getType();
500  }
501  else if(snap_type == MOST_BOUND_PARTICLE_SNAPHOT)
502  {
503  if(Sp->P[index].ID.is_previously_most_bound())
504  return Sp->P[index].getType();
505  else
506  return -1; // marks that particle has type that is not written out
507  }
508  else if(snap_type == MOST_BOUND_PARTICLE_SNAPHOT_REORDERED)
509  {
510  return Sp->P[index].getType();
511  }
512  else
513  Terminate("can't be");
514 }
515 
516 void snap_io::set_type_of_element(int index, int type)
517 {
518  if(type < NTYPES)
519  Sp->P[index].setType(type);
520 }
521 
522 void *snap_io::get_base_address_of_structure(enum arrays array, int index)
523 {
524  switch(array)
525  {
526  case A_SPHP:
527  return (void *)(Sp->SphP + index);
528  case A_P:
529  return (void *)(Sp->P + index);
530  case A_PS:
531  return (void *)(Sp->PS + index);
532 #ifdef LIGHTCONE_PARTICLES
533  case A_LC:
534  Terminate("a"); // return (void *) (Lp->P + index);
535 #endif
536 #ifdef LIGHTCONE_MASSMAPS
537  case A_MM:
538  Terminate("b"); // return (void *) (Lp->P + index);
539 #endif
540 #if defined(REARRANGE_OPTION) && defined(MERGERTREE)
541  case A_TT:
542  return (void *)(MergerTree->TreeTable + index);
543 #endif
544  default:
545  Terminate("we don't expect to get here");
546  }
547 
548  return NULL;
549 }
550 
559 void snap_io::write_snapshot(int num, mysnaptype loc_snap_type)
560 {
561 #ifdef DO_NOT_PRODUCE_BIG_OUTPUT
562  if(snap_type != MOST_BOUND_PARTICLE_SNAPHOT)
563  {
564  mpi_printf("\nSNAPSHOT: We skip writing snapshot file #%d @ time %g\n", num, All.Time);
565  return;
566  }
567 #endif
568 
569  snap_type = loc_snap_type;
570 
571  TIMER_START(CPU_SNAPSHOT);
572 
573  mpi_printf("\nSNAPSHOT: writing snapshot file #%d @ time %g ... \n", num, All.Time);
574 
575  double t0 = Logs.second();
577 
579 
580  int n_type[NTYPES];
581 
582  /* determine global and local particle numbers */
583  for(int n = 0; n < NTYPES; n++)
584  n_type[n] = 0;
585 
586  for(int n = 0; n < Sp->NumPart; n++)
587  {
588  int type = get_type_of_element(n);
589  if(type >= 0)
590  n_type[type]++;
591  }
592 
593  sumup_large_ints(NTYPES, n_type, ntot_type_all, Communicator);
594 
595  if(All.NumFilesPerSnapshot > 1)
596  {
597  if(ThisTask == 0)
598  {
599  char buf[MAXLEN_PATH_EXTRA];
600  sprintf(buf, "%s/snapdir_%03d", All.OutputDir, num);
601  mkdir(buf, 02755);
602  }
603  MPI_Barrier(Communicator);
604  }
605 
606  char buf[MAXLEN_PATH_EXTRA];
607  if(All.NumFilesPerSnapshot > 1)
608  sprintf(buf, "%s/snapdir_%03d/%s_%03d", All.OutputDir, num, All.SnapshotFileBase, num);
609  else
610  sprintf(buf, "%s%s_%03d", All.OutputDir, All.SnapshotFileBase, num);
611 
612  if(snap_type == MOST_BOUND_PARTICLE_SNAPHOT)
613  {
614  if(All.NumFilesPerSnapshot > 1)
615  sprintf(buf, "%s/snapdir_%03d/%s-prevmostboundonly_%03d", All.OutputDir, num, All.SnapshotFileBase, num);
616  else
617  sprintf(buf, "%s%s-prevmostboundonly_%03d", All.OutputDir, All.SnapshotFileBase, num);
618  }
619  else if(snap_type == MOST_BOUND_PARTICLE_SNAPHOT_REORDERED)
620  {
621  if(All.NumFilesPerSnapshot > 1)
622  sprintf(buf, "%s/snapdir_%03d/%s-prevmostboundonly-treeorder_%03d", All.OutputDir, num, All.SnapshotFileBase, num);
623  else
624  sprintf(buf, "%s%s-prevmostboundonly-treeorder_%03d", All.OutputDir, All.SnapshotFileBase, num);
625  }
626 
627  /* now write the files */
629 
630  long long byte_count = get_io_byte_count(), byte_count_all;
631  sumup_longs(1, &byte_count, &byte_count_all, Communicator);
632 
633  double t1 = Logs.second();
634 
635  mpi_printf("SNAPSHOT: done with writing snapshot. Took %g sec, total size %g MB, corresponds to effective I/O rate of %g MB/sec\n",
636  Logs.timediff(t0, t1), byte_count_all / (1024.0 * 1024.0), byte_count_all / (1024.0 * 1024.0) / Logs.timediff(t0, t1));
637 
639 
640  TIMER_STOP(CPU_SNAPSHOT);
641 }
642 
643 void snap_io::fill_file_header(int writeTask, int lastTask, long long *n_type, long long *ntot_type)
644 {
645  /* determine global and local particle numbers */
646  for(int n = 0; n < NTYPES + 1; n++)
647  n_type[n] = 0;
648 
649  for(int n = 0; n < Sp->NumPart; n++)
650  {
651  int type = get_type_of_element(n);
652  if(type >= 0)
653  n_type[type]++;
654  }
655 
656 #if defined(REARRANGE_OPTION) && defined(MERGERTREE)
658  n_type[NTYPES] = MergerTree->Ntrees;
659 #endif
660 
661  /* determine particle numbers of each type in file */
662  if(ThisTask == writeTask)
663  {
664  for(int n = 0; n < NTYPES + 1; n++)
665  ntot_type[n] = n_type[n];
666 
667  for(int task = writeTask + 1; task <= lastTask; task++)
668  {
669  long long nn[NTYPES + 1];
670  MPI_Recv(&nn[0], NTYPES + 1, MPI_LONG_LONG, task, TAG_LOCALN, Communicator, MPI_STATUS_IGNORE);
671  for(int n = 0; n < NTYPES + 1; n++)
672  ntot_type[n] += nn[n];
673  }
674 
675  for(int task = writeTask + 1; task <= lastTask; task++)
676  MPI_Send(&ntot_type[0], NTYPES + 1, MPI_LONG_LONG, task, TAG_N, Communicator);
677  }
678  else
679  {
680  MPI_Send(&n_type[0], NTYPES + 1, MPI_LONG_LONG, writeTask, TAG_LOCALN, Communicator);
681  MPI_Recv(&ntot_type[0], NTYPES + 1, MPI_LONG_LONG, writeTask, TAG_N, Communicator, MPI_STATUS_IGNORE);
682  }
683 
684  /* fill file header */
685 
686  for(int n = 0; n < NTYPES; n++)
687  {
688  header.npart[n] = ntot_type[n];
689  header.npartTotal[n] = ntot_type_all[n];
690  }
691 
692 #ifdef MERGERTREE
694  {
695  header.Ntrees = ntot_type[NTYPES];
696 #if defined(REARRANGE_OPTION) && defined(MERGERTREE)
697  header.TotNtrees = MergerTree->TotNtrees;
698 #endif
699  }
700  else
701  {
702  header.Ntrees = 0;
703  header.TotNtrees = 0;
704  }
705 #endif
706 
707  for(int n = 0; n < NTYPES; n++)
708  header.mass[n] = All.MassTable[n];
709 
710  header.time = All.Time;
711 
713  header.redshift = 1.0 / All.Time - 1;
714  else
715  header.redshift = 0;
716 
719 
720 #ifdef GADGET2_HEADER
721  for(int n = 0; n < NTYPES; n++)
722  header.npartTotalLowWord[n] = ntot_type_all[n];
723 
724  header.flag_sfr = 0;
725  header.flag_feedback = 0;
726  header.flag_cooling = 0;
727 
728 #ifdef COOLING
729  header.flag_cooling = 1;
730 #endif
731 
732 #ifdef STARFORMATION
733  header.flag_sfr = 1;
734  header.flag_feedback = 1;
735 #endif
736  header.Omega0 = All.Omega0;
737  header.OmegaLambda = All.OmegaLambda;
738 
739 #if !(defined(REARRANGE_OPTION) && defined(MERGERTREE))
740  header.HubbleParam = All.HubbleParam;
741  header.Hubble = All.Hubble;
742 #endif
743 
744 #ifdef OUTPUT_IN_DOUBLEPRECISION
745  header.flag_doubleprecision = 1;
746 #else
747  header.flag_doubleprecision = 0;
748 #endif
749 #endif
750 }
751 
752 void snap_io::read_file_header(const char *fname, int filenr, int readTask, int lastTask, long long *n_type, long long *ntot_type,
753  int *nstart)
754 {
755  n_type[NTYPES] = 0;
756  ntot_type[NTYPES] = 0;
757 
758 #ifdef GADGET2_HEADER
759  for(int i = 0; i < NTYPES_HEADER; i++)
760  if(header.npartTotalLowWord[i] > 0)
761  header.npartTotal[i] = header.npartTotalLowWord[i]; // + (((long long)header.npartTotalHighWord[i]) << 32);
762 #endif
763 
764  if(Sp->TotNumPart == 0)
765  {
766  if(header.num_files <= 1)
767  for(int i = 0; i < NTYPES; i++)
768  header.npartTotal[i] = header.npart[i];
769 
770  Sp->TotNumGas = header.npartTotal[0];
771  Sp->TotNumPart = 0;
772 
773  for(int i = 0; i < NTYPES; i++)
774  Sp->TotNumPart += header.npartTotal[i];
775 
776 #ifdef GADGET2_HEADER
777 #if defined(SECOND_ORDER_LPT_ICS)
778  if(header.flag_ic_info == FLAG_SECOND_ORDER_ICS)
779  {
780  mpi_printf("READIC: Second order ICs detected. Will complete them first before starting run.\n");
781  All.LptScalingfactor = header.lpt_scalingfactor;
782  }
783  else
784  Terminate("READIC: No second order ICs detected even though you activated SECOND_ORDER_LPT_ICS.\n");
785 #else
786  if(header.flag_ic_info == FLAG_SECOND_ORDER_ICS)
787  Terminate("Detected second order ICs but SECOND_ORDER_LPT_ICS is not enabled.\n");
788 #endif
789 #endif
790 
791 #ifdef GENERATE_GAS_IN_ICS
792  if(All.RestartFlag == RST_BEGIN)
793  {
794 #ifdef SPLIT_PARTICLE_TYPE
795  for(int i = 0; i < NTYPES; i++)
796  if((1 << i) & (SPLIT_PARTICLE_TYPE))
797  {
798  Sp->TotNumGas += header.npartTotal[i];
799  Sp->TotNumPart += header.npartTotal[i];
800  }
801 #else
802  Sp->TotNumGas += header.npartTotal[1];
803  Sp->TotNumPart += header.npartTotal[1];
804 #endif
805  }
806 #endif
807 
808  for(int i = 0; i < NTYPES; i++)
809  All.MassTable[i] = header.mass[i];
810 
812  All.Time = All.TimeBegin;
813  else
815 
817  }
818 
819  if(ThisTask == readTask)
820  {
821  if(filenr == 0 && nstart == NULL)
822  {
823  mpi_printf(
824  "\nREADIC: filenr=%d, '%s' contains:\n"
825  "READIC: Type 0 (gas): %8lld (tot=%15lld) masstab= %g\n",
826  filenr, fname, (long long)header.npart[0], (long long)header.npartTotal[0], All.MassTable[0]);
827 
828  for(int type = 1; type < NTYPES; type++)
829  {
830  mpi_printf("READIC: Type %d: %8lld (tot=%15lld) masstab= %g\n", type, (long long)header.npart[type],
831  (long long)header.npartTotal[type], All.MassTable[type]);
832  }
833  mpi_printf("\n");
834  }
835  }
836 
837  /* to collect the gas particles all at the beginning (in case several
838  snapshot files are read on the current CPU) we move the collisionless
839  particles such that a gap of the right size is created */
840 
841  long long nall = 0;
842  for(int type = 0; type < NTYPES; type++)
843  {
844  ntot_type[type] = header.npart[type];
845 
846  long long n_in_file = header.npart[type];
847  int ntask = lastTask - readTask + 1;
848  int n_for_this_task = n_in_file / ntask;
849  if((ThisTask - readTask) < (n_in_file % ntask))
850  n_for_this_task++;
851 
852  n_type[type] = n_for_this_task;
853 
854  nall += n_for_this_task;
855  }
856 
857 #ifdef MERGERTREE
859  ntot_type[NTYPES] = header.Ntrees;
860 #endif
861 
862  if(nstart)
863  {
864  memmove(static_cast<void *>(&Sp->P[Sp->NumGas + nall]), static_cast<void *>(&Sp->P[Sp->NumGas]),
865  (Sp->NumPart - Sp->NumGas) * sizeof(particle_data));
866 #ifndef OUTPUT_COORDINATES_AS_INTEGERS
867  memmove(&Ptmp[Sp->NumGas + nall], &Ptmp[Sp->NumGas], (Sp->NumPart - Sp->NumGas) * sizeof(ptmp_data));
868 #endif
869  *nstart = Sp->NumGas;
870  }
871 }
872 
881 {
882 #ifdef GADGET2_HEADER
883  write_vector_attribute(handle, "NumPart_ThisFile", header.npart, H5T_NATIVE_UINT, NTYPES);
884 #else
885  write_vector_attribute(handle, "NumPart_ThisFile", header.npart, H5T_NATIVE_UINT64, NTYPES);
886 #endif
887  write_vector_attribute(handle, "NumPart_Total", header.npartTotal, H5T_NATIVE_UINT64, NTYPES);
888  write_vector_attribute(handle, "MassTable", header.mass, H5T_NATIVE_DOUBLE, NTYPES);
889  write_scalar_attribute(handle, "Time", &header.time, H5T_NATIVE_DOUBLE);
890  write_scalar_attribute(handle, "Redshift", &header.redshift, H5T_NATIVE_DOUBLE);
891  write_scalar_attribute(handle, "BoxSize", &header.BoxSize, H5T_NATIVE_DOUBLE);
892  write_scalar_attribute(handle, "NumFilesPerSnapshot", &header.num_files, H5T_NATIVE_INT);
893  write_string_attribute(handle, "Git_commit", GIT_COMMIT);
894  write_string_attribute(handle, "Git_date", GIT_DATE);
895 
896 #if defined(REARRANGE_OPTION) && defined(MERGERTREE)
898  {
899  write_scalar_attribute(handle, "Ntrees_ThisFile", &header.Ntrees, H5T_NATIVE_UINT64);
900  write_scalar_attribute(handle, "Ntrees_Total", &header.TotNtrees, H5T_NATIVE_UINT64);
901  }
902 #endif
903 }
904 
909 void snap_io::read_header_fields(const char *fname)
910 {
911  for(int i = 0; i < NTYPES; i++)
912  {
913  header.npart[i] = 0;
914  header.npartTotal[i] = 0;
915  header.mass[i] = 0;
916  }
917 
918  int ntypes = NTYPES;
919 
920  hid_t hdf5_file = my_H5Fopen(fname, H5F_ACC_RDONLY, H5P_DEFAULT);
921  hid_t handle = my_H5Gopen(hdf5_file, "/Header");
922 
923  /* check if the file in question actually has this number of types */
924  hid_t hdf5_attribute = my_H5Aopen_name(handle, "NumPart_ThisFile");
925  hid_t space = H5Aget_space(hdf5_attribute);
926  hsize_t dims, len;
927  H5Sget_simple_extent_dims(space, &dims, &len);
928  H5Sclose(space);
929  if(len != (size_t)ntypes)
930  Terminate("Length of NumPart_ThisFile attribute (%d) does not match NTYPES(ICS) (%d)", (int)len, ntypes);
931  my_H5Aclose(hdf5_attribute, "NumPart_ThisFile");
932 
933  /* now read the header fields */
934 
935 #ifdef GADGET2_HEADER
936  read_vector_attribute(handle, "NumPart_ThisFile", header.npart, H5T_NATIVE_UINT, ntypes);
937 #else
938  read_vector_attribute(handle, "NumPart_ThisFile", header.npart, H5T_NATIVE_UINT64, ntypes);
939 #endif
940 
941  read_vector_attribute(handle, "NumPart_Total", header.npartTotal, H5T_NATIVE_UINT64, ntypes);
942 
943 #ifdef MERGERTREE
945  {
946  read_scalar_attribute(handle, "Ntrees_ThisFile", &header.Ntrees, H5T_NATIVE_UINT64);
947  read_scalar_attribute(handle, "Ntrees_Total", &header.TotNtrees, H5T_NATIVE_UINT64);
948  }
949 #endif
950 
951  read_vector_attribute(handle, "MassTable", header.mass, H5T_NATIVE_DOUBLE, ntypes);
952  read_scalar_attribute(handle, "Time", &header.time, H5T_NATIVE_DOUBLE);
953  read_scalar_attribute(handle, "Redshift", &header.redshift, H5T_NATIVE_DOUBLE);
954  read_scalar_attribute(handle, "BoxSize", &header.BoxSize, H5T_NATIVE_DOUBLE);
955  read_scalar_attribute(handle, "NumFilesPerSnapshot", &header.num_files, H5T_NATIVE_INT);
956 
957  my_H5Gclose(handle, "/Header");
958  my_H5Fclose(hdf5_file, fname);
959 }
960 
962 
963 void snap_io::set_filenr_in_header(int numfiles) { header.num_files = numfiles; }
964 
965 void snap_io::read_increase_numbers(int type, int n_for_this_task)
966 {
967  if(type < NTYPES)
968  Sp->NumPart += n_for_this_task;
969 
970  if(type == 0)
971  Sp->NumGas += n_for_this_task;
972 }
973 
974 void snap_io::get_datagroup_name(int type, char *buf)
975 {
976  if(type < NTYPES)
977  sprintf(buf, "/PartType%d", type);
978  else if(type == NTYPES)
979  sprintf(buf, "/TreeTable");
980  else
981  Terminate("wrong group");
982 }
global_data_all_processes All
Definition: main.cc:40
int N_DataGroups
Definition: io.h:139
void write_multiple_files(char *fname, int numfilesperdump, int append_flag=0, int chunk_size=0)
Definition: io.cc:700
void * header_buf
Definition: io.h:174
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
char info[100]
Definition: io.h:177
enum file_contents type_of_file
Definition: io.h:187
void read_files_driver(const char *fname, int rep, int numfiles)
Definition: io.cc:227
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
bool is_previously_most_bound(void)
Definition: idstorage.h:117
MyReal FacCoordToInt
Definition: intposconvert.h:38
void pos_to_intpos(T *posdiff, MyIntPosType *intpos)
MyReal RegionLen
Definition: intposconvert.h:39
MyReal RegionCenter[3]
Definition: intposconvert.h:43
MyReal RegionCorner[3]
Definition: intposconvert.h:42
MyReal FacIntToCoord
Definition: intposconvert.h:37
long long byte_count
void reset_io_byte_count(void)
long long get_io_byte_count(void)
double timediff(double t0, double t1)
Definition: logs.cc:488
double second(void)
Definition: logs.cc:471
void mpi_printf(const char *fmt,...)
Definition: setcomm.h:55
int ThisTask
Definition: setcomm.h:33
MPI_Comm Communicator
Definition: setcomm.h:31
subfind_data * PS
Definition: simparticles.h:63
long long TotNumPart
Definition: simparticles.h:46
long long TotNumGas
Definition: simparticles.h:47
void allocate_memory(void)
Definition: simparticles.h:273
particle_data * P
Definition: simparticles.h:54
sph_particle_data * SphP
Definition: simparticles.h:59
int get_filenr_from_header(void)
Definition: snap_io.cc:961
void read_header_fields(const char *fname)
This function reads the snapshot header in case of hdf5 files (i.e. format 3)
Definition: snap_io.cc:909
void * get_base_address_of_structure(enum arrays array, int index)
Definition: snap_io.cc:522
void read_increase_numbers(int type, int n_for_this_task)
Definition: snap_io.cc:965
void read_snapshot(int num, mysnaptype snap_type)
Definition: snap_io.cc:244
void fill_file_header(int writeTask, int lastTask, long long *nloc_part, long long *npart)
Definition: snap_io.cc:643
void read_file_header(const char *fname, int filenr, int readTask, int lastTask, long long *nloc_part, long long *npart, int *nstart)
Definition: snap_io.cc:752
void set_type_of_element(int index, int type)
Definition: snap_io.cc:516
void write_header_fields(hid_t)
Write the fields contained in the header group of the HDF5 snapshot file.
Definition: snap_io.cc:880
void get_datagroup_name(int grnr, char *gname)
Definition: snap_io.cc:974
void read_ic(const char *fname)
This function reads initial conditions that are in on of the default file formats of Gadget.
Definition: snap_io.cc:289
void set_filenr_in_header(int)
Definition: snap_io.cc:963
void init_basic(simparticles *Sp_ptr)
Function for field registering.
Definition: snap_io.cc:50
io_header header
Definition: snap_io.h:125
int get_type_of_element(int index)
Definition: snap_io.cc:495
void write_snapshot(int num, mysnaptype snap_type)
Save snapshot to disk.
Definition: snap_io.cc:559
#define ALLOC_TOLERANCE
Definition: constants.h:74
#define SOLAR_MASS
Definition: constants.h:119
#define MAXLEN_PATH_EXTRA
Definition: constants.h:301
#define GAMMA
Definition: constants.h:99
#define SEC_PER_YEAR
Definition: constants.h:128
#define NTYPES
Definition: constants.h:308
#define MAX_REAL_NUMBER
Definition: constants.h:94
mysnaptype
Definition: dtypes.h:305
@ MOST_BOUND_PARTICLE_SNAPHOT
Definition: dtypes.h:307
@ NORMAL_SNAPSHOT
Definition: dtypes.h:306
@ MOST_BOUND_PARTICLE_SNAPHOT_REORDERED
Definition: dtypes.h:308
uint32_t MyIntPosType
Definition: dtypes.h:35
@ RST_CREATEICS
Definition: dtypes.h:319
@ RST_FOF
Definition: dtypes.h:316
@ RST_BEGIN
Definition: dtypes.h:313
@ RST_RESUME
Definition: dtypes.h:314
#define BITS_FOR_POSITIONS
Definition: dtypes.h:37
void write_vector_attribute(hid_t handle, const char *attr_name, const void *buf, hid_t mem_type_id, int length)
Definition: hdf5_util.cc:636
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 write_scalar_attribute(hid_t handle, const char *attr_name, const void *buf, hid_t mem_type_id)
Definition: hdf5_util.cc:621
void read_scalar_attribute(hid_t handle, const char *attr_name, void *buf, hid_t mem_type_id)
Definition: hdf5_util.cc:589
void write_string_attribute(hid_t handle, const char *attr_name, const char *buf)
Definition: hdf5_util.cc:654
@ READ_IF_PRESENT
Definition: io.h:124
@ SKIP_ON_READ
Definition: io.h:125
@ ALL_TYPES
Definition: io.h:103
@ AGE_BLOCK
Definition: io.h:105
@ Z_BLOCK
Definition: io.h:106
@ MASS_BLOCK
Definition: io.h:104
@ GAS_ONLY
Definition: io.h:100
@ TREETABLE
Definition: io.h:116
#define FLAG_SECOND_ORDER_ICS
Definition: io.h:272
arrays
Definition: io.h:30
@ A_NONE
Definition: io.h:31
@ A_MM
Definition: io.h:46
@ A_SPHP
Definition: io.h:32
@ A_P
Definition: io.h:33
@ A_TT
Definition: io.h:42
@ A_PS
Definition: io.h:34
@ A_LC
Definition: io.h:45
@ FILE_MY_INTPOS_TYPE
Definition: io.h:72
@ 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 NTYPES_HEADER
Definition: io.h:268
@ MEM_INT
Definition: io.h:79
@ MEM_MY_FLOAT
Definition: io.h:85
@ 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
logs Logs
Definition: main.cc:43
#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
#define TAG_LOCALN
Definition: mpi_utils.h:52
void sumup_large_ints(int n, int *src, long long *res, MPI_Comm comm)
#define TAG_N
Definition: mpi_utils.h:25
void sumup_longs(int n, long long *src, long long *res, MPI_Comm comm)
memory Mem
Definition: main.cc:44
expr pow(half base, half exp)
Definition: half.hpp:2803
double UnitVelocity_in_cm_per_s
Definition: allvars.h:119
enum restart_options RestartFlag
Definition: allvars.h:68
double MassTable[NTYPES]
Definition: allvars.h:269
integertime Ti_Current
Definition: allvars.h:188
char OutputDir[MAXLEN_PATH]
Definition: allvars.h:272
integertime Ti_lastoutput
Definition: allvars.h:190
char SnapshotFileBase[MAXLEN_PATH]
Definition: allvars.h:272
void set_cosmo_factors_for_current_time(void)
Definition: allvars.cc:20
MyIDStorage ID
Definition: particle_data.h:70
unsigned char getType(void)
void setType(unsigned char type)
MyIntPosType IntPos[3]
Definition: particle_data.h:53
long long TotNtrees
Definition: snap_io.h:120
double mass[NTYPES_HEADER]
Definition: snap_io.h:112
long long npartTotal[NTYPES_HEADER]
Definition: snap_io.h:110
long long npart[NTYPES_HEADER]
Definition: snap_io.h:109
long long Ntrees
Definition: snap_io.h:119
const char * GIT_COMMIT
const char * GIT_DATE