GADGET-4
kicks.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 <math.h>
15 #include <mpi.h>
16 #include <stdio.h>
17 #include <stdlib.h>
18 #include <string.h>
19 
20 #include "../data/allvars.h"
21 #include "../data/dtypes.h"
22 #include "../domain/domain.h"
23 #include "../gravtree/gravtree.h"
24 #include "../logs/logs.h"
25 #include "../logs/timer.h"
26 #include "../main/main.h"
27 #include "../main/simulation.h"
28 #include "../mpi_utils/mpi_utils.h"
29 #include "../ngbtree/ngbtree.h"
30 #include "../system/system.h"
31 #include "../time_integration/driftfac.h"
32 #include "../time_integration/timestep.h"
33 
41 {
42  TIMER_START(CPU_DRIFTS);
43 
45 
46 #if defined(PMGRID) && !defined(TREEPM_NOTIMESPLIT)
47  if(All.PM_Ti_endstep == All.Ti_Current) /* need to do long-range kick */
48  {
49  /* first, determine the new PM timestep */
50  integertime ti_step = Sp.get_timestep_pm();
51 
52  All.PM_Ti_begstep = All.PM_Ti_endstep;
53  All.PM_Ti_endstep = All.PM_Ti_begstep + ti_step;
54 
55  integertime tstart = All.PM_Ti_begstep;
56  integertime tend = tstart + ti_step / 2;
57 
58  double dt_gravkick;
59 
61  dt_gravkick = Driftfac.get_gravkick_factor(tstart, tend);
62  else
63  dt_gravkick = (tend - tstart) * All.Timebase_interval;
64 
65  for(int i = 0; i < Sp.NumPart; i++)
66  {
67  for(int j = 0; j < 3; j++)
68  Sp.P[i].Vel[j] += Sp.P[i].GravPM[j] * dt_gravkick;
69  }
70  }
71 #endif
72 
75 
76 #ifdef FORCE_EQUAL_TIMESTEPS
78 #endif
79 
80 #ifdef HIERARCHICAL_GRAVITY
81  /* First, move all active particles to the highest allowed timestep for this synchronization time.
82  * They will then cascade down to smaller timesteps as needed.
83  */
84  for(int i = 0; i < Sp.TimeBinsGravity.NActiveParticles; i++)
85  {
86  int target = Sp.TimeBinsGravity.ActiveParticleList[i];
88  int binold = Sp.P[target].TimeBinGrav;
89 
90  Sp.TimeBinsGravity.timebin_move_particle(target, binold, bin);
91  Sp.P[target].TimeBinGrav = bin;
92  }
93 
94  long long Previous_GlobalNActiveGravity = Sp.TimeBinsGravity.GlobalNActiveParticles;
95 
96  double dt_gravsum = 0;
97 
98  int bin_highest_occupied = 0;
99 
100  /* go over all timebins */
101  for(int timebin = All.HighestSynchronizedTimeBin; timebin >= 0; timebin--)
102  {
104 
107 
108  if(Sp.TimeBinsGravity.GlobalNActiveParticles == 0) /* we are done at this point */
109  break;
110 
111  /* calculate gravity for all active particles */
112  if(Sp.TimeBinsGravity.GlobalNActiveParticles != Previous_GlobalNActiveGravity)
113  {
114  TIMER_STOP(CPU_DRIFTS);
115 
117 
118  TIMER_START(CPU_DRIFTS);
119  }
120 
121  /* now check whether the current timestep should be reduced */
122 
123  int nfine = 0;
124  for(int i = 0; i < Sp.TimeBinsGravity.NActiveParticles; i++)
125  {
126  int target = Sp.TimeBinsGravity.ActiveParticleList[i];
127  int binold = Sp.P[target].TimeBinGrav;
128 
129  if(Sp.test_if_grav_timestep_is_too_large(target, binold))
130  nfine++;
131  }
132 
133  long long nfine_tot;
134  sumup_large_ints(1, &nfine, &nfine_tot, Communicator);
135 
136  int push_down_flag = 0;
138  nfine_tot > 0.33 * Sp.TimeBinsGravity.GlobalNActiveParticles)
139  {
140  mpi_printf(
141  "KICKS: We reduce the highest occupied timestep by pushing %lld particles on timebin=%d down in timestep (fraction "
142  "wanting lower bin is: %g)\n",
143  Sp.TimeBinsGravity.GlobalNActiveParticles - nfine_tot, timebin,
144  ((double)nfine_tot) / Sp.TimeBinsGravity.GlobalNActiveParticles);
145 
146  push_down_flag = 1;
147  }
148 
149  for(int i = 0; i < Sp.TimeBinsGravity.NActiveParticles; i++)
150  {
151  int target = Sp.TimeBinsGravity.ActiveParticleList[i];
152  int binold = Sp.P[target].TimeBinGrav;
153 
154  if(push_down_flag || Sp.test_if_grav_timestep_is_too_large(target, binold))
155  {
156  int bin = binold - 1;
157 
158  if(bin == 0)
159  {
160  Sp.print_particle_info(target);
161  Terminate("timestep too small");
162  }
163 
164  Sp.TimeBinsGravity.timebin_move_particle(target, binold, bin);
165  Sp.P[target].TimeBinGrav = bin;
166  }
167  else if(binold > bin_highest_occupied)
168  bin_highest_occupied = binold;
169  }
170 
171  if(All.HighestOccupiedGravTimeBin == 0) /* this will only be the case in the very first step */
172  {
173  MPI_Allreduce(&bin_highest_occupied, &All.HighestOccupiedGravTimeBin, 1, MPI_INT, MPI_MAX, Communicator);
175  mpi_printf("KICKS: Special Start-up: All.HighestOccupiedGravTimeBin=%d\n", All.HighestOccupiedGravTimeBin);
176  }
177 
179  {
180  integertime ti_step = timebin ? (((integertime)1) << timebin) : 0;
181 
182  integertime tstart = All.Ti_begstep[timebin]; /* beginning of step */
183  integertime tend = tstart + ti_step / 2; /* midpoint of step */
184 
185  double dt_gravkick;
186 
188  dt_gravkick = Driftfac.get_gravkick_factor(tstart, tend);
189  else
190  dt_gravkick = (tend - tstart) * All.Timebase_interval;
191 
192  double dt_save = dt_gravkick;
193 
194  if(timebin < All.HighestSynchronizedTimeBin)
195  {
196  ti_step = (timebin + 1) ? (((integertime)1) << (timebin + 1)) : 0;
197 
198  tstart = All.Ti_begstep[timebin + 1]; /* beginning of step */
199  tend = tstart + ti_step / 2; /* midpoint of step */
200 
202  dt_gravkick -= Driftfac.get_gravkick_factor(tstart, tend);
203  else
204  dt_gravkick -= (tend - tstart) * All.Timebase_interval;
205  }
206 
207  dt_gravsum += dt_gravkick;
208 
209  mpi_printf("KICKS: 1st gravity for hierarchical timebin=%d: %lld particles dt_gravkick=%g %g %g\n", timebin,
210  Sp.TimeBinsGravity.GlobalNActiveParticles, dt_gravkick, dt_gravsum, dt_save);
211 
212  for(int i = 0; i < Sp.TimeBinsGravity.NActiveParticles; i++)
213  {
214  int target = Sp.TimeBinsGravity.ActiveParticleList[i];
215 
216  for(int j = 0; j < 3; j++)
217  Sp.P[target].Vel[j] += Sp.P[target].GravAccel[j] * dt_gravkick;
218  }
219  }
220 
221  Previous_GlobalNActiveGravity = Sp.TimeBinsGravity.GlobalNActiveParticles;
222  }
223 
224 #else
225 
226  mpi_printf("KICKS: 1st gravity for highest active timebin=%d: particles %lld\n", All.HighestActiveTimeBin,
228 
229  for(int i = 0; i < Sp.TimeBinsGravity.NActiveParticles; i++)
230  {
231  int target = Sp.TimeBinsGravity.ActiveParticleList[i];
232 
233 #ifndef FORCE_EQUAL_TIMESTEPS
234  integertime ti_step = Sp.get_timestep_grav(target);
235  int timebin;
236 
237  Sp.timebins_get_bin_and_do_validity_checks(ti_step, &timebin, Sp.P[target].TimeBinGrav);
238 
239  ti_step = timebin ? (((integertime)1) << timebin) : 0;
240 
241  Sp.TimeBinsGravity.timebin_move_particle(target, Sp.P[target].TimeBinGrav, timebin);
242  Sp.P[target].TimeBinGrav = timebin;
243 #else
244  int timebin = Sp.P[target].TimeBinGrav;
245  integertime ti_step = timebin ? (((integertime)1) << timebin) : 0;
246 #endif
247 
248  integertime tstart = All.Ti_begstep[timebin]; /* beginning of step */
249  integertime tend = tstart + ti_step / 2; /* midpoint of step */
250 
251  double dt_gravkick;
252 
254  dt_gravkick = Driftfac.get_gravkick_factor(tstart, tend);
255  else
256  dt_gravkick = (tend - tstart) * All.Timebase_interval;
257 
258  for(int j = 0; j < 3; j++)
259  Sp.P[target].Vel[j] += Sp.P[target].GravAccel[j] * dt_gravkick;
260  }
261 
262 #endif
263 
264  TIMER_STOP(CPU_DRIFTS);
265 }
266 
275 {
276  TIMER_START(CPU_DRIFTS);
277 
278  char fullmark[8];
279 
281  sprintf(fullmark, "(*)");
282  else
283  fullmark[0] = 0;
284 
285  if(ThisTask == 0)
286  {
287  fprintf(Logs.FdTimings, "\nStep%s: %d, t: %g, dt: %g, highest active timebin: %d (lowest active: %d, highest occupied: %d)\n",
290 
291  fprintf(Logs.FdDensity, "\nStep%s: %d, t: %g, dt: %g, highest active timebin: %d (lowest active: %d, highest occupied: %d)\n",
294 
295  fprintf(Logs.FdHydro, "\nStep%s: %d, t: %g, dt: %g, highest active timebin: %d (lowest active: %d, highest occupied: %d)\n",
298  }
299 
300  double dt_gravkick;
301 
302 #ifdef HIERARCHICAL_GRAVITY
303  /* go over all timebins, in inverse sequence so that we end up getting the cumulative force at the end */
304  for(int timebin = 0; timebin <= All.HighestActiveTimeBin; timebin++)
305  {
306  if(Sp.TimeBinSynchronized[timebin])
307  {
308  /* need to make all timebins below the current one active */
311 
313  {
314  /* calculate gravity for all active particles */
315 
316  TIMER_STOP(CPU_DRIFTS);
317 
319 
320  TIMER_START(CPU_DRIFTS);
321 
322  mpi_printf("KICKS: 2nd gravity for hierarchical timebin=%d: particles %lld\n", timebin,
324 
325  integertime ti_step = timebin ? (((integertime)1) << timebin) : 0;
326 
327  integertime tend = All.Ti_begstep[timebin]; /* end of step (Note: All.Ti_begstep[] has already been advanced for the next
328  step at this point) */
329  integertime tstart = tend - ti_step / 2; /* midpoint of step */
330 
332  dt_gravkick = Driftfac.get_gravkick_factor(tstart, tend);
333  else
334  dt_gravkick = (tend - tstart) * All.Timebase_interval;
335 
336  if(timebin < All.HighestActiveTimeBin)
337  {
338  ti_step = (timebin + 1) ? (((integertime)1) << (timebin + 1)) : 0;
339 
340  tend = All.Ti_begstep[timebin + 1]; /* end of step (Note: All.Ti_begstep[] has already been advanced for the next
341  step at this point) */
342  tstart = tend - ti_step / 2; /* midpoint of step */
343 
345  dt_gravkick -= Driftfac.get_gravkick_factor(tstart, tend);
346  else
347  dt_gravkick -= (tend - tstart) * All.Timebase_interval;
348  }
349 
350  for(int i = 0; i < Sp.TimeBinsGravity.NActiveParticles; i++)
351  {
352  int target = Sp.TimeBinsGravity.ActiveParticleList[i];
353 
354  for(int j = 0; j < 3; j++)
355  Sp.P[target].Vel[j] += Sp.P[target].GravAccel[j] * dt_gravkick;
356 
357  if(Sp.P[target].getType() == 0 && All.HighestOccupiedGravTimeBin == timebin)
358  {
359  for(int j = 0; j < 3; j++)
360  {
361  Sp.SphP[target].VelPred[j] = Sp.P[target].Vel[j];
362  Sp.SphP[target].FullGravAccel[j] = Sp.P[target].GravAccel[j];
363  }
364  }
365  }
366  }
367  }
368  }
369 #else
370 
373 
375  {
376  TIMER_STOP(CPU_DRIFTS);
377 
378  /* calculate gravity for all active particles */
380 
381  TIMER_START(CPU_DRIFTS);
382 
383  mpi_printf("KICKS: 2nd gravity for highest active timebin=%d: particles %lld\n", All.HighestActiveTimeBin,
385 
386  for(int i = 0; i < Sp.TimeBinsGravity.NActiveParticles; i++)
387  {
388  int target = Sp.TimeBinsGravity.ActiveParticleList[i];
389  int timebin = Sp.P[target].TimeBinGrav;
390 
391  integertime ti_step = (timebin) ? (((integertime)1) << (timebin)) : 0;
392  integertime tend = All.Ti_Current;
393  integertime tstart = tend - ti_step / 2; /* midpoint of step */
394 
396  dt_gravkick = Driftfac.get_gravkick_factor(tstart, tend);
397  else
398  dt_gravkick = (tend - tstart) * All.Timebase_interval;
399 
400  for(int j = 0; j < 3; j++)
401  Sp.P[target].Vel[j] += Sp.P[target].GravAccel[j] * dt_gravkick;
402 
403  if(Sp.P[target].getType() == 0)
404  {
405  for(int j = 0; j < 3; j++)
406  Sp.SphP[target].VelPred[j] = Sp.P[target].Vel[j];
407  }
408  }
409  }
410 
411 #endif
412 
413 #if defined(PMGRID) && defined(PERIODIC) && !defined(TREEPM_NOTIMESPLIT)
414  if(All.PM_Ti_endstep == All.Ti_Current) /* need to do long-range kick */
415  {
416  TIMER_STOP(CPU_DRIFTS);
417 
419 
420  TIMER_START(CPU_DRIFTS);
421 
422  integertime ti_step = All.PM_Ti_endstep - All.PM_Ti_begstep;
423  integertime tstart = All.PM_Ti_begstep + ti_step / 2;
424  integertime tend = tstart + ti_step / 2;
425 
427  dt_gravkick = Driftfac.get_gravkick_factor(tstart, tend);
428  else
429  dt_gravkick = (tend - tstart) * All.Timebase_interval;
430 
431  for(int i = 0; i < Sp.NumPart; i++)
432  for(int j = 0; j < 3; j++)
433  Sp.P[i].Vel[j] += Sp.P[i].GravPM[j] * dt_gravkick;
434 
435  for(int i = 0; i < Sp.NumGas; i++)
436  if(Sp.P[i].getType() == 0)
437  for(int j = 0; j < 3; j++)
438  Sp.SphP[i].VelPred[j] = Sp.P[i].Vel[j];
439 
441  }
442 #else
444 #endif
445 
446  TIMER_STOP(CPU_DRIFTS);
447 }
448 
450 {
451  if(All.NumCurrentTiStep == 0) /* special domain decomposition now that we know the timebins of both gravity and hydro */
452  {
454 
455  NgbTree.treefree();
456 
459 
461  NgbTree.treebuild(Sp.NumGas, NULL);
462  }
463 
464  /* now we can calculate the hydro forces */
465  hydro_force(FIRST_HALF_STEP); /* computes hydrodynamical accelerations and rate of chnange of entropy,
466  and applies this where appropriate directly (half step kicks) */
467 }
468 
470 {
471  /* now we can calculate the hydro forces */
472  hydro_force(SECOND_HALF_STEP); /* computes hydrodynamical accelerations and change rate of entropy,
473  and applies this where appropriate directly (half step kicks) */
474 }
475 
480 void sim::hydro_force(int step_indicator)
481 {
483  return;
484 
485  /* Create list of targets. */
486  int *targetlist = (int *)Mem.mymalloc("targetlist", Sp.NumGas * sizeof(int));
487 
488  struct old_hydro_accel
489  {
490  MyFloat HydroAccel[3];
491  };
492  old_hydro_accel *Old = NULL;
493 
494  if(step_indicator == SECOND_HALF_STEP)
495  Old = (old_hydro_accel *)Mem.mymalloc("Old", Sp.TimeBinsHydro.NActiveParticles * sizeof(old_hydro_accel));
496 
497  int Nhydroforces = 0;
498 
499  for(int i = 0; i < Sp.TimeBinsHydro.NActiveParticles; i++)
500  {
501  int target = Sp.TimeBinsHydro.ActiveParticleList[i];
502 
503  if(target < 0 || target >= Sp.NumGas)
504  Terminate("target=%d i=%d\n", target, i);
505 
506  if(step_indicator == SECOND_HALF_STEP)
507  {
508  Old[i].HydroAccel[0] = Sp.SphP[target].HydroAccel[0];
509  Old[i].HydroAccel[1] = Sp.SphP[target].HydroAccel[1];
510  Old[i].HydroAccel[2] = Sp.SphP[target].HydroAccel[2];
511  }
512 
513  targetlist[Nhydroforces++] = target;
514  }
515 
516 #ifdef REUSE_HYDRO_ACCELERATIONS_FROM_PREVIOUS_STEP
517  /* in this case, the forces for the first hydro step are simply kept */
518  if(step_indicator != FIRST_HALF_STEP)
519 #endif
520  {
521  NgbTree.hydro_forces_determine(Nhydroforces, targetlist);
522  }
523 
524  /* let's now do the hydrodynamical kicks */
525  for(int i = 0; i < Nhydroforces; i++)
526  {
527  int target = Sp.TimeBinsHydro.ActiveParticleList[i];
528 
529  int timebin = Sp.P[target].getTimeBinHydro();
530 
531  integertime tstart, tend, ti_step = timebin ? (((integertime)1) << timebin) : 0;
532 
533  if(step_indicator == SECOND_HALF_STEP)
534  {
535  tend = All.Ti_Current;
536  tstart = tend - ti_step / 2; /* midpoint of step */
537  }
538  else
539  {
540  tstart = All.Ti_Current;
541  tend = tstart + ti_step / 2; /* midpoint of step */
542  }
543 
544  double dt_hydrokick, dt_entr = (tend - tstart) * All.Timebase_interval;
545 
547  dt_hydrokick = Driftfac.get_hydrokick_factor(tstart, tend);
548  else
549  dt_hydrokick = dt_entr;
550 
551  Sp.SphP[target].Entropy += Sp.SphP[target].DtEntropy * dt_entr;
552 
553  Sp.P[target].Vel[0] += Sp.SphP[target].HydroAccel[0] * dt_hydrokick;
554  Sp.P[target].Vel[1] += Sp.SphP[target].HydroAccel[1] * dt_hydrokick;
555  Sp.P[target].Vel[2] += Sp.SphP[target].HydroAccel[2] * dt_hydrokick;
556 
557  if(step_indicator == SECOND_HALF_STEP)
558  {
559  Sp.SphP[target].EntropyPred = Sp.SphP[target].Entropy;
560 #ifdef PRESSURE_ENTROPY_SPH
561  Sp.SphP[target].EntropyToInvGammaPred = pow(Sp.SphP[target].EntropyPred, 1.0 / GAMMA);
562 #endif
564 
565  Sp.SphP[target].VelPred[0] += (Sp.SphP[target].HydroAccel[0] - Old[i].HydroAccel[0]) * dt_hydrokick;
566  Sp.SphP[target].VelPred[1] += (Sp.SphP[target].HydroAccel[1] - Old[i].HydroAccel[1]) * dt_hydrokick;
567  Sp.SphP[target].VelPred[2] += (Sp.SphP[target].HydroAccel[2] - Old[i].HydroAccel[2]) * dt_hydrokick;
568 
569  /* note: if there is no gravity, we should instead set VelPred = Vel (if this is not done anymore in the gravity
570  * routine)
571  */
572  }
573  }
574 
575  if(step_indicator == SECOND_HALF_STEP)
576  Mem.myfree(Old);
577 
578  Mem.myfree(targetlist);
579 }
global_data_all_processes All
Definition: main.cc:40
void domain_free(void)
Definition: domain.cc:179
void domain_decomposition(domain_options mode)
Definition: domain.cc:51
double get_gravkick_factor(integertime time0, integertime time1)
Definition: driftfac.cc:109
double get_hydrokick_factor(integertime time0, integertime time1)
Definition: driftfac.cc:150
FILE * FdHydro
Definition: logs.h:49
FILE * FdTimings
Definition: logs.h:47
FILE * FdDensity
Definition: logs.h:48
void mpi_printf(const char *fmt,...)
Definition: setcomm.h:55
int ThisTask
Definition: setcomm.h:33
MPI_Comm Communicator
Definition: setcomm.h:31
void compute_grav_accelerations(int timebin)
This routine computes the gravitational accelerations for all active particles.
Definition: gravity.cc:41
domain< simparticles > Domain
Definition: simulation.h:58
void do_hydro_step_second_half(void)
Definition: kicks.cc:469
void hydro_force(int step_indicator)
Definition: kicks.cc:480
void do_hydro_step_first_half(void)
Definition: kicks.cc:449
void gravity_long_range_force(void)
void find_timesteps_and_do_gravity_step_first_half(void)
performs the first half step kick operator for the gravity
Definition: kicks.cc:40
sph NgbTree
Definition: simulation.h:68
void gravity_set_oldacc(int timebin)
Definition: gravity.cc:248
void find_global_timesteps(void)
void do_gravity_step_second_half(void)
performs the second gravity half step kick operator
Definition: kicks.cc:274
simparticles Sp
Definition: simulation.h:56
integertime get_timestep_grav(int p)
Definition: timestep.cc:177
long long TotNumPart
Definition: simparticles.h:46
int test_if_grav_timestep_is_too_large(int p, int bin)
Definition: timestep.cc:158
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
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
int TimeBinSynchronized[TIMEBINS]
Definition: simparticles.h:203
void mark_active_timebins(void)
Definition: predict.cc:155
void hydro_forces_determine(int ntarget, int *targetlist)
Definition: hydra.cc:271
void treeallocate(int max_partindex, partset *Pptr, domain< partset > *Dptr)
Definition: tree.cc:776
void treefree(void)
Definition: tree.cc:1228
int treebuild(int ninsert, int *indexlist)
Definition: tree.cc:52
#define FIRST_HALF_STEP
Definition: constants.h:26
#define SECOND_HALF_STEP
Definition: constants.h:27
#define GAMMA
Definition: constants.h:99
int integertime
Definition: constants.h:331
@ STANDARD
Definition: domain.h:23
float MyFloat
Definition: dtypes.h:86
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
driftfac Driftfac
Definition: main.cc:41
void sumup_large_ints(int n, int *src, long long *res, MPI_Comm comm)
memory Mem
Definition: main.cc:44
expr pow(half base, half exp)
Definition: half.hpp:2803
int NActiveParticles
Definition: timestep.h:18
int * ActiveParticleList
Definition: timestep.h:20
void timebin_make_list_of_active_particles_up_to_timebin(int timebin)
Definition: timestep.cc:663
void timebin_add_particles_of_timebin_to_list_of_active_particles(int timebin)
Definition: timestep.cc:670
long long GlobalNActiveParticles
Definition: timestep.h:19
void timebin_move_particle(int p, int timeBin_old, int timeBin_new)
Definition: timestep.cc:538
integertime Ti_Current
Definition: allvars.h:188
integertime Ti_begstep[TIMEBINS]
Definition: allvars.h:210
void set_cosmo_factors_for_current_time(void)
Definition: allvars.cc:20
unsigned char getTimeBinHydro(void)
MyFloat Vel[3]
Definition: particle_data.h:54
signed char TimeBinGrav
Definition: particle_data.h:71
unsigned char getType(void)
vector< MyFloat > GravAccel
Definition: particle_data.h:55
void set_thermodynamic_variables(void)