GADGET-4
run.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 <ctype.h>
15 #include <math.h>
16 #include <mpi.h>
17 #include <stdio.h>
18 #include <stdlib.h>
19 #include <string.h>
20 #include <unistd.h>
21 
22 #include "../cooling_sfr/cooling.h"
23 #include "../data/allvars.h"
24 #include "../data/dtypes.h"
25 #include "../data/mymalloc.h"
26 #include "../domain/domain.h"
27 #include "../gravtree/gravtree.h"
28 #include "../io/io.h"
29 #include "../io/snap_io.h"
30 #include "../lightcone/lightcone_massmap_io.h"
31 #include "../lightcone/lightcone_particle_io.h"
32 #include "../logs/logs.h"
33 #include "../main/main.h"
34 #include "../main/simulation.h"
35 #include "../ngbtree/ngbtree.h"
36 #include "../sort/parallel_sort.h"
37 #include "../system/system.h"
38 
49 void sim::run(void)
50 {
51 #if defined(NGENIC_TEST) && defined(PERIODIC) && defined(PMGRID)
52  snap_io Snap(&Sp, Communicator, All.SnapFormat); /* get an I/O object */
53  Snap.write_snapshot(All.SnapshotFileCount, NORMAL_SNAPSHOT); /* write snapshot file */
54 #if defined(POWERSPEC_ON_OUTPUT)
55  PM.calculate_power_spectra(All.SnapshotFileCount);
56 #endif
57  return;
58 #endif
59 
60  while(1) /* main loop over synchronization points */
61  {
62  /* store old time for logging purposes */
63  All.TimeOld = All.Time;
64 
65  /* determine the next synchronization time at which we have active particles */
66  integertime ti_next_kick_global = Sp.find_next_sync_point();
67 
68 #ifdef OUTPUT_NON_SYNCHRONIZED_ALLOWED
69  while(ti_next_kick_global > All.Ti_nextoutput && All.Ti_nextoutput >= 0)
70  {
74 
77  }
78 #endif
79 
80  All.Ti_Current = ti_next_kick_global;
84 
85 #ifdef LIGHTCONE
86 #ifdef LIGHTCONE_PARTICLES
87  mpi_printf("LIGHTCONE_PARTICLES: Lp.NumPart=%d Checked %d box replicas out of list of length %d\n", Lp.NumPart,
88  LightCone.NumLastCheck, LightCone.NumBoxes);
89 #endif
90 #ifdef LIGHTCONE_MASSMAPS
91  mpi_printf("LIGHTCONE_MASSMAPS: Mp.NumPart=%d \n", Mp.NumPart);
92 #endif
93 #if defined(LIGHTCONE_MASSMAPS) || defined(LIGHTCONE_PARTICLES_GROUPS)
94  Sp.drift_all_particles(); // we do this here to be able to avoid large buffer sizes, if needed multiple binning operations are
95  // done
96 #endif
97 #endif
98 
99  /* mark the timebins that are active on this step */
101 
102  /* create lists with the particles that are synchronized on this step */
104 
105  /* produce some further log messages */
107 
108  /* call functions that update certain 'extra' physics settings to new current time */
109  set_non_standard_physics_for_current_time();
110 
111  /* for sufficiently large steps, carry out a new domain decomposition */
113  {
114  NgbTree.treefree();
116 
118 
119 #ifdef LIGHTCONE_PARTICLES
120  LightCone.lightcone_clear_boxlist(All.Time);
121 #endif
122 
123 #ifdef DEBUG_MD5
124  Logs.log_debug_md5("C");
125 #endif
127 
128 #ifdef DEBUG_MD5
129  Logs.log_debug_md5("D");
130 #endif
131 
133  NgbTree.treebuild(Sp.NumGas, NULL);
134  }
135 
136  /* compute SPH densities and smoothing lengths for active SPH particles, and optionally those
137  * accessed passively. This also creates the list of active hydro particles at this
138  * synchronization point, which is stored in the list TimeBinsHydro.ActiveParticleList[].
139  * This list is reused for the subsequent second and first hydro step. */
141 
142  /* if particles have increased their smoothing lengths, this is recorded in parent tree nodes */
144 
145  /* hydro-forces, second half-step. This will also update the predicted velocities/entropies with the new current ones */
147 
148  /* this does the closing gravity half-step for the timebins that end at the current synchronization point */
150 
151  /* do any extra physics, in a Strang-split way, for the timesteps that are finished */
152  calculate_non_standard_physics_end_of_step();
153 
154 #ifdef DEBUG_MD5
155  Logs.log_debug_md5("A");
156 #endif
157 
159 
161 
162 #ifdef DEBUG_MD5
163  Logs.log_debug_md5("BEFORE SNAP");
164 #endif
166 
167 #ifdef DEBUG_MD5
168  Logs.log_debug_md5("AFTER SNAP");
169 #endif
170 
171  if(All.Ti_Current >= TIMEBASE) /* did we reached the final time? */
172  {
173  mpi_printf("\nFinal time=%g reached. Simulation ends.\n", All.TimeMax);
174 
175  /* make a snapshot at the final time in case none has been produced at this time yet */
177  {
180  }
181 
182  break;
183  }
184 
185  /* kicks particles by half a gravity step */
187 
188 #ifdef DEBUG_MD5
189  Logs.log_debug_md5("B");
190 #endif
191 
192  /* Find new hydro timesteps. This will not change the set of active hydro particles at this synchronization point,
193  * but it can change how they are distributed over timebins. */
195 
196  /* compute hydro-forces and apply momentum changes to interacting particle pairs for first half-steps */
198 
199  /* update the neighbor tree with the new velocities */
201 
202  /* output some CPU usage log-info (accounts for everything needed up to complete the previous timestep) */
204 
205 #ifdef STOP_AFTER_STEP
206  if(All.NumCurrentTiStep == STOP_AFTER_STEP)
207  {
208  mpi_printf("RUN: We have reached the timestep specified with STOP_AFTER_STEP and therefore stop.");
209  endrun();
210  }
211 #endif
212 
214 
215  /* Check whether we should write a restart file */
216  if(check_for_interruption_of_run())
217  return;
218  }
219 
220  restart Restart{Communicator};
221  Restart.write(this); /* write a restart file at final time - can be used to continue simulation beyond final time */
222 
223  Logs.write_cpu_log(); /* output final cpu measurements */
224 }
225 
230 void sim::set_non_standard_physics_for_current_time(void)
231 {
232 #ifdef COOLING
233  CoolSfr.IonizeParams(); /* set UV background for the current time */
234 #endif
235 }
236 
243 void sim::calculate_non_standard_physics_end_of_step(void)
244 {
245 #ifdef COOLING
246 #ifdef STARFORMATION
247  CoolSfr.sfr_create_star_particles(&Sp);
248  CoolSfr.cooling_and_starformation(&Sp);
249 #else
250  CoolSfr.cooling_only(&Sp);
251 #endif
252 #endif
253 
254 #ifdef MEASURE_TOTAL_MOMENTUM
256 #endif
257 }
258 
268 int sim::check_for_interruption_of_run(void)
269 {
270  /* Check whether we need to interrupt the run */
271  int stopflag = 0;
272  if(ThisTask == 0)
273  {
274  FILE *fd;
275  char stopfname[MAXLEN_PATH_EXTRA];
276 
277  sprintf(stopfname, "%sstop", All.OutputDir);
278  if((fd = fopen(stopfname, "r"))) /* Is the stop-file present? If yes, interrupt the run. */
279  {
280  fclose(fd);
281  printf("stop-file detected. stopping.\n");
282  stopflag = 1;
283  unlink(stopfname);
284  }
285 
286  sprintf(stopfname, "%srestart", All.OutputDir);
287  if((fd = fopen(stopfname, "r"))) /* Is the restart-file present? If yes, write a user-requested restart file. */
288  {
289  fclose(fd);
290  printf("restart-file detected. writing restart files.\n");
291  stopflag = 3;
292  unlink(stopfname);
293  }
294 
295  if(Logs.CPUThisRun > 0.85 * All.TimeLimitCPU) /* are we running out of CPU-time ? If yes, interrupt run. */
296  {
297  printf("reaching time-limit. stopping.\n");
298  stopflag = 2;
299  }
300  }
301 
302  MPI_Bcast(&stopflag, 1, MPI_INT, 0, Communicator);
303 
304  if(stopflag)
305  {
306  restart Restart{Communicator};
307  Restart.write(this); /* write restart file */
308  MPI_Barrier(Communicator);
309 
310  if(stopflag == 3)
311  return 0;
312 
313  if(stopflag == 2 && ThisTask == 0)
314  {
315  FILE *fd;
316  char contfname[MAXLEN_PATH_EXTRA];
317  sprintf(contfname, "%scont", All.OutputDir);
318  if((fd = fopen(contfname, "w")))
319  fclose(fd);
320  }
321  return 1;
322  }
323 
324  /* is it time to write a regular restart-file? (for security) */
325  if(ThisTask == 0)
326  {
328  {
330  stopflag = 3;
331  }
332  else
333  stopflag = 0;
334  }
335 
336  MPI_Bcast(&stopflag, 1, MPI_INT, 0, Communicator);
337 
338  if(stopflag == 3)
339  {
340  restart Restart{Communicator};
341  Restart.write(this); /* write an occasional restart file */
342  stopflag = 0;
343  }
344  return 0;
345 }
346 
352 integertime sim::find_next_outputtime(integertime ti_curr)
353 {
354  integertime ti;
355  integertime ti_next = -1;
356  double time;
357 
359 
360  if(All.OutputListOn)
361  {
362  for(int i = 0; i < All.OutputListLength; i++)
363  {
364  time = All.OutputListTimes[i];
365 
366  if(time >= All.TimeBegin && time <= All.TimeMax)
367  {
369  ti = (integertime)(log(time / All.TimeBegin) / All.Timebase_interval);
370  else
371  ti = (integertime)((time - All.TimeBegin) / All.Timebase_interval);
372 
373 #ifndef OUTPUT_NON_SYNCHRONIZED_ALLOWED
374  /* We will now modify 'ti' to map it to the closest available output time according to the specified MaxSizeTimestep.
375  * The real output time may hence deviate by +/- 0.5*MaxSizeTimestep from the desired output time.
376  */
377 
378  /* first, determine maximum output interval based on All.MaxSizeTimestep */
380 
381  /* make it a power 2 subdivision */
382  integertime ti_min = TIMEBASE;
383  while(ti_min > timax)
384  ti_min >>= 1;
385  timax = ti_min;
386 
387  double multiplier = ti / ((double)timax);
388 
389  /* now round this to the nearest multiple of timax */
390  ti = ((integertime)(multiplier + 0.5)) * timax;
391 #endif
392 
393  if(ti >= ti_curr)
394  {
395  if(ti_next == -1)
396  {
397  ti_next = ti;
399  }
400 
401  if(ti_next > ti)
402  {
403  ti_next = ti;
405  }
406  }
407  }
408  }
409  }
410  else
411  {
413  {
414  if(All.TimeBetSnapshot <= 1.0)
415  Terminate("TimeBetSnapshot > 1.0 required for your simulation.\n");
416  }
417  else
418  {
419  if(All.TimeBetSnapshot <= 0.0)
420  Terminate("TimeBetSnapshot > 0.0 required for your simulation.\n");
421  }
422  time = All.TimeOfFirstSnapshot;
423  int iter = 0;
424 
425  while(time < All.TimeBegin)
426  {
428  time *= All.TimeBetSnapshot;
429  else
430  time += All.TimeBetSnapshot;
431 
432  iter++;
433 
434  if(iter > 1000000)
435  Terminate("Can't determine next output time.\n");
436  }
437 
438  while(time <= All.TimeMax)
439  {
441  ti = (integertime)(log(time / All.TimeBegin) / All.Timebase_interval);
442  else
443  ti = (integertime)((time - All.TimeBegin) / All.Timebase_interval);
444 
445 #ifndef OUTPUT_NON_SYNCHRONIZED_ALLOWED
446  /* We will now modify 'ti' to map it to the closest available output time according to the specified MaxSizeTimestep.
447  * The real output time may hence deviate by +/- 0.5*MaxSizeTimestep from the desired output time.
448  */
449 
450  /* first, determine maximum output interval based on All.MaxSizeTimestep */
452 
453  /* make it a power 2 subdivision */
454  integertime ti_min = TIMEBASE;
455  while(ti_min > timax)
456  ti_min >>= 1;
457  timax = ti_min;
458 
459  double multiplier = ti / ((double)timax);
460 
461  /* now round this to the nearest multiple of timax */
462  ti = ((integertime)(multiplier + 0.5)) * timax;
463 #endif
464 
465  if(ti >= ti_curr)
466  {
467  ti_next = ti;
468  break;
469  }
470 
472  time *= All.TimeBetSnapshot;
473  else
474  time += All.TimeBetSnapshot;
475 
476  iter++;
477 
478  if(iter > MAXITER)
479  Terminate("Can't determine next output time.\n");
480  }
481  }
482 
483  if(ti_next == -1)
484  {
485  ti_next = 2 * TIMEBASE; /* this will prevent any further output */
486 
487  mpi_printf("\nSNAPSHOT: There is no valid time for a further snapshot file.\n");
488  }
489  else
490  {
491  double next = All.get_absolutetime_from_integertime(ti_next);
492 
493  mpi_printf("\nSNAPSHOT: Setting next time for snapshot file to Time_next= %g (DumpFlag=%d)\n\n", next, All.DumpFlag_nextoutput);
494  }
495 
496  return ti_next;
497 }
498 
508 {
509 #if defined(LIGHTCONE_MASSMAPS)
510  /* we may do this on partial timesteps since for massmaps we always drift all particles, i.e. the lightcone is complete up to
511  * All.Time */
512  LightCone.lightcone_massmap_flush(1);
513 #endif
514 
515 #ifndef OUTPUT_NON_SYNCHRONIZED_ALLOWED
516  if(All.HighestActiveTimeBin == All.HighestOccupiedTimeBin) /* allow only top-level synchronization points */
517 #endif
519  {
520  for(int i = 0; i < Sp.NumPart; i++)
521  if(Sp.P[i].Ti_Current != All.Ti_Current)
522  Terminate("P[i].Ti_Current != All.Ti_Current");
523 
524 #if defined(STARFORMATION) && defined(FOF)
525  // do an extra domain decomposition here to make sure that there are no new stars among the block of gas particles
526  NgbTree.treefree();
530  NgbTree.treebuild(Sp.NumGas, NULL);
531 #endif
532 
533 #ifndef OUTPUT_NON_SYNCHRONIZED_ALLOWED
534  NgbTree.treefree();
537 #endif
538 
539 #ifdef FOF
540  mpi_printf("\nFOF: We shall first compute a group catalog for this snapshot file\n");
541 
542  /* this structure will hold auxiliary information for each particle, needed only during group finding */
543  Sp.PS = (subfind_data *)Mem.mymalloc_movable(&Sp.PS, "PS", Sp.MaxPart * sizeof(subfind_data));
544  memset(Sp.PS, 0, Sp.MaxPart * sizeof(subfind_data));
545 
546  /* First, we save the original location of the particles, in order to be able to revert to this layout later on */
547  for(int i = 0; i < Sp.NumPart; i++)
548  {
549  Sp.PS[i].OriginTask = ThisTask;
550  Sp.PS[i].OriginIndex = i;
551  }
552 
553  fof<simparticles> FoF{Communicator, &Sp, &Domain};
554 
555  FoF.fof_fof(All.SnapshotFileCount, "fof", "groups", 0);
556 
557 #if defined(MERGERTREE) && defined(SUBFIND)
558  MergerTree.CurrTotNsubhalos = FoF.TotNsubhalos;
559  MergerTree.CurrNsubhalos = FoF.Nsubhalos;
560 
561  MergerTree.mergertree_determine_descendants_on_the_fly(All.SnapshotFileCount);
562 
563  MergerTree.PrevTotNsubhalos = FoF.TotNsubhalos;
564  MergerTree.PrevNsubhalos = FoF.Nsubhalos;
565 
566  for(int n = 0; n < Sp.NumPart; n++)
567  {
568  Sp.P[n].PrevSubhaloNr = Sp.PS[n].SubhaloNr;
569  Sp.P[n].PrevSizeOfSubhalo = Sp.PS[n].SizeOfSubhalo;
570  Sp.P[n].PrevRankInSubhalo = Sp.PS[n].RankInSubhalo;
571 
572  if(Sp.P[n].PrevSubhaloNr.get() >= MergerTree.PrevTotNsubhalos && Sp.P[n].PrevSubhaloNr.get() != HALONR_MAX)
573  Terminate("Sp.P[n].PrevSubhaloNr=%lld MergerTree.PrevTotNsubhalos=%lld\n", (long long)Sp.P[n].PrevSubhaloNr.get(),
574  (long long)MergerTree.PrevTotNsubhalos);
575 
576  if(Sp.P[n].PrevSizeOfSubhalo.get() > 0 && Sp.P[n].PrevSubhaloNr.get() == HALONR_MAX)
577  Terminate("Sp.P[n].PrevSizeOfSubhalo=%d Sp.P[n].PrevSubhaloNr=%lld\n", (int)Sp.P[n].PrevSizeOfSubhalo.get(),
578  (long long)Sp.P[n].PrevSubhaloNr.get());
579  }
580 #endif
581 #endif
582 
584  {
585  snap_io Snap(&Sp, Communicator, All.SnapFormat); /* get an I/O object */
586  Snap.write_snapshot(All.SnapshotFileCount, NORMAL_SNAPSHOT); /* write snapshot file */
587  }
588 
589 #ifdef SUBFIND_ORPHAN_TREATMENT
590  {
592  Snap.write_snapshot(All.SnapshotFileCount, MOST_BOUND_PARTICLE_SNAPHOT); /* write special snapshot file */
593  }
594 #endif
595 
596 #ifdef FOF
597  /* now revert from output order to the original order */
598  for(int n = 0; n < Sp.NumPart; n++)
599  {
600  Sp.PS[n].TargetTask = Sp.PS[n].OriginTask;
601  Sp.PS[n].TargetIndex = Sp.PS[n].OriginIndex;
602  }
603 
604  TIMER_START(CPU_FOF);
605 
607 
608  TIMER_STOP(CPU_FOF);
609 
610  Mem.myfree(Sp.PS);
611 #endif
612 
613 #if defined(POWERSPEC_ON_OUTPUT) && defined(PERIODIC) && defined(PMGRID)
614  PM.calculate_power_spectra(All.SnapshotFileCount);
615 #endif
616 
618  All.Ti_nextoutput = find_next_outputtime(All.Ti_Current + 1);
619 
620 #ifndef OUTPUT_NON_SYNCHRONIZED_ALLOWED
623 
624  /* we need to reconstruct the timebins here. Even though the particles are in the same place again,
625  * it could have happened that Sp.P was reduced in size temporarily below NumPart on a certain task,
626  * in which case the timebin data may have become invalid.
627  */
629 
631  NgbTree.treebuild(Sp.NumGas, NULL);
632 #endif
633  }
634 
635 #if defined(LIGHTCONE_PARTICLES)
636  if(Lp.TestIfAboveFillFactor(std::min<int>(Lp.MaxPart, Sp.MaxPart)))
637  {
638 #if defined(LIGHTCONE_PARTICLES_GROUPS) && defined(FOF)
639  /* do this only on full timesteps if groups are calculated on lightcone */
641  {
642  mpi_printf("\nLIGHTCONE_PARTICLES_GROUPS: We shall first compute a group catalogue for the lightcone particles\n");
643 
644  /* assign unique IDs to Lp particles */
645 
646  int *numlist = (int *)Mem.mymalloc("numlist", Lp.NumPart * sizeof(int));
647 
648  MPI_Allgather(&Lp.NumPart, 1, MPI_INT, numlist, 1, MPI_INT, Communicator);
649 
650  long long newID = 1;
651  for(int i = 0; i < ThisTask; i++)
652  newID += numlist[i];
653 
654  for(int i = 0; i < Lp.NumPart; i++)
655  Lp.P[i].ID.set(newID++);
656 
657  Mem.myfree(numlist);
658 
659  domain<lcparticles> LcDomain(Communicator, &Lp);
660 
661  LcDomain.domain_decomposition(STANDARD);
662 
663  /* this structure will hold auxiliary information for each particle, needed only during group finding */
664  Lp.PS = (subfind_data *)Mem.mymalloc_movable(&Lp.PS, "PS", Lp.MaxPart * sizeof(subfind_data));
665  memset(Lp.PS, 0, Lp.MaxPart * sizeof(subfind_data));
666 
667  /* First, we save the original location of the particles, in order to be able to revert to this layout later on */
668  for(int i = 0; i < Lp.NumPart; i++)
669  {
670  Lp.PS[i].OriginTask = ThisTask;
671  Lp.PS[i].OriginIndex = i;
672  }
673 
674  fof<lcparticles> FoF{Communicator, &Lp, &LcDomain};
675 
676  double inner_distance = Driftfac.get_comoving_distance_for_scalefactor(All.Time);
677 
678  FoF.fof_fof(All.LightconeFileCount, "lc_fof", "lc_groups", inner_distance);
679 
680 #endif
681 
682  {
683 #ifdef MERGERTREE
684  MergerTree.Ntrees = 0;
685  lightcone_particle_io Lcone(&Lp, &LightCone, &MergerTree, Communicator, All.SnapFormat); /* get an I/O object */
686 #else
687  lightcone_particle_io Lcone(&Lp, &LightCone, Communicator, All.SnapFormat); /* get an I/O object */
688 #endif
689  long long NumLP_tot = Lp.NumPart;
690  MPI_Allreduce(MPI_IN_PLACE, &NumLP_tot, 1, MPI_LONG_LONG, MPI_SUM, Communicator);
691  mpi_printf("\nLIGHTCONE: writing particle lightcone conesnap files #%d ... (NumLP_tot = %lld)\n", All.LightconeFileCount,
692  NumLP_tot);
693 
694  for(int i = 0; i < Lp.NumPart; i++)
695  {
696  double pos[3];
697  Lp.signedintpos_to_pos((MySignedIntPosType *)Lp.P[i].IntPos, pos);
698  vec2pix_ring(LIGHTCONE_ORDER_NSIDE, pos, &Lp.P[i].ipnest);
699  }
700 
701 #if !defined(LIGHTCONE_PARTICLES_GROUPS)
702  /* let's now sort the lightcone_particle_data according to healpix pixel number */
703  mycxxsort_parallel(Lp.P, Lp.P + Lp.NumPart, Lp.compare_ipnest, Communicator);
704 #endif
705 
706  for(int conenr = 0; conenr < LightCone.Nlightcones; conenr++)
707  Lcone.lightcone_save(All.LightconeFileCount, conenr, false);
708 
709  mpi_printf("LIGHTCONE: done with writing files.\n");
710 
711  All.LightconeFileCount++;
712  /* let's put this in brackets so that the object's destructor will be called already here */
713  }
714 
715 #if defined(LIGHTCONE_PARTICLES_GROUPS) && defined(FOF)
716  /* now revert from output order to the original order */
717  for(int n = 0; n < Lp.NumPart; n++)
718  {
719  Lp.PS[n].TargetTask = Lp.PS[n].OriginTask;
720  Lp.PS[n].TargetIndex = Lp.PS[n].OriginIndex;
721  }
722 
723  TIMER_START(CPU_FOF);
724 
726 
727  Mem.myfree(Lp.PS);
728 
729  LcDomain.domain_free();
730 
731  TIMER_STOP(CPU_FOF);
732 
733  int ncount[2] = {0, 0};
734  long long nsum[2];
735  for(int n = 0; n < Lp.NumPart; n++)
736  {
737  if(Lp.P[n].getFlagSaveDistance())
738  {
739  Lp.P[n--] = Lp.P[--Lp.NumPart];
740  ncount[0]++;
741  }
742  else
743  {
744  ncount[1]++;
745  }
746  }
747 
748  sumup_large_ints(2, ncount, nsum, Communicator);
749  mpi_printf("LIGHTCONE_PARTICLES_GROUPS: We could store %lld particles from the buffer, but had to keep %lld\n", nsum[0],
750  nsum[1]);
751  }
752 #else
753  Lp.NumPart = 0;
754 #endif
755 
756  if(Lp.MaxPart > LIGHTCONE_ALLOC_FAC * Sp.MaxPart + 1 && Lp.NumPart < LIGHTCONE_ALLOC_FAC * Sp.MaxPart)
757  Lp.reallocate_memory_maxpart(LIGHTCONE_ALLOC_FAC * Sp.MaxPart);
758  }
759 #endif
760 }
global_data_all_processes All
Definition: main.cc:40
Definition: domain.h:31
void domain_free(void)
Definition: domain.cc:179
void particle_exchange_based_on_PS(MPI_Comm Communicator)
void domain_decomposition(domain_options mode)
Definition: domain.cc:51
double get_comoving_distance_for_scalefactor(double ascale)
Definition: driftfac.cc:199
int flush_everything(void)
Definition: logs.cc:129
void write_cpu_log(void)
Write the FdBalance and FdCPU files.
Definition: logs.cc:330
double CPUThisRun
Definition: logs.h:43
void output_log_messages(void)
Write the FdInfo and FdTimeBin files.
Definition: logs.cc:175
void compute_statistics(void)
Computes new global statistics if needed (done by energy_statistics())
Definition: logs.cc:501
void log_debug_md5(const char *msg)
void compute_total_momentum(void)
void update_maxhsml(void)
void update_velocities(void)
void mpi_printf(const char *fmt,...)
Definition: setcomm.h:55
int ThisTask
Definition: setcomm.h:33
MPI_Comm Communicator
Definition: setcomm.h:31
void run(void)
Definition: run.cc:49
domain< simparticles > Domain
Definition: simulation.h:58
void do_hydro_step_second_half(void)
Definition: kicks.cc:469
void do_hydro_step_first_half(void)
Definition: kicks.cc:449
void find_hydro_timesteps(void)
Definition: timestep.cc:39
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 do_gravity_step_second_half(void)
performs the second gravity half step kick operator
Definition: kicks.cc:274
void create_snapshot_if_desired(void)
Check if a snapshot should be saved.
Definition: run.cc:507
void endrun(void)
This function aborts the simulations.
Definition: begrun.cc:395
simparticles Sp
Definition: simulation.h:56
subfind_data * PS
Definition: simparticles.h:63
void reconstruct_timebins(void)
Definition: predict.cc:43
TimeBinData TimeBinsGravity
Definition: simparticles.h:205
particle_data * P
Definition: simparticles.h:54
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
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
void write_snapshot(int num, mysnaptype snap_type)
Save snapshot to disk.
Definition: snap_io.cc:559
void compute_densities(void)
Definition: density.cc:410
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 MAXLEN_PATH_EXTRA
Definition: constants.h:301
#define MAXITER
Definition: constants.h:305
int integertime
Definition: constants.h:331
#define TIMEBASE
Definition: constants.h:333
#define LIGHTCONE_ALLOC_FAC
Definition: constants.h:62
@ STANDARD
Definition: domain.h:23
@ MOST_BOUND_PARTICLE_SNAPHOT
Definition: dtypes.h:307
@ NORMAL_SNAPSHOT
Definition: dtypes.h:306
int32_t MySignedIntPosType
Definition: dtypes.h:36
#define HALONR_MAX
Definition: idstorage.h:20
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 log(half arg)
Definition: half.hpp:2745
double mycxxsort_parallel(T *begin, T *end, Comp comp, MPI_Comm comm)
void timebins_allocate(void)
Definition: timestep.cc:457
void timebins_free(void)
Definition: timestep.cc:472
double get_absolutetime_from_integertime(integertime ti)
Definition: allvars.h:200
integertime Ti_Current
Definition: allvars.h:188
int SmallestTimeBinWithDomainDecomposition
Definition: allvars.h:160
char OutputDir[MAXLEN_PATH]
Definition: allvars.h:272
integertime Ti_lastoutput
Definition: allvars.h:190
integertime Ti_nextoutput
Definition: allvars.h:189
double OutputListTimes[MAXLEN_OUTPUTLIST]
Definition: allvars.h:275
void set_cosmo_factors_for_current_time(void)
Definition: allvars.cc:20
char OutputListFlag[MAXLEN_OUTPUTLIST]
Definition: allvars.h:276
std::atomic< integertime > Ti_Current
Definition: particle_data.h:60