GADGET-4
predict.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 "../data/intposconvert.h"
23 #include "../data/simparticles.h"
24 #include "../lightcone/lightcone.h"
25 #include "../logs/logs.h"
26 #include "../logs/timer.h"
27 #include "../main/main.h"
28 #include "../main/simulation.h"
29 #include "../mpi_utils/mpi_utils.h"
30 #include "../system/system.h"
31 #include "../time_integration/driftfac.h"
32 #include "../time_integration/timestep.h"
33 
34 /*
35  * It counts the number of particles in each timebin and updates the
36  * linked lists containing the particles of each time bin. Afterwards the
37  * linked list of active particles is updated by make_list_of_active_particles().
38  *
39  * The linked lists for each timebin are stored in #FirstInTimeBin[], #LastInTimeBin[],
40  * #PrevInTimeBin[] and #NextInTimeBin[]. The counters of particles per timebin are
41  * #TimeBinCount and #TimeBinCountSph.
42  */
44 {
45  TIMER_START(CPU_TIMELINE);
46 
47  for(int bin = 0; bin < TIMEBINS; bin++)
48  {
49  TimeBinsHydro.TimeBinCount[bin] = 0;
50  TimeBinsHydro.FirstInTimeBin[bin] = -1;
51  TimeBinsHydro.LastInTimeBin[bin] = -1;
52 
56 
57 #ifdef STARFORMATION
58  TimeBinSfr[bin] = 0;
59 #endif
60  }
61 
62  for(int i = 0; i < NumPart; i++)
63  {
64  int bin = P[i].TimeBinGrav;
65 
66  if(TimeBinsGravity.TimeBinCount[bin] > 0)
67  {
72  }
73  else
74  {
77  }
78 
80 
81  if(P[i].getType() == 0)
82  {
83  bin = P[i].getTimeBinHydro();
84 
85  if(TimeBinsHydro.TimeBinCount[bin] > 0)
86  {
91  }
92  else
93  {
96  }
97 
99 
100 #ifdef STARFORMATION
101  TimeBinSfr[bin] += SphP[i].Sfr;
102 #endif
103  }
104  }
105 
107 
108  TIMER_STOP(CPU_TIMELINE);
109 }
110 
121 {
122  /* find the next kick time */
123  integertime ti_next_kick = TIMEBASE;
124 
125  for(int n = 0; n < TIMEBINS; n++)
126  {
128  {
129  integertime ti_next_for_bin;
130  if(n > 0)
131  {
132  integertime dt_bin = (((integertime)1) << n);
133  ti_next_for_bin = (All.Ti_Current / dt_bin) * dt_bin + dt_bin; /* next kick time for this timebin */
134  }
135  else
136  {
137  ti_next_for_bin = All.Ti_Current;
138  }
139 
140  if(ti_next_for_bin < ti_next_kick)
141  ti_next_kick = ti_next_for_bin;
142  }
143  }
144 
145  integertime ti_next_kick_global;
146 #ifdef ENLARGE_DYNAMIC_RANGE_IN_TIME
147  minimum_large_ints(1, &ti_next_kick, &ti_next_kick_global, Communicator);
148 #else
149  MPI_Allreduce(&ti_next_kick, &ti_next_kick_global, 1, MPI_INT, MPI_MIN, Communicator);
150 #endif
151 
152  return ti_next_kick_global;
153 }
154 
156 {
157  int lowest_active_bin = TIMEBINS, highest_active_bin = 0;
158  int lowest_occupied_bin = TIMEBINS, highest_occupied_bin = 0;
159  int lowest_occupied_gravity_bin = TIMEBINS, highest_occupied_gravity_bin = 0;
160  int highest_synchronized_bin = 0;
161  int nsynchronized_gravity = 0, nsynchronized_hydro = 0;
162 
163  /* mark the bins that will be synchronized/active */
164 
165  for(int n = 0; n < TIMEBINS; n++)
166  {
168  {
169  if(highest_occupied_gravity_bin < n)
170  highest_occupied_gravity_bin = n;
171 
172  if(lowest_occupied_gravity_bin > n)
173  lowest_occupied_gravity_bin = n;
174  }
175 
177 
178  if(active)
179  {
180  if(highest_occupied_bin < n)
181  highest_occupied_bin = n;
182 
183  if(lowest_occupied_bin > n)
184  lowest_occupied_bin = n;
185  }
186 
187  integertime dt_bin = (((integertime)1) << n);
188 
189  if((All.Ti_Current % dt_bin) == 0)
190  {
191  TimeBinSynchronized[n] = 1;
193 
194  nsynchronized_gravity += TimeBinsGravity.TimeBinCount[n];
195  nsynchronized_hydro += TimeBinsHydro.TimeBinCount[n];
196 
197  if(highest_synchronized_bin < n)
198  highest_synchronized_bin = n;
199 
200  if(active)
201  {
202  if(highest_active_bin < n)
203  highest_active_bin = n;
204 
205  if(lowest_active_bin > n)
206  lowest_active_bin = n;
207  }
208  }
209  else
210  TimeBinSynchronized[n] = 0;
211  }
212 
213  int lowest_in[3], lowest_out[3];
214  lowest_in[0] = lowest_occupied_bin;
215  lowest_in[1] = lowest_occupied_gravity_bin;
216  lowest_in[2] = lowest_active_bin;
217  MPI_Allreduce(lowest_in, lowest_out, 3, MPI_INT, MPI_MIN, Communicator);
218  All.LowestOccupiedTimeBin = lowest_out[0];
219  All.LowestOccupiedGravTimeBin = lowest_out[1];
220  All.LowestActiveTimeBin = lowest_out[2];
221 
222  int highest_in[4], highest_out[4];
223  highest_in[0] = highest_occupied_bin;
224  highest_in[1] = highest_occupied_gravity_bin;
225  highest_in[2] = highest_active_bin;
226  highest_in[3] = highest_synchronized_bin;
227  MPI_Allreduce(highest_in, highest_out, 4, MPI_INT, MPI_MAX, Communicator);
228  All.HighestOccupiedTimeBin = highest_out[0];
229  All.HighestOccupiedGravTimeBin = highest_out[1];
230  All.HighestActiveTimeBin = highest_out[2];
231  All.HighestSynchronizedTimeBin = highest_out[3];
232 
233  /* note: the lowest synchronized bin is always 1 */
234 
235  int input_ints[2 + 2 * TIMEBINS];
236  long long output_longs[2 + 2 * TIMEBINS];
237 
238  input_ints[0] = nsynchronized_hydro;
239  input_ints[1] = nsynchronized_gravity;
240  memcpy(input_ints + 2, TimeBinsGravity.TimeBinCount, TIMEBINS * sizeof(int));
241  memcpy(input_ints + 2 + TIMEBINS, TimeBinsHydro.TimeBinCount, TIMEBINS * sizeof(int));
242 
243  sumup_large_ints(2 + 2 * TIMEBINS, input_ints, output_longs, Communicator);
244 
245  All.GlobalNSynchronizedHydro = output_longs[0];
246  All.GlobalNSynchronizedGravity = output_longs[1];
247  long long *tot_count_grav = output_longs + 2;
248  long long *tot_count_sph = output_longs + 2 + TIMEBINS;
249 
250  long long tot_grav = 0, tot_sph = 0;
251 
252  for(int n = 0; n < TIMEBINS; n++)
253  {
254  tot_grav += tot_count_grav[n];
255  tot_sph += tot_count_sph[n];
256 
257  if(n > 0)
258  {
259  tot_count_grav[n] += tot_count_grav[n - 1];
260  tot_count_sph[n] += tot_count_sph[n - 1];
261  }
262  }
263 
265 
266  for(int n = All.HighestOccupiedTimeBin; n >= All.LowestOccupiedTimeBin; n--)
267  {
268  if(tot_count_grav[n] > All.ActivePartFracForNewDomainDecomp * tot_grav ||
269  tot_count_sph[n] > All.ActivePartFracForNewDomainDecomp * tot_sph)
271  }
272 }
273 
275 {
276  TIMER_START(CPU_DRIFTS);
277 
278  for(int i = 0; i < NumPart; i++)
279  {
280 #ifdef LIGHTCONE_MASSMAPS
281  int flag = drift_particle(&P[i], &SphP[i], All.Ti_Current);
282 
283  if(flag)
284  {
285  MPI_Allreduce(MPI_IN_PLACE, &flag, 1, MPI_INT, MPI_MAX, Communicator);
286  LightCone->lightcone_massmap_flush(0);
287  }
288 #else
289  drift_particle(&P[i], &SphP[i], All.Ti_Current);
290 #endif
291  }
292 
293 #ifdef LIGHTCONE_MASSMAPS
294  int flag = 0;
295  do
296  {
297  flag = 0;
298  MPI_Allreduce(MPI_IN_PLACE, &flag, 1, MPI_INT, MPI_MAX, Communicator);
299  if(flag)
300  LightCone->lightcone_massmap_flush(0);
301  }
302  while(flag);
303 #endif
304 
305  TIMER_STOP(CPU_DRIFTS);
306 }
307 
313 int simparticles::drift_particle(particle_data *P, sph_particle_data *SphP, integertime time1, bool ignore_light_cone)
314 {
315  int buffer_full_flag = 0;
316 #ifndef LEAN
317  while(P->access.test_and_set(std::memory_order_acquire))
318  {
319  // acquire spin lock
320  }
321 #endif
322 
323  integertime time0 = P->Ti_Current.load(std::memory_order_acquire);
324 
325  if(time1 == time0)
326  {
327 #ifndef LEAN
328  P->access.clear(std::memory_order_release);
329 #endif
330  return buffer_full_flag;
331  }
332 
333  if(time1 < time0)
334  Terminate("no prediction into past allowed: time0=%lld time1=%lld\n", (long long)time0, (long long)time1);
335 
336  double dt_drift;
337 
339  dt_drift = Driftfac.get_drift_factor(time0, time1);
340  else
341  dt_drift = (time1 - time0) * All.Timebase_interval;
342 
343 #ifdef LIGHTCONE
344  if(ignore_light_cone == false)
345  buffer_full_flag = LightCone->lightcone_test_for_particle_addition(P, time0, time1, dt_drift);
346 #endif
347 
348  double posdiff[3];
349  for(int j = 0; j < 3; j++)
350  posdiff[j] = P->Vel[j] * dt_drift;
351 
352  MyIntPosType delta[3];
353  pos_to_signedintpos(posdiff, (MySignedIntPosType *)delta);
354 
355  for(int j = 0; j < 3; j++)
356  P->IntPos[j] += delta[j];
357 
358  constrain_intpos(P->IntPos); /* will only do something if we have a stretched box */
359 
360  if(P->getType() == 0)
361  {
362  double dt_hydrokick, dt_entr, dt_gravkick;
363 
365  {
366  dt_entr = (time1 - time0) * All.Timebase_interval;
367  dt_hydrokick = Driftfac.get_hydrokick_factor(time0, time1);
368  dt_gravkick = Driftfac.get_gravkick_factor(time0, time1);
369  }
370  else
371  dt_gravkick = dt_entr = dt_hydrokick = dt_drift;
372 
373  for(int j = 0; j < 3; j++)
374  {
375  SphP->VelPred[j] += SphP->HydroAccel[j] * dt_hydrokick;
376 #ifdef HIERARCHICAL_GRAVITY
377  SphP->VelPred[j] += SphP->FullGravAccel[j] * dt_gravkick;
378 #else
379  SphP->VelPred[j] += P->GravAccel[j] * dt_gravkick;
380 #endif
381 #if defined(PMGRID) && !defined(TREEPM_NOTIMESPLIT)
382  SphP->VelPred[j] += P->GravPM[j] * dt_gravkick;
383 #endif
384  }
385 
386  SphP->EntropyPred += SphP->DtEntropy * dt_entr;
387 
388  SphP->Density += SphP->DtDensity * dt_drift;
389 
390  SphP->Hsml += SphP->DtHsml * dt_drift;
391 
392 #ifdef PRESSURE_ENTROPY_SPH
393  SphP->PressureSphDensity += SphP->DtPressureSphDensity * dt_drift;
394 
395  SphP->EntropyToInvGammaPred = pow(SphP->EntropyPred, 1.0 / GAMMA);
396 #endif
397 
399  }
400 
401  P->Ti_Current = time1;
402 
403 #ifndef LEAN
404  P->access.clear(std::memory_order_release);
405 #endif
406 
407  return buffer_full_flag;
408 }
409 
411 {
412  TIMER_START(CPU_DRIFTS);
413 
415 
416  for(int n = 0; n < TIMEBINS; n++)
417  {
418  if(TimeBinSynchronized[n])
419  {
420  for(int i = TimeBinsHydro.FirstInTimeBin[n]; i >= 0; i = TimeBinsHydro.NextInTimeBin[i])
421  {
422  if(P[i].getType() == 0)
423  {
424  if(P[i].getTimeBinHydro() != n)
425  Terminate("P[i].TimeBinHydro=%d != timebin=%d", P[i].getTimeBinHydro(), n);
426 
427  if(P[i].Ti_Current.load(std::memory_order_acquire) != All.Ti_Current)
428  drift_particle(&P[i], &SphP[i], All.Ti_Current);
429 
431  }
432  }
433  }
434  }
435 
437 
438  for(int n = 0; n < TIMEBINS; n++)
439  {
440  if(TimeBinSynchronized[n])
441  {
442  for(int i = TimeBinsGravity.FirstInTimeBin[n]; i >= 0; i = TimeBinsGravity.NextInTimeBin[i])
443  {
444  if(P[i].Ti_Current.load(std::memory_order_acquire) != All.Ti_Current)
445  drift_particle(&P[i], &SphP[i], All.Ti_Current);
446 
448  }
449  }
450  }
451 
453  long long out[2];
454 
455  sumup_large_ints(2, in, out, Communicator);
456 
459 
460  TIMER_STOP(CPU_DRIFTS);
461 }
global_data_all_processes All
Definition: main.cc:40
double get_drift_factor(integertime time0, integertime time1)
Definition: driftfac.cc:68
double get_gravkick_factor(integertime time0, integertime time1)
Definition: driftfac.cc:109
double get_hydrokick_factor(integertime time0, integertime time1)
Definition: driftfac.cc:150
MySignedIntPosType pos_to_signedintpos(T posdiff)
void constrain_intpos(MyIntPosType *pos)
Definition: intposconvert.h:68
MPI_Comm Communicator
Definition: setcomm.h:31
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
void reconstruct_timebins(void)
Definition: predict.cc:43
TimeBinData TimeBinsGravity
Definition: simparticles.h:205
particle_data * P
Definition: simparticles.h:54
sph_particle_data * SphP
Definition: simparticles.h:59
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
void make_list_of_active_particles(void)
Definition: predict.cc:410
void drift_all_particles(void)
Definition: predict.cc:274
void mark_active_timebins(void)
Definition: predict.cc:155
#define TIMEBINS
Definition: constants.h:332
#define GAMMA
Definition: constants.h:99
int integertime
Definition: constants.h:331
#define TIMEBASE
Definition: constants.h:333
uint32_t MyIntPosType
Definition: dtypes.h:35
int32_t MySignedIntPosType
Definition: dtypes.h:36
#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)
void minimum_large_ints(int n, long long *src, long long *res, MPI_Comm comm)
expr pow(half base, half exp)
Definition: half.hpp:2803
int NActiveParticles
Definition: timestep.h:18
int * ActiveParticleList
Definition: timestep.h:20
int TimeBinCount[TIMEBINS]
Definition: timestep.h:21
int FirstInTimeBin[TIMEBINS]
Definition: timestep.h:23
long long GlobalNActiveParticles
Definition: timestep.h:19
int * PrevInTimeBin
Definition: timestep.h:26
int * NextInTimeBin
Definition: timestep.h:25
int LastInTimeBin[TIMEBINS]
Definition: timestep.h:24
long long GlobalNSynchronizedHydro
Definition: allvars.h:90
long long GlobalNSynchronizedGravity
Definition: allvars.h:91
integertime Ti_Current
Definition: allvars.h:188
integertime Ti_begstep[TIMEBINS]
Definition: allvars.h:210
int SmallestTimeBinWithDomainDecomposition
Definition: allvars.h:160
double ActivePartFracForNewDomainDecomp
Definition: allvars.h:161
std::atomic_flag access
Definition: particle_data.h:88
std::atomic< integertime > Ti_Current
Definition: particle_data.h:60
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
MyIntPosType IntPos[3]
Definition: particle_data.h:53
void set_thermodynamic_variables(void)