GADGET-4
timestep.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 "../cooling_sfr/cooling.h"
21 #include "../data/allvars.h"
22 #include "../data/dtypes.h"
23 #include "../data/intposconvert.h"
24 #include "../data/mymalloc.h"
25 #include "../data/simparticles.h"
26 #include "../logs/logs.h"
27 #include "../logs/timer.h"
28 #include "../main/simulation.h"
29 #include "../system/system.h"
30 #include "../time_integration/driftfac.h"
31 #include "../time_integration/timestep.h"
32 
40 {
41 #ifndef FORCE_EQUAL_TIMESTEPS
42 
44 
46 
47  TIMER_START(CPU_TIMELINE);
48 
50 
51  TIMER_STOP(CPU_TIMELINE);
52 
53 #endif
54 }
55 
57 {
58  /* Now assign new timesteps for the hydro particles that are synchronized */
59  for(int i = 0; i < TimeBinsHydro.NActiveParticles; i++)
60  {
61  int target = TimeBinsHydro.ActiveParticleList[i];
62  if(P[target].getType() != 0)
63  continue;
64 
65  if(TimeBinSynchronized[P[target].getTimeBinHydro()])
66  {
67  integertime ti_step = get_timestep_hydro(target);
68 
69  int bin;
70  timebins_get_bin_and_do_validity_checks(ti_step, &bin, P[target].getTimeBinHydro());
71 
72 #ifdef SELFGRAVITY
73  /* we enforce that the hydro timestep is nested inside the gravity step */
74  if(bin > P[target].TimeBinGrav)
75  bin = P[target].TimeBinGrav;
76 #endif
77 
78  TimeBinsHydro.timebin_move_particle(target, P[target].getTimeBinHydro(), bin);
79 
80  P[target].setTimeBinHydro(bin);
81  }
82  }
83 }
84 
85 #ifdef FORCE_EQUAL_TIMESTEPS
87 {
89 
91 
92  TIMER_START(CPU_TIMELINE);
93 
94  integertime globTimeStep = TIMEBASE;
95 
96 #if defined(PMGRID) && !defined(TREEPM_NOTIMESPLIT)
97  globTimeStep = Sp.get_timestep_pm();
98 #endif
99 
100 #if defined(SELFGRAVITY) || defined(EXTERNALGRAVITY)
101  for(int idx = 0; idx < Sp.TimeBinsGravity.NActiveParticles; idx++)
102  {
104 
105  integertime ti_step = Sp.get_timestep_grav(i);
106  if(ti_step < globTimeStep)
107  globTimeStep = ti_step;
108  }
109 #endif
110 
111  for(int idx = 0; idx < Sp.TimeBinsHydro.NActiveParticles; idx++)
112  {
113  int i = Sp.TimeBinsHydro.ActiveParticleList[idx];
114  if(Sp.P[i].getType() != 0)
115  continue;
116 
117  integertime ti_step = Sp.get_timestep_hydro(i);
118  if(ti_step < globTimeStep)
119  globTimeStep = ti_step;
120  }
121 
122 #ifdef ENLARGE_DYNAMIC_RANGE_IN_TIME
123  minimum_large_ints(1, &globTimeStep, &All.GlobalTimeStep, Communicator);
124 #else
125  MPI_Allreduce(&globTimeStep, &All.GlobalTimeStep, 1, MPI_INT, MPI_MIN, Communicator);
126 #endif
127 
128  for(int idx = 0; idx < Sp.TimeBinsGravity.NActiveParticles; idx++)
129  {
130  int target = Sp.TimeBinsGravity.ActiveParticleList[idx];
131 
132  int bin;
133 
134  Sp.timebins_get_bin_and_do_validity_checks(All.GlobalTimeStep, &bin, Sp.P[target].TimeBinGrav);
135  Sp.TimeBinsGravity.timebin_move_particle(target, Sp.P[target].TimeBinGrav, bin);
136  Sp.P[target].TimeBinGrav = bin;
137  }
138 
139  for(int idx = 0; idx < Sp.TimeBinsHydro.NActiveParticles; idx++)
140  {
141  int target = Sp.TimeBinsHydro.ActiveParticleList[idx];
142  if(Sp.P[target].getType() != 0)
143  continue;
144 
145  int bin;
146 
147  Sp.timebins_get_bin_and_do_validity_checks(All.GlobalTimeStep, &bin, Sp.P[target].getTimeBinHydro());
148  Sp.TimeBinsHydro.timebin_move_particle(target, Sp.P[target].getTimeBinHydro(), bin);
149 #ifndef LEAN
150  Sp.P[target].TimeBinHydro = bin;
151 #endif
152  }
153 
154  TIMER_STOP(CPU_TIMELINE);
155 }
156 #endif
157 
159 {
160  integertime ti_step_bin = bin ? (((integertime)1) << bin) : 0;
161 
162  integertime ti_step = get_timestep_grav(p);
163 
164  if(ti_step < ti_step_bin)
165  return 1;
166  else
167  return 0;
168 }
169 
178 {
179  double ax = All.cf_a2inv * P[p].GravAccel[0];
180  double ay = All.cf_a2inv * P[p].GravAccel[1];
181  double az = All.cf_a2inv * P[p].GravAccel[2];
182 
183 #if defined(PMGRID) && !defined(TREEPM_NOTIMESPLIT)
184  ax += All.cf_a2inv * P[p].GravPM[0];
185  ay += All.cf_a2inv * P[p].GravPM[1];
186  az += All.cf_a2inv * P[p].GravPM[2];
187 #endif
188 
189  double ac = sqrt(ax * ax + ay * ay + az * az); /* this is now the physical acceleration */
190 
191  if(ac == 0)
192  ac = MIN_FLOAT_NUMBER;
193 
194  /* determine the "kinematic" timestep dt_grav, in physical units */
195 #if NSOFTCLASSES > 1
196  double dt_grav = sqrt(2 * All.ErrTolIntAccuracy * All.cf_atime * All.SofteningTable[P[p].getSofteningClass()] / ac);
197 #else
198  double dt_grav = sqrt(2 * All.ErrTolIntAccuracy * All.cf_atime * All.SofteningTable[0] / ac);
199 #endif
200 
201  double dt = dt_grav;
202 
203  /* convert the physical timestep to dloga if needed. Note: If comoving integration has not been selected,
204  All.cf_hubble_a=1.
205  */
206  dt *= All.cf_hubble_a;
207 
208  if(dt >= All.MaxSizeTimestep)
209  dt = All.MaxSizeTimestep;
210 
211 #if defined(SELFGRAVITY) && defined(PMGRID) && !defined(TREEPM_NOTIMESPLIT)
212  if(dt >= All.DtDisplacement)
213  dt = All.DtDisplacement;
214 #endif
215 
216  if(dt < All.MinSizeTimestep)
217  {
218  dt = All.MinSizeTimestep;
219 #ifndef NO_STOP_BELOW_MINTIMESTEP
220  Terminate(
221  "Timestep wants to be below the limit MinSizeTimestep=%g\n"
222  "Part-ID=%lld task=%d type=%d dtgrav=%g ac=%g soft=%g\n",
223  All.MinSizeTimestep, (long long)P[p].ID.get(), ThisTask, P[p].getType(), dt, ac,
225 #endif
226  }
227 
228  integertime ti_step = (integertime)(dt / All.Timebase_interval);
229 
230  if(!(ti_step > 0 && ti_step < TIMEBASE))
231  {
232  double pos[3];
233  intpos_to_pos(P[p].IntPos, pos); /* converts the integer coordinates to floating point */
234 
235  Terminate(
236  "\nError: A timestep of size zero was assigned on the integer timeline!\n"
237  "We better stop.\n"
238  "Task=%d Part-ID=%lld type=%d dt_grav=%g dt=%g tibase=%g ac=%g xyz=(%g|%g|%g) vel=(%g|%g|%g) tree=(%g|%g|%g) mass=%g\n\n",
239  ThisTask, (long long)P[p].ID.get(), P[p].getType(), dt_grav, dt, All.Timebase_interval, ac, pos[0], pos[1], pos[2],
240  P[p].Vel[0], P[p].Vel[1], P[p].Vel[2], P[p].GravAccel[0], P[p].GravAccel[1], P[p].GravAccel[2], P[p].getMass());
241 
242  myflush(stdout);
243  Terminate("integer timestep reached zero");
244  }
245 
246  return ti_step;
247 }
248 
249 #if defined(PMGRID) && !defined(TREEPM_NOTIMESPLIT)
250 integertime simparticles::get_timestep_pm(void)
251 {
252  integertime ti_step = TIMEBASE;
253  while(ti_step > (All.DtDisplacement / All.Timebase_interval))
254  ti_step >>= 1;
255 
256  if(ti_step > (All.PM_Ti_endstep - All.PM_Ti_begstep)) /* PM-timestep wants to increase */
257  {
258  int bin = get_timestep_bin(ti_step);
259  int binold = get_timestep_bin(All.PM_Ti_endstep - All.PM_Ti_begstep);
260 
261  while(TimeBinSynchronized[bin] == 0 && bin > binold) /* make sure the new step is synchronized */
262  bin--;
263 
264  ti_step = bin ? (((integertime)1) << bin) : 0;
265  }
266 
267  if(All.Ti_Current == TIMEBASE) /* we here finish the last timestep. */
268  ti_step = 0;
269 
270  return ti_step;
271 }
272 #endif
273 
275 {
276  if(P[p].getType() != 0)
277  Terminate("P[p].getType() != 0");
278 
279  double ax = All.cf_afac2 * SphP[p].HydroAccel[0];
280  double ay = All.cf_afac2 * SphP[p].HydroAccel[1];
281  double az = All.cf_afac2 * SphP[p].HydroAccel[2];
282 
283  ax += All.cf_a2inv * P[p].GravAccel[0];
284  ay += All.cf_a2inv * P[p].GravAccel[1];
285  az += All.cf_a2inv * P[p].GravAccel[2];
286 
287 #if defined(PMGRID) && !defined(TREEPM_NOTIMESPLIT)
288  ax += All.cf_a2inv * P[p].GravPM[0];
289  ay += All.cf_a2inv * P[p].GravPM[1];
290  az += All.cf_a2inv * P[p].GravPM[2];
291 #endif
292 
293  double ac = sqrt(ax * ax + ay * ay + az * az); /* this is now the physical acceleration */
294 
295  if(ac == 0)
296  ac = MIN_FLOAT_NUMBER;
297 
298  /* determine the "kinematic" timestep dt_grav, in physical units */
299  double dt_kin = sqrt(2 * All.ErrTolIntAccuracy * All.cf_atime * SphP[p].Hsml / ac);
300 
301  /* calculate local Courant timestep and treebased maximum timestep in physical units */
302 
303  double dt_courant = (All.cf_atime / All.cf_afac3) * All.CourantFac * 2.0 * SphP[p].Hsml / (SphP[p].MaxSignalVel + MIN_FLOAT_NUMBER);
304 
305  double dt_treebased = (All.cf_atime / All.cf_afac3) * SphP[p].CurrentMaxTiStep;
306 
307  /* calculate a timestep that restricts the rate at which the smoothing length may change,
308  * in physical units
309  */
310  double dt_hsml = All.cf_atime2 * All.CourantFac * fabs(SphP[p].Hsml / (SphP[p].DtHsml + MIN_FLOAT_NUMBER));
311 
312  /* now take the smallest of these four criteria */
313  double dt = dt_kin;
314  if(dt > dt_courant)
315  dt = dt_courant;
316  if(dt > dt_treebased)
317  dt = dt_treebased;
318  if(dt > dt_hsml)
319  dt = dt_hsml;
320 
321 #ifdef STARFORMATION
322  if(P[p].getType() == 0) /* to protect using a particle that has been turned into a star */
323  {
324  if(SphP[p].Sfr > 0)
325  {
326  double dt_sfr = 0.1 * P[p].getMass() / SphP[p].Sfr;
327  if(dt_sfr < dt)
328  dt = dt_sfr;
329  }
330  }
331 #endif
332 
333  /* convert the physical timestep to dloga in the cosmological case.
334  * Note: If comoving integration has not been selected, All.cf_hubble_a = 1.0.
335  */
336  dt *= All.cf_hubble_a;
337 
338  if(dt >= All.MaxSizeTimestep)
339  dt = All.MaxSizeTimestep;
340 
341 #if defined(PMGRID) && !defined(TREEPM_NOTIMESPLIT)
342  if(dt >= All.DtDisplacement)
343  dt = All.DtDisplacement;
344 #endif
345 
346  if(dt < All.MinSizeTimestep)
347  {
348  if(P[p].getType() == 0)
349  Terminate(
350  "Timestep wants to be below the limit MinSizeTimestep=%g\n"
351  "Part-ID=%lld task=%d dtkin=%g dtcourant=%g ac=%g\n",
352  All.MinSizeTimestep, (long long)P[p].ID.get(), ThisTask, dt_kin * All.cf_hubble_a, dt_courant * All.cf_hubble_a, ac);
353  dt = All.MinSizeTimestep;
354  }
355 
356  integertime ti_step = (integertime)(dt / All.Timebase_interval);
357 
358  if(!(ti_step > 0 && ti_step < TIMEBASE))
359  {
360  double pos[3];
361  intpos_to_pos(P[p].IntPos, pos); /* converts the integer coordinates to floating point */
362 
363  Terminate(
364  "\nError: A timestep of size zero was assigned on the integer timeline!\n"
365  "We better stop.\n"
366  "Task=%d Part-ID=%lld type=%d dt=%g dtc=%g dt_kin=%g dt_treebased=%g dt_hsml=%g tibase=%g ti_step=%d ac=%g xyz=(%g|%g|%g) "
367  "vel=(%g|%g|%g) "
368  "tree=(%g|%g|%g) mass=%g All.cf_hubble_a=%g\n\n",
369  ThisTask, (long long)P[p].ID.get(), P[p].getType(), dt, dt_courant, dt_kin, dt_treebased, dt_hsml, All.Timebase_interval,
370  (int)ti_step, ac, pos[0], pos[1], pos[2], P[p].Vel[0], P[p].Vel[1], P[p].Vel[2], P[p].GravAccel[0], P[p].GravAccel[1],
371  P[p].GravAccel[2], P[p].getMass(), All.cf_hubble_a);
372  }
373 
374  return ti_step;
375 }
376 
377 #if defined(PMGRID) && defined(PERIODIC) && !defined(TREEPM_NOTIMESPLIT)
378 void simparticles::find_long_range_step_constraint(void)
379 {
380  double dtmin = MAX_DOUBLE_NUMBER;
381 
382  for(int p = 0; p < NumPart; p++)
383  {
384  if(P[p].getType() == 0)
385  continue;
386 
387  /* calculate acceleration */
388  double ax = All.cf_a2inv * P[p].GravPM[0];
389  double ay = All.cf_a2inv * P[p].GravPM[1];
390  double az = All.cf_a2inv * P[p].GravPM[2];
391 
392  double ac = sqrt(ax * ax + ay * ay + az * az); /* this is now the physical acceleration */
393 
394  if(ac < 1.0e-30)
395  ac = 1.0e-30;
396 
397 #if NSOFTCLASSES > 1
398  double dt = sqrt(2 * All.ErrTolIntAccuracy * All.cf_atime * All.ForceSoftening[P[p].getSofteningClass()] / 2.8 / ac);
399 #else
400  double dt = sqrt(2 * All.ErrTolIntAccuracy * All.cf_atime * All.ForceSoftening[0] / 2.8 / ac);
401 #endif
402  dt *= All.cf_hubble_a;
403 
404  if(dt < dtmin)
405  dtmin = dt;
406  }
407 
408  dtmin *= 2.0; /* move it one timebin higher to prevent being too conservative */
409 
410  MPI_Allreduce(&dtmin, &All.DtDisplacement, 1, MPI_DOUBLE, MPI_MIN, Communicator);
411 
412  mpi_printf("TIMESTEPS: displacement time constraint: %g (%g)\n", All.DtDisplacement, All.MaxSizeTimestep);
413 
414  if(All.DtDisplacement > All.MaxSizeTimestep)
415  All.DtDisplacement = All.MaxSizeTimestep;
416 }
417 #endif
418 
420 {
421  int bin = -1;
422 
423  if(ti_step == 0)
424  return 0;
425 
426  if(ti_step == 1)
427  Terminate("time-step of integer size 1 not allowed\n");
428 
429  while(ti_step)
430  {
431  bin++;
432  ti_step >>= 1;
433  }
434 
435  return bin;
436 }
437 
438 void TimeBinData::timebins_init(const char *name, int *maxPart)
439 {
440  NActiveParticles = 0;
441  ActiveParticleList = 0;
442 
443  for(int i = 0; i < TIMEBINS; i++)
444  {
445  FirstInTimeBin[i] = -1;
446  LastInTimeBin[i] = -1;
447  }
448 
449  NextInTimeBin = 0;
450  PrevInTimeBin = 0;
451 
452  strncpy(Name, name, 99);
453  Name[99] = 0;
454  MaxPart = maxPart;
455 }
456 
458 {
459  char Identifier[200];
460  Identifier[199] = 0;
461 
462  snprintf(Identifier, 199, "NextActiveParticle%s", Name);
463  ActiveParticleList = (int *)Mem.mymalloc_movable(&ActiveParticleList, Identifier, *(MaxPart) * sizeof(int));
464 
465  snprintf(Identifier, 199, "NextInTimeBin%s", Name);
466  NextInTimeBin = (int *)Mem.mymalloc_movable(&NextInTimeBin, Identifier, *(MaxPart) * sizeof(int));
467 
468  snprintf(Identifier, 199, "PrevInTimeBin%s", Name);
469  PrevInTimeBin = (int *)Mem.mymalloc_movable(&PrevInTimeBin, Identifier, *(MaxPart) * sizeof(int));
470 }
471 
473 {
474  Mem.myfree_movable(PrevInTimeBin);
475  Mem.myfree_movable(NextInTimeBin);
476  Mem.myfree_movable(ActiveParticleList);
477 
478  PrevInTimeBin = NULL;
479  NextInTimeBin = NULL;
480  ActiveParticleList = NULL;
481 }
482 
484 {
485  if(ActiveParticleList != NULL)
486  {
487  ActiveParticleList = (int *)Mem.myrealloc_movable(ActiveParticleList, *(MaxPart) * sizeof(int));
488  NextInTimeBin = (int *)Mem.myrealloc_movable(NextInTimeBin, *(MaxPart) * sizeof(int));
489  PrevInTimeBin = (int *)Mem.myrealloc_movable(PrevInTimeBin, *(MaxPart) * sizeof(int));
490  }
491 }
492 
494 {
495  /* make it a power 2 subdivision */
496  integertime ti_min = TIMEBASE;
497  while(ti_min > ti_step)
498  ti_min >>= 1;
499  ti_step = ti_min;
500 
501  /* get timestep bin */
502  int bin = -1;
503 
504  if(ti_step == 0)
505  bin = 0;
506 
507  if(ti_step == 1)
508  Terminate("time-step of integer size 1 not allowed\n");
509 
510  while(ti_step)
511  {
512  bin++;
513  ti_step >>= 1;
514  }
515 
516  if(bin > bin_old) /* timestep wants to increase */
517  {
518  while(TimeBinSynchronized[bin] == 0 && bin > bin_old) /* make sure the new step is synchronized */
519  bin--;
520 
521  ti_step = bin ? (((integertime)1) << bin) : 0;
522  }
523 
524  if(All.Ti_Current >= TIMEBASE) /* we here finish the last timestep. */
525  {
526  ti_step = 0;
527  bin = 0;
528  }
529 
530  if((TIMEBASE - All.Ti_Current) < ti_step) /* check that we don't run beyond the end */
531  {
532  Terminate("we are beyond the end of the timeline"); /* should not happen */
533  }
534 
535  *bin_new = bin;
536 }
537 
538 void TimeBinData::timebin_move_particle(int p, int timeBin_old, int timeBin_new)
539 {
540  if(timeBin_old == timeBin_new)
541  return;
542 
543  TimeBinCount[timeBin_old]--;
544 
545  int prev = PrevInTimeBin[p];
546  int next = NextInTimeBin[p];
547 
548  if(FirstInTimeBin[timeBin_old] == p)
549  FirstInTimeBin[timeBin_old] = next;
550  if(LastInTimeBin[timeBin_old] == p)
551  LastInTimeBin[timeBin_old] = prev;
552  if(prev >= 0)
553  NextInTimeBin[prev] = next;
554  if(next >= 0)
555  PrevInTimeBin[next] = prev;
556 
557  if(TimeBinCount[timeBin_new] > 0)
558  {
559  PrevInTimeBin[p] = LastInTimeBin[timeBin_new];
560  NextInTimeBin[LastInTimeBin[timeBin_new]] = p;
561  NextInTimeBin[p] = -1;
562  LastInTimeBin[timeBin_new] = p;
563  }
564  else
565  {
566  FirstInTimeBin[timeBin_new] = LastInTimeBin[timeBin_new] = p;
567  PrevInTimeBin[p] = NextInTimeBin[p] = -1;
568  }
569 
570  TimeBinCount[timeBin_new]++;
571 }
572 
574 {
575  int p = ActiveParticleList[idx];
576  ActiveParticleList[idx] = -1;
577 
578  TimeBinCount[bin]--;
579 
580  if(p >= 0)
581  {
582  int prev = PrevInTimeBin[p];
583  int next = NextInTimeBin[p];
584 
585  if(prev >= 0)
586  NextInTimeBin[prev] = next;
587  if(next >= 0)
588  PrevInTimeBin[next] = prev;
589 
590  if(FirstInTimeBin[bin] == p)
591  FirstInTimeBin[bin] = next;
592  if(LastInTimeBin[bin] == p)
593  LastInTimeBin[bin] = prev;
594  }
595 }
596 
597 /* insert a particle into the timebin struct behind another already existing particle */
598 void TimeBinData::timebin_add_particle(int i_new, int i_old, int timeBin, int addToListOfActiveParticles)
599 {
600  TimeBinCount[timeBin]++;
601 
602  if(i_old < 0)
603  {
604  /* if we don't have an existing particle to add if after, let's take the last one in this timebin */
605  i_old = LastInTimeBin[timeBin];
606 
607  if(i_old < 0)
608  {
609  /* the timebin is empty at the moment, so just add the new particle */
610  FirstInTimeBin[timeBin] = i_new;
611  LastInTimeBin[timeBin] = i_new;
612  NextInTimeBin[i_new] = -1;
613  PrevInTimeBin[i_new] = -1;
614  }
615  }
616 
617  if(i_old >= 0)
618  {
619  /* otherwise we added it already */
620  PrevInTimeBin[i_new] = i_old;
621  NextInTimeBin[i_new] = NextInTimeBin[i_old];
622  if(NextInTimeBin[i_old] >= 0)
623  PrevInTimeBin[NextInTimeBin[i_old]] = i_new;
624  NextInTimeBin[i_old] = i_new;
625  if(LastInTimeBin[timeBin] == i_old)
626  LastInTimeBin[timeBin] = i_new;
627  }
628 
629  if(addToListOfActiveParticles)
630  {
633  }
634 }
635 
637 {
638  for(int idx = 0; idx < TimeBinsGravity.NActiveParticles; idx++)
639  {
640  int i = TimeBinsGravity.ActiveParticleList[idx];
641  if(i < 0)
642  continue;
643 
644  if(P[i].ID.get() == 0 && P[i].getMass() == 0)
645  {
646  TimeBinsGravity.timebin_remove_particle(idx, P[i].TimeBinGrav);
647  }
648  }
649 
650  for(int idx = 0; idx < TimeBinsHydro.NActiveParticles; idx++)
651  {
652  int i = TimeBinsHydro.ActiveParticleList[idx];
653  if(i < 0)
654  continue;
655 
656  if(P[i].ID.get() == 0 && P[i].getMass() == 0 && P[i].getType() == 0)
657  {
658  TimeBinsHydro.timebin_remove_particle(idx, P[i].getTimeBinHydro());
659  }
660  }
661 }
662 
664 {
665  NActiveParticles = 0;
666  for(int tbin = timebin; tbin >= 0; tbin--)
668 }
669 
671 {
672  for(int i = FirstInTimeBin[timebin]; i >= 0; i = NextInTimeBin[i])
673  {
676  }
677 }
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
void find_hydro_timesteps(void)
Definition: timestep.cc:39
sph NgbTree
Definition: simulation.h:68
void find_global_timesteps(void)
simparticles Sp
Definition: simulation.h:56
int get_timestep_bin(integertime ti_step)
Definition: timestep.cc:419
void timebin_cleanup_list_of_active_particles(void)
Definition: timestep.cc:636
integertime get_timestep_grav(int p)
Definition: timestep.cc:177
void assign_hydro_timesteps(void)
Definition: timestep.cc:56
integertime get_timestep_hydro(int p)
Definition: timestep.cc:274
int test_if_grav_timestep_is_too_large(int p, int bin)
Definition: timestep.cc:158
TimeBinData TimeBinsGravity
Definition: simparticles.h:205
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 tree_based_timesteps(void)
#define TIMEBINS
Definition: constants.h:332
#define MIN_FLOAT_NUMBER
Definition: constants.h:80
int integertime
Definition: constants.h:331
#define MAX_DOUBLE_NUMBER
Definition: constants.h:81
#define TIMEBASE
Definition: constants.h:333
#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
void minimum_large_ints(int n, long long *src, long long *res, MPI_Comm comm)
memory Mem
Definition: main.cc:44
half fabs(half arg)
Definition: half.hpp:2627
expr sqrt(half arg)
Definition: half.hpp:2777
void timebins_reallocate(void)
Definition: timestep.cc:483
void timebins_init(const char *name, int *MaxPart)
Definition: timestep.cc:438
int NActiveParticles
Definition: timestep.h:18
int * ActiveParticleList
Definition: timestep.h:20
void timebin_remove_particle(int idx, int bin)
Definition: timestep.cc:573
void timebin_make_list_of_active_particles_up_to_timebin(int timebin)
Definition: timestep.cc:663
char Name[100]
Definition: timestep.h:27
int TimeBinCount[TIMEBINS]
Definition: timestep.h:21
void timebin_add_particles_of_timebin_to_list_of_active_particles(int timebin)
Definition: timestep.cc:670
int FirstInTimeBin[TIMEBINS]
Definition: timestep.h:23
void timebins_allocate(void)
Definition: timestep.cc:457
int * PrevInTimeBin
Definition: timestep.h:26
int * MaxPart
Definition: timestep.h:28
void timebin_move_particle(int p, int timeBin_old, int timeBin_new)
Definition: timestep.cc:538
int * NextInTimeBin
Definition: timestep.h:25
int LastInTimeBin[TIMEBINS]
Definition: timestep.h:24
void timebin_add_particle(int i_new, int i_old, int timeBin, int addToListOfActiveParticles)
Definition: timestep.cc:598
void timebins_free(void)
Definition: timestep.cc:472
integertime Ti_Current
Definition: allvars.h:188
double ForceSoftening[NSOFTCLASSES+NSOFTCLASSES_HYDRO+2]
Definition: allvars.h:258
void set_cosmo_factors_for_current_time(void)
Definition: allvars.cc:20
double SofteningTable[NSOFTCLASSES+NSOFTCLASSES_HYDRO]
Definition: allvars.h:256
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 setTimeBinHydro(unsigned char bin)
signed char TimeBinHydro
Definition: particle_data.h:73
vector< MyFloat > GravAccel
Definition: particle_data.h:55
void myflush(FILE *fstream)
Definition: system.cc:125