GADGET-4
simparticles.h
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 #ifndef SIMPART_H
13 #define SIMPART_H
14 
15 #include <math.h>
16 
17 #include "../data/allvars.h"
18 #include "../data/constants.h"
19 #include "../data/dtypes.h"
20 #include "../data/intposconvert.h"
21 #include "../data/macros.h"
22 #include "../data/mymalloc.h"
23 #include "../data/particle_data.h"
24 #include "../data/sph_particle_data.h"
25 #include "../main/main.h"
26 #include "../mpi_utils/mpi_utils.h"
27 #include "../mpi_utils/setcomm.h"
28 #include "../system/system.h"
29 #include "../time_integration/timestep.h"
30 
31 #ifdef LIGHTCONE
32 class lightcone;
33 #endif
34 
35 class simparticles : public intposconvert, public setcomm
36 {
37  public:
38  simparticles(MPI_Comm comm) : setcomm(comm) {}
39 
40  int NumPart;
41  int NumGas;
43  int MaxPart;
44  int MaxPartSph;
46  long long TotNumPart;
47  long long TotNumGas;
50 
56  /* the following struture holds data that is stored for each SPH particle in addition to the collisionless
57  * variables.
58  */
61  unsigned short int MarkerValue;
62 
64 
65  inline void copy_particle(particle_data *Ptarget, particle_data *Psource)
66  {
67  // we do this ugly trick here because the atomic_flag in particle_data has an implicitly deleted copy operator...
68  // but we know what we are doing, and this is the easiest way at the moment to work around this in our case unnecessary protection
69  memcpy(static_cast<void *>(Ptarget), static_cast<void *>(Psource), sizeof(particle_data));
70  }
71 
72  static bool inline compare_IDs(const MyIDType &a, const MyIDType &b) { return a < b; }
73 
74 #if defined(LIGHTCONE_PARTICLES_GROUPS) && defined(FOF)
75  double *DistanceOrigin;
76 #endif
77 
78 #ifdef SUBFIND_ORPHAN_TREATMENT
79  idstoredata IdStore;
80  static inline bool compare_SpP_ID(const particle_data &a, const particle_data &b) { return a.ID.get() < b.ID.get(); }
81 #endif
82 
83 #ifdef LIGHTCONE
84  lightcone *LightCone;
85 #endif
86 
87 #ifdef FOF
88  MyIDStorage *MinID;
89  int *Len; // this is here enough in 32bit because only the group segments on the local processor are treated
90  int *Head, *Next, *Tail, *MinIDTask;
91  MyFloat *fof_nearest_distance;
92  MyFloat *fof_nearest_hsml;
93 
94  struct bit_flags
95  {
96  unsigned char Nonlocal : 2, MinIDChanged : 2, Marked : 2;
97  } * Flags;
98 
99  double LinkL;
100 
101  inline void link_two_particles(int target, int j)
102  {
103  if(Head[target] != Head[j]) /* only if not yet linked */
104  {
105  int p, s;
106  if(Len[Head[target]] > Len[Head[j]]) /* p group is longer */
107  {
108  p = target;
109  s = j;
110  }
111  else
112  {
113  p = j;
114  s = target;
115  }
116  Next[Tail[Head[p]]] = Head[s];
117 
118  Tail[Head[p]] = Tail[Head[s]];
119 
120  Len[Head[p]] += Len[Head[s]];
121 
122  if(MinID[Head[s]].get() < MinID[Head[p]].get())
123  {
124  MinID[Head[p]] = MinID[Head[s]];
125  MinIDTask[Head[p]] = MinIDTask[Head[s]];
126  }
127 
128  int ss = Head[s];
129  do
130  Head[ss] = Head[p];
131  while((ss = Next[ss]) >= 0);
132  }
133  }
134 
135 #ifdef SUBFIND
136  struct nearest_r2_data
137  {
138  double dist[2];
139  } * R2Loc;
140 
141 #endif
142 #endif
143 
144 #ifdef PMGRID
145  double Asmth[2], Rcut[2];
146 #endif
147 
148 #if defined(PMGRID) && (!defined(PERIODIC) || defined(PLACEHIGHRESREGION))
149  double TotalMeshSize[2]; /* this is in integer space but should be double here to protect against overflows */
150  MySignedIntPosType Corner[2][3];
151  MySignedIntPosType Xmintot[2][3], Xmaxtot[2][3];
152  MyIntPosType MeshSize[2][3];
153  MyIntPosType Left[2][3];
154  MyIntPosType OldMeshSize[2];
155  MyIntPosType ReferenceIntPos[2][3];
156 #endif
157 
158 #if defined(RANDOMIZE_DOMAINCENTER_TYPES) || defined(PLACEHIGHRESREGION)
159  MyIntPosType PlacingMask;
160  MyIntPosType PlacingBlocksize;
161 #endif
162 
163 #ifdef PLACEHIGHRESREGION
164  inline int check_high_res_overlap(MyIntPosType *center, MyIntPosType halflen)
165  {
166  MyIntPosType intleft[3] = {center[0] - halflen - ReferenceIntPos[HIGH_MESH][0],
167  center[1] - halflen - ReferenceIntPos[HIGH_MESH][1],
168  center[2] - halflen - ReferenceIntPos[HIGH_MESH][2]};
169 
170  MyIntPosType intright[3] = {center[0] + halflen - ReferenceIntPos[HIGH_MESH][0],
171  center[1] + halflen - ReferenceIntPos[HIGH_MESH][1],
172  center[2] + halflen - ReferenceIntPos[HIGH_MESH][2]};
173 
174  MySignedIntPosType *left = (MySignedIntPosType *)intleft;
175  MySignedIntPosType *right = (MySignedIntPosType *)intright;
176 
177  if(right[0] <= Xmintot[HIGH_MESH][0] || left[0] >= Xmaxtot[HIGH_MESH][0] || right[1] <= Xmintot[HIGH_MESH][1] ||
178  left[1] >= Xmaxtot[HIGH_MESH][1] || right[2] <= Xmintot[HIGH_MESH][2] || left[2] >= Xmaxtot[HIGH_MESH][2])
179  return FLAG_OUTSIDE;
180  else if(right[0] <= Xmaxtot[HIGH_MESH][0] && left[0] >= Xmintot[HIGH_MESH][0] && right[1] <= Xmaxtot[HIGH_MESH][1] &&
181  left[1] >= Xmintot[HIGH_MESH][1] && right[2] <= Xmaxtot[HIGH_MESH][2] && left[2] >= Xmintot[HIGH_MESH][2])
182  return FLAG_INSIDE;
183  else
184  return FLAG_BOUNDARYOVERLAP;
185  }
186 
187  inline int check_high_res_point_location(MyIntPosType *intpos)
188  {
189  MyIntPosType relpos[3] = {intpos[0] - ReferenceIntPos[HIGH_MESH][0], intpos[1] - ReferenceIntPos[HIGH_MESH][1],
190  intpos[2] - ReferenceIntPos[HIGH_MESH][2]};
191 
192  MySignedIntPosType *pos = (MySignedIntPosType *)relpos;
193 
194  if(pos[0] < Xmintot[HIGH_MESH][0] || pos[0] >= Xmaxtot[HIGH_MESH][0] || pos[1] < Xmintot[HIGH_MESH][1] ||
195  pos[1] >= Xmaxtot[HIGH_MESH][1] || pos[2] < Xmintot[HIGH_MESH][2] || pos[2] >= Xmaxtot[HIGH_MESH][2])
196  return FLAG_OUTSIDE;
197  else
198  return FLAG_INSIDE;
199  }
200 
201 #endif
202 
206 
207  int nsource;
208  int *indexlist;
209 
210 #ifdef STARFORMATION
211  double TimeBinSfr[TIMEBINS];
212 #endif
213 
214  inline int getTimeBinSynchronized(int bin) { return TimeBinSynchronized[bin]; }
215 
216 #ifdef REARRANGE_OPTION
217  static bool compare_TreeID_ID(const particle_data &a, const particle_data &b)
218  {
219  if(a.TreeID < b.TreeID)
220  return true;
221 
222  if(a.TreeID > b.TreeID)
223  return false;
224 
225  return a.ID.get() < b.ID.get();
226  }
227 
228  static bool compare_ID(const particle_data &a, const particle_data &b) { return a.ID.get() < b.ID.get(); }
229 #endif
230 
231  inline MyFloat get_DtHsml(int i) { return SphP[i].DtHsml; }
232 
233  inline MyFloat get_Csnd(int i) { return SphP[i].Csnd; }
234 
235  inline MyFloat get_OldAcc(int i) { return P[i].OldAcc; }
236 
237  /* sets the internal energy per unit mass of particle i from its entropy */
238  inline double get_utherm_from_entropy(int i)
239  {
240 #ifdef ISOTHERM_EQS
241  return SphP[i].Entropy;
242 #else
243  double fact_entropy_to_u = pow(SphP[i].Density * All.cf_a3inv, GAMMA_MINUS1) / GAMMA_MINUS1;
244  return SphP[i].Entropy * fact_entropy_to_u;
245 #endif
246  }
247 
248  /* sets the entropy of particle i from its internal energy per unit mass */
249  inline void set_entropy_from_utherm(double utherm, int i)
250  {
251  double fact_u_to_entropy = GAMMA_MINUS1 / pow(SphP[i].Density * All.cf_a3inv, GAMMA_MINUS1);
252  SphP[i].Entropy = utherm * fact_u_to_entropy;
253  SphP[i].EntropyPred = SphP[i].Entropy;
254 
255 #ifdef PRESSURE_ENTROPY_SPH
256  SphP[i].EntropyToInvGammaPred = pow(SphP[i].EntropyPred, 1.0 / GAMMA);
257 #endif
258  }
259 
261  {
263 
264  for(int i = 0; i < NumPart; i++)
266  }
267 
268  /* This routine allocates memory for
269  * particle storage, both the collisionless and the SPH particles.
270  * The memory for the ordered binary tree of the timeline
271  * is also allocated.
272  */
273  void allocate_memory(void)
274  {
275  /* Note: P and SphP are initialized to zero */
276  P = (particle_data *)Mem.mymalloc_movable_clear(&P, "P", MaxPart * sizeof(particle_data));
277  SphP = (sph_particle_data *)Mem.mymalloc_movable_clear(&SphP, "SphP", MaxPartSph * sizeof(sph_particle_data));
278 
281  }
282 
283  void free_memory(void)
284  {
287 
288  Mem.myfree(SphP);
289  Mem.myfree(P);
290  }
291 
292  void reallocate_memory_maxpart(int maxpartNew)
293  {
294  mpi_printf("ALLOCATE: Changing to MaxPart = %d\n", maxpartNew);
295 
296  P = (particle_data *)Mem.myrealloc_movable(P, maxpartNew * sizeof(particle_data));
297  if(maxpartNew > MaxPart)
298  memset(((char *)P) + MaxPart * sizeof(particle_data), 0, (maxpartNew - MaxPart) * sizeof(particle_data));
299  MaxPart = maxpartNew;
300 
302  }
303 
304  void reallocate_memory_maxpartsph(int maxpartsphNew)
305  {
306  mpi_printf("ALLOCATE: Changing to MaxPartSph = %d\n", maxpartsphNew);
307 
308  SphP = (sph_particle_data *)Mem.myrealloc_movable(SphP, maxpartsphNew * sizeof(sph_particle_data));
309  if(maxpartsphNew > MaxPartSph)
310  memset(((char *)SphP) + MaxPartSph * sizeof(sph_particle_data), 0, (maxpartsphNew - MaxPartSph) * sizeof(sph_particle_data));
311  MaxPartSph = maxpartsphNew;
312 
314  }
315 
321  void dump_particles(void)
322  {
323  FILE *fd;
324  char buffer[200];
325  sprintf(buffer, "particles_%d.dat", ThisTask);
326  if((fd = fopen(buffer, "w")))
327  {
328  fwrite(&NumPart, 1, sizeof(int), fd);
329  for(int i = 0; i < NumPart; i++)
330  fwrite(&P[i].IntPos[0], 3, sizeof(MyIntPosType), fd);
331  for(int i = 0; i < NumPart; i++)
332  fwrite(&P[i].Vel[0], 3, sizeof(MyFloat), fd);
333  for(int i = 0; i < NumPart; i++)
334  fwrite(&P[i].ID, 1, sizeof(MyIDStorage), fd);
335  fclose(fd);
336  }
337  }
338 
344  {
345  MyReal pos[3];
346  intpos_to_pos(P[i].IntPos, pos); /* converts the integer coordinates to floating point */
347 
348  printf("Task=%d, ID=%llu, Type=%d, TimeBinGrav=%d, TimeBinHydro=%d, Mass=%g, pos=%g|%g|%g, vel=%g|%g|%g, OldAcc=%g\n", ThisTask,
349  (unsigned long long)P[i].ID.get(), P[i].getType(), P[i].TimeBinGrav, P[i].getTimeBinHydro(), P[i].getMass(), pos[0], pos[1],
350  pos[2], P[i].Vel[0], P[i].Vel[1], P[i].Vel[2], P[i].OldAcc);
351 #if defined(PMGRID) && defined(PERIODIC) && !defined(TREEPM_NOTIMESPLIT)
352  printf("GravAccel=%g|%g|%g, GravPM=%g|%g|%g, Soft=%g, SoftClass=%d\n", P[i].GravAccel[0], P[i].GravAccel[1], P[i].GravAccel[2],
353  P[i].GravPM[0], P[i].GravPM[1], P[i].GravPM[2], All.ForceSoftening[P[i].getSofteningClass()], P[i].getSofteningClass());
354 #else
355 #ifndef LEAN
356  printf("GravAccel=%g|%g|%g, Soft=%g, SoftType=%d\n", P[i].GravAccel[0], P[i].GravAccel[1], P[i].GravAccel[2],
358 #endif
359 #endif
360 
361  if(P[i].getType() == 0)
362  {
363  printf("rho=%g, hsml=%g, entr=%g, csnd=%g\n", SphP[i].Density, SphP[i].Hsml, SphP[i].Entropy, SphP[i].get_sound_speed());
364  printf("ID=%llu SphP[p].CurrentMaxTiStep=%g\n", (unsigned long long)P[i].ID.get(), SphP[i].CurrentMaxTiStep);
365  }
366 
367  myflush(stdout);
368  }
369 
375  {
376  for(int i = 0; i < NumPart; i++)
377  if(P[i].ID.get() == ID)
379  }
380 
381  public:
382  inline int get_active_index(int idx)
383  {
384 #ifdef HIERARCHICAL_GRAVITY
386 #else
387  return idx;
388 #endif
389  }
390 
391  void reconstruct_timebins(void);
393  void mark_active_timebins(void);
394  void drift_all_particles(void);
395  int drift_particle(particle_data *P, sph_particle_data *SphP, integertime time1, bool ignore_light_cone = false);
399 
400 #if defined(PMGRID) && !defined(TREEPM_NOTIMESPLIT)
401  integertime get_timestep_pm(void);
402 #endif
403 
404 #if defined(PMGRID) && defined(PERIODIC) && !defined(TREEPM_NOTIMESPLIT)
405  void find_long_range_step_constraint(void);
406 #endif
407 
408  void timebins_get_bin_and_do_validity_checks(integertime ti_step, int *bin_new, int bin_old);
409 
410  void assign_hydro_timesteps(void);
412 
413  int test_if_grav_timestep_is_too_large(int p, int bin);
414  int get_timestep_bin(integertime ti_step);
415 
416 #ifdef ADAPTIVE_HYDRO_SOFTENING
417  int get_softeningtype_for_hydro_particle(int i)
418  {
419  double soft = All.GasSoftFactor * SphP[i].Hsml;
420 
421  if(soft <= All.ForceSoftening[NSOFTCLASSES])
422  return NSOFTCLASSES;
423 
424  int k = 0.5 + log(soft / All.ForceSoftening[NSOFTCLASSES]) / log(All.AdaptiveHydroSofteningSpacing);
425  if(k >= NSOFTCLASSES_HYDRO)
426  k = NSOFTCLASSES_HYDRO - 1;
427 
428  return NSOFTCLASSES + k;
429  }
430 #endif
431 
432 #ifdef INDIVIDUAL_GRAVITY_SOFTENING
433 
434  public:
435 #if(INDIVIDUAL_GRAVITY_SOFTENING) & 2
436 #error "INDIVIDUAL_GRAVITY_SOFTENING may not include particle type 1 which is used as a reference point"
437 #endif
438 
439 #if((INDIVIDUAL_GRAVITY_SOFTENING)&1) && defined(ADAPTIVE_HYDRO_SOFTENING)
440 #error "INDIVIDUAL_GRAVITY_SOFTENING may not include particle type 0 when ADAPTIVE_HYDRO_SOFTENING is used"
441 #endif
442 
443  int get_softening_type_from_mass(double mass)
444  {
445  int min_type = -1;
446  double eps = get_desired_softening_from_mass(mass);
447  double min_dln = MAX_FLOAT_NUMBER;
448 
449  for(int i = 0; i < NSOFTCLASSES; i++)
450  {
451  if(All.ForceSoftening[i] > 0)
452  {
453  double dln = fabs(log(eps) - log(All.ForceSoftening[i]));
454 
455  if(dln < min_dln)
456  {
457  min_dln = dln;
458  min_type = i;
459  }
460  }
461  }
462 
463  if(min_type < 0)
464  Terminate("min_type < 0");
465 
466  return min_type;
467  }
468 
474  double get_desired_softening_from_mass(double mass)
475  {
476  return All.ForceSoftening[All.SofteningClassOfPartType[1]] * pow(mass / All.AvgType1Mass, 1.0 / 3);
477  }
478 
483  void init_individual_softenings(void)
484  {
485  int ndm = 0;
486  double mass = 0, masstot, massmin = MAX_DOUBLE_NUMBER, massmax = 0;
487  long long ndmtot;
488 
489  for(int i = 0; i < NumPart; i++)
490  if(P[i].getType() == 1)
491  {
492  ndm++;
493  mass += P[i].getMass();
494 
495  if(massmin > P[i].getMass())
496  massmin = P[i].getMass();
497 
498  if(massmax < P[i].getMass())
499  massmax = P[i].getMass();
500  }
501 
502  sumup_large_ints(1, &ndm, &ndmtot, Communicator);
503  MPI_Allreduce(&mass, &masstot, 1, MPI_DOUBLE, MPI_SUM, Communicator);
504 
505  MPI_Allreduce(MPI_IN_PLACE, &massmin, 1, MPI_DOUBLE, MPI_MIN, Communicator);
506  MPI_Allreduce(MPI_IN_PLACE, &massmax, 1, MPI_DOUBLE, MPI_MAX, Communicator);
507 
508  All.AvgType1Mass = masstot / ndmtot;
509 
510  mpi_printf("INIT: AvgType1Mass = %g (min=%g max=%g) Ndm1tot=%lld\n", All.AvgType1Mass, massmin, massmax, ndmtot);
511 
512  if(massmax > 1.00001 * massmin)
513  Terminate("Strange: Should use constant mass type-1 particles if INDIVIDUAL_GRAVITY_SOFTENING is used\n");
514 
516  {
517  double rhomean_dm = (All.Omega0 - All.OmegaBaryon) * (3 * All.Hubble * All.Hubble / (8 * M_PI * All.G));
518 
519  mpi_printf("INIT: For this AvgType1Mass, the mean particle spacing is %g and the assigned softening is %g\n",
520  pow(All.AvgType1Mass / rhomean_dm, 1.0 / 3), All.SofteningTable[All.SofteningClassOfPartType[1]]);
521  }
522 
523  for(int i = 0; i < NumPart; i++)
524  if(((1 << P[i].getType()) & (INDIVIDUAL_GRAVITY_SOFTENING)))
525  P[i].setSofteningClass(get_softening_type_from_mass(P[i].getMass()));
526  }
527 #endif
528 };
529 
530 #endif
global_data_all_processes All
Definition: main.cc:40
MyIDType get(void) const
Definition: idstorage.h:87
void intpos_to_pos(MyIntPosType *intpos, T *posdiff)
void mpi_printf(const char *fmt,...)
Definition: setcomm.h:55
int ThisTask
Definition: setcomm.h:33
MPI_Comm Communicator
Definition: setcomm.h:31
int get_timestep_bin(integertime ti_step)
Definition: timestep.cc:419
void copy_particle(particle_data *Ptarget, particle_data *Psource)
Definition: simparticles.h:65
unsigned short int MarkerValue
Definition: simparticles.h:61
int get_active_index(int idx)
Definition: simparticles.h:382
MyFloat get_DtHsml(int i)
Definition: simparticles.h:231
subfind_data * PS
Definition: simparticles.h:63
void timebin_cleanup_list_of_active_particles(void)
Definition: timestep.cc:636
void fill_active_gravity_list_with_all_particles(void)
Definition: simparticles.h:260
integertime get_timestep_grav(int p)
Definition: timestep.cc:177
void assign_hydro_timesteps(void)
Definition: timestep.cc:56
void print_particle_info_from_ID(MyIDType ID)
Print information relative to a particle / cell to standard output given its ID. *.
Definition: simparticles.h:374
int drift_particle(particle_data *P, sph_particle_data *SphP, integertime time1, bool ignore_light_cone=false)
This function drifts a particle i to time1.
Definition: predict.cc:313
long long TotNumPart
Definition: simparticles.h:46
void reallocate_memory_maxpartsph(int maxpartsphNew)
Definition: simparticles.h:304
integertime get_timestep_hydro(int p)
Definition: timestep.cc:274
void reconstruct_timebins(void)
Definition: predict.cc:43
long long TotNumGas
Definition: simparticles.h:47
double get_utherm_from_entropy(int i)
Definition: simparticles.h:238
MyFloat get_OldAcc(int i)
Definition: simparticles.h:235
void reallocate_memory_maxpart(int maxpartNew)
Definition: simparticles.h:292
void allocate_memory(void)
Definition: simparticles.h:273
simparticles(MPI_Comm comm)
Definition: simparticles.h:38
MyFloat get_Csnd(int i)
Definition: simparticles.h:233
int test_if_grav_timestep_is_too_large(int p, int bin)
Definition: timestep.cc:158
void free_memory(void)
Definition: simparticles.h:283
void set_entropy_from_utherm(double utherm, int i)
Definition: simparticles.h:249
TimeBinData TimeBinsGravity
Definition: simparticles.h:205
void print_particle_info(int i)
Print information relative to a particle / cell to standard output.
Definition: simparticles.h:343
void dump_particles(void)
Definition: simparticles.h:321
particle_data * P
Definition: simparticles.h:54
sph_particle_data * SphP
Definition: simparticles.h:59
void timebins_get_bin_and_do_validity_checks(integertime ti_step, int *bin_new, int bin_old)
Definition: timestep.cc:493
TimeBinData TimeBinsHydro
Definition: simparticles.h:204
integertime find_next_sync_point(void)
This function finds the next synchronization point of the system. (i.e. the earliest point of time an...
Definition: predict.cc:120
int TimeBinSynchronized[TIMEBINS]
Definition: simparticles.h:203
particle_data pdata
Definition: simparticles.h:49
int getTimeBinSynchronized(int bin)
Definition: simparticles.h:214
void make_list_of_active_particles(void)
Definition: predict.cc:410
static bool compare_IDs(const MyIDType &a, const MyIDType &b)
Definition: simparticles.h:72
void drift_all_particles(void)
Definition: predict.cc:274
void mark_active_timebins(void)
Definition: predict.cc:155
#define FLAG_OUTSIDE
Definition: constants.h:29
#define FLAG_BOUNDARYOVERLAP
Definition: constants.h:31
#define FLAG_INSIDE
Definition: constants.h:30
#define HIGH_MESH
Definition: constants.h:34
#define TIMEBINS
Definition: constants.h:332
#define NSOFTCLASSES_HYDRO
Definition: constants.h:321
#define NSOFTCLASSES
Definition: constants.h:312
#define GAMMA_MINUS1
Definition: constants.h:110
#define GAMMA
Definition: constants.h:99
int integertime
Definition: constants.h:331
#define MAX_DOUBLE_NUMBER
Definition: constants.h:81
#define MAX_FLOAT_NUMBER
Definition: constants.h:79
#define M_PI
Definition: constants.h:56
float MyFloat
Definition: dtypes.h:86
uint32_t MyIntPosType
Definition: dtypes.h:35
int32_t MySignedIntPosType
Definition: dtypes.h:36
double MyReal
Definition: dtypes.h:82
unsigned int MyIDType
Definition: dtypes.h:68
#define Terminate(...)
Definition: macros.h:19
void sumup_large_ints(int n, int *src, long long *res, MPI_Comm comm)
memory Mem
Definition: main.cc:44
half fabs(half arg)
Definition: half.hpp:2627
expr log(half arg)
Definition: half.hpp:2745
expr pow(half base, half exp)
Definition: half.hpp:2803
void timebins_reallocate(void)
Definition: timestep.cc:483
int NActiveParticles
Definition: timestep.h:18
int * ActiveParticleList
Definition: timestep.h:20
void timebins_allocate(void)
Definition: timestep.cc:457
void timebins_free(void)
Definition: timestep.cc:472
int SofteningClassOfPartType[NTYPES]
Definition: allvars.h:250
double ForceSoftening[NSOFTCLASSES+NSOFTCLASSES_HYDRO+2]
Definition: allvars.h:258
double SofteningTable[NSOFTCLASSES+NSOFTCLASSES_HYDRO]
Definition: allvars.h:256
void setSofteningClass(unsigned char softclass)
MyDouble getMass(void)
unsigned char getSofteningClass(void)
unsigned char getTimeBinHydro(void)
MyFloat Vel[3]
Definition: particle_data.h:54
MyIDStorage ID
Definition: particle_data.h:70
signed char TimeBinGrav
Definition: particle_data.h:71
unsigned char getType(void)
void myflush(FILE *fstream)
Definition: system.cc:125