GADGET-4
init.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 <mpi.h>
15#include <algorithm>
16#include <cmath>
17#include <cstdio>
18#include <cstdlib>
19#include <cstring>
20
21#include "../cooling_sfr/cooling.h"
22#include "../data/allvars.h"
23#include "../data/dtypes.h"
24#include "../data/mymalloc.h"
25#include "../domain/domain.h"
26#include "../fof/fof.h"
27#include "../gravtree/gravtree.h"
28#include "../io/io.h"
29#include "../io/snap_io.h"
30#include "../logs/timer.h"
31#include "../main/main.h"
32#include "../main/simulation.h"
33#include "../mpi_utils/mpi_utils.h"
34#include "../ngbtree/ngbtree.h"
35#include "../ngenic/ngenic.h"
36#include "../pm/pm.h"
37#include "../sort/parallel_sort.h"
38#include "../subfind/subfind_readid_io.h"
39#include "../system/system.h"
40#include "../time_integration/timestep.h"
41
42using namespace std;
43
51void sim::init(int RestartSnapNum)
52{
53#ifdef NGENIC
55 {
56 Ngenic.ngenic_displace_particles();
57
59 {
60 double fac = 1 / sqrt(All.cf_a3inv);
61 for(int i = 0; i < Sp.NumPart; i++)
62 for(int k = 0; k < 3; k++)
63 Sp.P[i].Vel[k] *= fac;
64
65 strcat(All.SnapshotFileBase, "_ics");
66 mpi_printf("Start writing file %s\nRestartSnapNum %d\n", All.SnapshotFileBase, 0);
67 snap_io Snap(&Sp, Communicator, All.SnapFormat); /* get an I/O object */
69 endrun();
70 }
71 }
72#else
74 {
75 Terminate("Compile with option NGENIC to create cosmological initial conditions");
76 }
77#endif
78
79#ifdef LIGHTCONE_PARTICLES
80 Lp.MaxPart = LIGHTCONE_ALLOC_FAC * Sp.MaxPart;
81 Lp.NumPart = 0;
82 Lp.allocate_memory();
83#endif
84
85#ifdef LIGHTCONE_MASSMAPS
86 LightCone.Mp->Npix = nside2npix(All.LightConeMassMapsNside);
87
88 subdivide_evenly(LightCone.Mp->Npix, NTask, ThisTask, &LightCone.Mp->FirstPix, &LightCone.Mp->NpixLoc);
89
91 Mp.NumPart = 0;
92 Mp.allocate_memory();
93 LightCone.MassMap = (double *)Mem.mymalloc_movable_clear(&LightCone.MassMap, "MassMap", LightCone.Mp->NpixLoc * sizeof(double));
94#endif
95
96 /* this makes sure that masses are initialized in the case that the mass-block
97 is empty for this particle type */
98
99 for(int i = 0; i < Sp.NumPart; i++)
100 if(All.MassTable[Sp.P[i].getType()] != 0)
101 {
102#ifndef LEAN
103 Sp.P[i].setMass(All.MassTable[Sp.P[i].getType()]);
104#else
105 All.PartMass = All.MassTable[Sp.P[i].getType()];
106#endif
107 }
108
109#if NSOFTCLASSES > 1
110 for(int i = 0; i < Sp.NumPart; i++)
112#endif
113
114#ifdef GENERATE_GAS_IN_ICS
116 {
117 /* determine maximum ID */
118 MyIDType maxid = 0;
119 for(int i = 0; i < Sp.NumPart; i++)
120 if(Sp.P[i].ID.get() > maxid)
121 maxid = Sp.P[i].ID.get();
122
123 MyIDType *tmp = (MyIDType *)Mem.mymalloc("tmp", NTask * sizeof(MyIDType));
124 MPI_Allgather(&maxid, sizeof(MyIDType), MPI_BYTE, tmp, sizeof(MyIDType), MPI_BYTE, Communicator);
125
126 for(int i = 0; i < NTask; i++)
127 if(tmp[i] > maxid)
128 maxid = tmp[i];
129
130 Mem.myfree(tmp);
131
132 int count = 0;
133 for(int i = 0; i < Sp.NumPart; i++)
134#ifdef SPLIT_PARTICLE_TYPE
135 if((1 << Sp.P[i].getType()) & (SPLIT_PARTICLE_TYPE))
136#else
137 if(Sp.P[i].getType() == 1)
138#endif
139 count++;
140
141 int *numpart_list = (int *)Mem.mymalloc("numpart_list", NTask * sizeof(int));
142 MPI_Allgather(&count, 1, MPI_INT, numpart_list, 1, MPI_INT, Communicator);
143
144 maxid++;
145
146 for(int i = 0; i < ThisTask; i++)
147 maxid += numpart_list[i];
148
149 Mem.myfree(numpart_list);
150
151 Domain.domain_resize_storage(count + Sp.NumPart, count, 0);
152
153 memmove(Sp.P + count, Sp.P, sizeof(particle_data) * Sp.NumPart);
154
155 Sp.NumPart += count;
156 Sp.NumGas += count;
157
158 if(Sp.NumGas > Sp.MaxPartSph)
159 Terminate("Task=%d ends up getting more SPH particles (%d) than allowed (%d)\n", ThisTask, Sp.NumGas, Sp.MaxPartSph);
160
161 if(Sp.NumPart > Sp.MaxPart)
162 Terminate("Task=%d ends up getting more particles (%d) than allowed (%d)\n", ThisTask, Sp.NumPart, Sp.MaxPart);
163
164 double fac = All.OmegaBaryon / All.Omega0;
165 double rho = All.Omega0 * 3 * All.Hubble * All.Hubble / (8 * M_PI * All.G);
166
167 int j = 0;
168 for(int i = count; i < Sp.NumPart; i++)
169#ifdef SPLIT_PARTICLE_TYPE
170 if((1 << Sp.P[i].getType()) & (SPLIT_PARTICLE_TYPE))
171#else
172 if(Sp.P[i].getType() == 1)
173#endif
174 {
175 double d = pow(Sp.P[i].getMass() / rho, 1.0 / 3);
176 double a = 0.5 * All.OmegaBaryon / All.Omega0 * d;
177 double b = 0.5 * (All.Omega0 - All.OmegaBaryon) / All.Omega0 * d;
178
179 MyIntPosType delta_a[3];
180 double aa[3] = {a, a, a};
182
183 MyIntPosType delta_b[3];
184 double bb[3] = {b, b, b};
186
187 Sp.P[j] = Sp.P[i];
188
189 Sp.P[j].setMass(Sp.P[j].getMass() * fac);
190 Sp.P[i].setMass(Sp.P[i].getMass() * (1 - fac));
191
192#ifdef SECOND_ORDER_LPT_ICS
193 Sp.P[j].OldAcc *= fac;
194 Sp.P[i].OldAcc *= (1 - fac);
195#endif
196
197 Sp.P[j].setType(0);
198#if NSOFTCLASSES > 1
200#endif
201 Sp.P[j].ID.set(maxid++);
202 Sp.P[i].IntPos[0] += delta_a[0];
203 Sp.P[i].IntPos[1] += delta_a[1];
204 Sp.P[i].IntPos[2] += delta_a[2];
205 Sp.P[j].IntPos[0] -= delta_b[0];
206 Sp.P[j].IntPos[1] -= delta_b[1];
207 Sp.P[j].IntPos[2] -= delta_b[2];
208
209 j++;
210 }
211
212 All.MassTable[0] = 0;
213
214#ifdef SPLIT_PARTICLE_TYPE
215 for(int i = 1; i < NTYPES; i++)
216 if((1 << i) & (SPLIT_PARTICLE_TYPE))
217 All.MassTable[i] *= (1 - fac);
218#else
219 All.MassTable[1] *= (1 - fac);
220#endif
221
222 mpi_printf("\nGENERATE_GAS_IN_ICS: Generated gas particles from DM particle distribution. TotNumGas=%lld\n\n", Sp.TotNumGas);
223 }
224#endif
225
226#ifdef STARFORMATION
228 {
229 if(All.MassTable[STAR_TYPE] == 0 && All.MassTable[0] > 0)
230 {
231 All.MassTable[0] = 0;
232 }
233 }
234#endif
235
236 double u_init = (1.0 / GAMMA_MINUS1) * (BOLTZMANN / PROTONMASS) * All.InitGasTemp;
237 u_init *= All.UnitMass_in_g / All.UnitEnergy_in_cgs; /* unit conversion */
238
239 double molecular_weight;
240 if(All.InitGasTemp > 1.0e4) /* assuming FULL ionization */
241 molecular_weight = 4 / (8 - 5 * (1 - HYDROGEN_MASSFRAC));
242 else /* assuming NEUTRAL GAS */
243 molecular_weight = 4 / (1 + 3 * HYDROGEN_MASSFRAC);
244
245 u_init /= molecular_weight;
246
247 All.InitGasU = u_init;
248
250 {
251 if(All.InitGasTemp > 0)
252 {
253 for(int i = 0; i < Sp.NumGas; i++)
254 {
255 if(ThisTask == 0 && i == 0 && Sp.SphP[i].Entropy == 0)
256 mpi_printf("READIC: Initializing u from InitGasTemp !\n");
257
258 if(Sp.SphP[i].Entropy == 0)
259 Sp.SphP[i].Entropy = All.InitGasU;
260 /* Note: the coversion to entropy will be done in the function init(),
261 after the densities have been computed */
262 }
263 }
264 }
265
266 for(int i = 0; i < Sp.NumGas; i++)
267 Sp.SphP[i].Entropy = std::max<double>(All.MinEgySpec, Sp.SphP[i].Entropy);
268
269#ifdef COOLING
270 CoolSfr.IonizeParams();
271#endif
272
274 {
276 All.Ti_Current = 0;
277 }
278 else
279 {
281 All.Ti_Current = 0;
282 }
283
285
286 All.NumCurrentTiStep = 0; /* setup some counters */
288
290 {
291 if(RestartSnapNum < 0)
292 All.SnapshotFileCount = atoi(All.InitCondFile + strlen(All.InitCondFile) - 3) + 1;
293 else
294 All.SnapshotFileCount = RestartSnapNum + 1;
295 }
296
299 All.TotNumDensity = 0;
300 All.TotNumHydro = 0;
301
302 All.TopNodeAllocFactor = 0.08;
303 All.TreeAllocFactor = 0.3;
305
307
308#if defined(EVALPOTENTIAL) && defined(PMGRID) && defined(PERIODIC)
309 double mass_sum = 0;
310
311 for(int i = 0; i < Sp.NumPart; i++)
312 mass_sum += Sp.P[i].getMass();
313
314 MPI_Allreduce(&mass_sum, &All.TotalMass, 1, MPI_DOUBLE, MPI_SUM, Communicator);
315#endif
316
317 if(All.ComovingIntegrationOn) /* change to new velocity variable */
318 {
319 double afac = sqrt(All.Time) * All.Time;
320
321 for(int i = 0; i < Sp.NumPart; i++)
322 {
323 for(int j = 0; j < 3; j++)
324 Sp.P[i].Vel[j] *= afac; /* for dm/gas particles, p = a^2 xdot */
325 }
326 }
327
328 for(int i = 0; i < TIMEBINS; i++)
329 All.Ti_begstep[i] = 0;
330
331#if defined(PMGRID) && !defined(TREEPM_NOTIMESPLIT)
332 All.PM_Ti_endstep = All.PM_Ti_begstep = 0;
333#endif
334
335 for(int i = 0; i < Sp.NumPart; i++) /* start-up initialization with non-zero values where required */
336 {
337#ifndef LEAN
338 Sp.P[i].access.clear();
339#endif
340#ifdef MERGERTREE
341 Sp.P[i].PrevSubhaloNr.set(HALONR_MAX);
342#endif
343 }
344
345 for(int i = 0; i < TIMEBINS; i++)
346 Sp.TimeBinSynchronized[i] = 1;
347
348 Sp.reconstruct_timebins();
349
350 for(int i = 0; i < Sp.NumGas; i++) /* initialize sph_properties with non-zero values where required */
351 {
352 Sp.SphP[i].EntropyPred = Sp.SphP[i].Entropy;
353
354 for(int j = 0; j < 3; j++)
355 Sp.SphP[i].VelPred[j] = Sp.P[i].Vel[j];
356
358 {
359#ifdef COOLING
360 Sp.SphP[i].Ne = 1.0;
361#endif
362 }
363
364#ifdef PRESSURE_ENTROPY_SPH
365 Sp.SphP[i].EntropyToInvGammaPred = pow(Sp.SphP[i].EntropyPred, 1.0 / GAMMA);
366#endif
367
368#ifdef TIMEDEP_ART_VISC
369#ifdef HIGH_ART_VISC_START
370 Sp.SphP[i].Alpha = All.ArtBulkViscConst;
371#else
372 Sp.SphP[i].Alpha = All.AlphaMin;
373#endif
374#endif
375 }
376
377#ifdef RECREATE_UNIQUE_IDS
378 recreate_unique_ids();
379#endif
380
381 test_id_uniqueness();
382
383 Domain.domain_decomposition(STANDARD); /* do initial domain decomposition (gives equal numbers of particles) */
384
385 GravTree.set_softenings();
386
387#ifdef ADAPTIVE_HYDRO_SOFTENING
388 mpi_printf("INIT: Adaptive hydro softening, minimum gravitational softening for SPH particles: %g\n",
389 All.MinimumComovingHydroSoftening);
390 mpi_printf("INIT: Adaptive hydro softening, maximum gravitational softening for SPH particles: %g\n",
391 All.MinimumComovingHydroSoftening * pow(All.AdaptiveHydroSofteningSpacing, NSOFTCLASSES_HYDRO - 1));
392 mpi_printf("INIT: Adaptive hydro softening, number of softening values: %d\n", NSOFTCLASSES_HYDRO);
393#endif
394
395#ifdef INDIVIDUAL_GRAVITY_SOFTENING
396 Sp.init_individual_softenings();
397#endif
398
399 if(All.RestartFlag == RST_FOF)
400 {
401#ifdef FOF
402
403#if defined(SUBFIND) && defined(MERGERTREE)
404 // we are reading the previous subhalo catalogue, if available, to assign the previous subhalo length to particles
405 MergerTree.get_previous_size_of_subhalo_for_each_particle(RestartSnapNum - 1);
406#endif
407
408 Sp.PS = (subfind_data *)Mem.mymalloc_movable(&Sp.PS, "PS", Sp.MaxPart * sizeof(subfind_data));
409 memset(Sp.PS, 0, Sp.MaxPart * sizeof(subfind_data));
410
411 for(int i = 0; i < Sp.NumGas; i++)
412 {
413 if(ThisTask == 0 && i == 0)
414 printf("INIT: Converting u -> entropy All.cf_a3inv=%g\n", All.cf_a3inv);
415
416 Sp.SphP[i].Entropy = GAMMA_MINUS1 * Sp.SphP[i].Entropy / pow(Sp.SphP[i].Density * All.cf_a3inv, GAMMA_MINUS1);
417 }
418
419 /* First, we save the original location of the particles, in order to be able to revert to this layout later on */
420 for(int i = 0; i < Sp.NumPart; i++)
421 {
422 Sp.PS[i].OriginTask = ThisTask;
423 Sp.PS[i].OriginIndex = i;
424 }
425 fof<simparticles> FoF{Communicator, &Sp, &Domain};
426 FoF.fof_fof(RestartSnapNum, "fof", "groups", 0);
427
428 {
430 snap_io Snap(&Sp, Communicator, All.SnapFormat); /* get an I/O object */
431 Snap.write_snapshot(RestartSnapNum, NORMAL_SNAPSHOT);
432 }
433
434#ifdef SUBFIND_ORPHAN_TREATMENT
435 {
436 /* now read the IDs of the most bound particles of a previously existing special dump, with the idea being to
437 * be able also on postprocessing to construct a cumulative list of all particles that used to be a most particles
438 * in any previous dump. This can be accomplished by computing the group catalogues consecutively from 000 to the last
439 * snapshoyt number
440 */
441
442 if(RestartSnapNum > 0)
443 {
444 subreadid_io SnapIDread(&Sp.IdStore, Communicator, All.SnapFormat);
445 SnapIDread.previously_bound_read_snap_ids(RestartSnapNum - 1);
446
447 FoF.subfind_match_ids_of_previously_most_bound_ids(&Sp);
448 }
449
450 snap_io Snap(&Sp, Communicator, All.SnapFormat);
451 Snap.write_snapshot(RestartSnapNum, MOST_BOUND_PARTICLE_SNAPHOT); /* write special snapshot file */
452 }
453#endif
454
455#endif
456 endrun();
457 }
458
459 /* build neighbor tree */
460 NgbTree.treeallocate(Sp.NumGas, &Sp, &Domain);
461 NgbTree.treebuild(Sp.NumGas, NULL);
462
464 {
465#if defined(PMGRID) && defined(PERIODIC)
466
467 PM.calculate_power_spectra(RestartSnapNum);
468#else
469 mpi_printf("\nThis option (Power Spectrum) only works for PERIODIC and enabled PMGRID.\n\n");
470#endif
471 endrun();
472 }
473
474 All.Ti_Current = 0;
475
476 setup_smoothinglengths();
477
478 /* at this point, the entropy variable actually contains the
479 * internal energy, read in from the initial conditions file.
480 * Once the density has been computed, we can convert to entropy.
481 */
482#ifdef PRESSURE_ENTROPY_SPH
483#ifndef INITIAL_CONDITIONS_CONTAIN_ENTROPY
484 NgbTree.setup_entropy_to_invgamma();
485#endif
486#endif
487
488 double mass = 0;
489 for(int i = 0; i < Sp.NumGas; i++)
490 {
491#ifndef INITIAL_CONDITIONS_CONTAIN_ENTROPY
492
493 if(ThisTask == 0 && i == 0)
494 printf("INIT: Converting u -> entropy\n");
495
496#if !defined(PRESSURE_ENTROPY_SPH) && !defined(ISOTHERM_EQS)
497 Sp.SphP[i].Entropy = GAMMA_MINUS1 * Sp.SphP[i].Entropy / pow(Sp.SphP[i].Density * All.cf_a3inv, GAMMA_MINUS1);
498#endif
499 Sp.SphP[i].EntropyPred = Sp.SphP[i].Entropy;
500
501#endif
502 /* The predicted entropy values have been already set for all SPH formulation */
503 /* so it should be ok computing pressure and csound now */
504 Sp.SphP[i].set_thermodynamic_variables();
505
506 mass += Sp.P[i].getMass();
507 }
508
510 {
511#ifdef PERIODIC
514 {
515 /* can't do this check when not all particles are loaded */
516 check_omega();
517 }
518 else
519 {
520 mpi_printf("Skipping Omega check since not all particles are loaded\n");
521 }
522#endif
523 }
524
525#ifdef STARFORMATION
526 /* initialize absolute masses in materials */
527 for(int i = 0; i < Sp.NumGas; i++)
528 {
529 Sp.SphP[i].Metallicity = Sp.P[i].Metallicity; // set above
530
531 Sp.SphP[i].MassMetallicity = Sp.SphP[i].Metallicity * Sp.P[i].getMass();
532 }
533#endif
534
535 // tree_based_timesteps_set_initialmaxtistep();
536
537#ifdef DEBUG_MD5
538 Logs.log_debug_md5("AFTER-INIT");
539#endif
540
541 return;
542}
543
544#ifdef PERIODIC
550void sim::check_omega(void)
551{
552 double mass = 0;
553
554 for(int i = 0; i < Sp.NumPart; i++)
555 mass += Sp.P[i].getMass();
556
557 double masstot;
558 MPI_Allreduce(&mass, &masstot, 1, MPI_DOUBLE, MPI_SUM, Communicator);
559
560 double omega = masstot * (LONG_Z * LONG_Y * LONG_Z) / (All.BoxSize * All.BoxSize * All.BoxSize) /
561 (3 * All.Hubble * All.Hubble / (8 * M_PI * All.G));
562 if(fabs(omega - All.Omega0) > 1.0e-2)
563 {
565 "\n\nI've found something odd!\nThe mass content accounts only for Omega=%g,\nbut you specified Omega=%g in the "
566 "parameterfile.\n\nI better stop.\n",
567 omega, All.Omega0);
568 endrun();
569 }
570}
571#endif
572
581void sim::setup_smoothinglengths(void)
582{
584
585 for(int i = 0; i < Sp.NumGas; i++)
587
589
591 {
592 mpi_printf("INIT: Setup smoothing lengths.\n");
593
596
597 for(int i = 0; i < Sp.NumGas; i++)
598 {
599 int no = GravTree.Father[i];
600
601 while(10 * All.DesNumNgb * Sp.P[i].getMass() > GravTree.get_nodep(no)->mass)
602 {
603 int p = GravTree.get_nodep(no)->father;
604
605 if(p < 0)
606 break;
607
608 no = p;
609 }
610
611 double len;
612 if(GravTree.get_nodep(no)->level > 0)
614 else
615 len = Sp.RegionLen;
616
617 Sp.SphP[i].Hsml = pow(3.0 / (4 * M_PI) * All.DesNumNgb * Sp.P[i].getMass() / GravTree.get_nodep(no)->mass, 1.0 / 3) * len;
618 }
619
621 for(int i = 0; i < Sp.NumGas; i++)
623
625
626#ifdef PRESSURE_ENTROPY_SPH
627 for(int i = 0; i < Sp.NumGas; i++)
628 Sp.SphP[i].PressureSphDensity = Sp.SphP[i].Density;
629#endif
630
632 }
633}
634
635void sim::recreate_unique_ids(void)
636{
637 mpi_printf("INIT: Setting new unique IDs.\n");
638
639 int *numpart_list = (int *)Mem.mymalloc("numpart_list", NTask * sizeof(int));
640
641 MPI_Allgather(&Sp.NumPart, 1, MPI_INT, numpart_list, 1, MPI_INT, Communicator);
642
643 MyIDType id = 1;
644
645 for(int i = 0; i < ThisTask; i++)
646 id += numpart_list[i];
647
648 for(int i = 0; i < Sp.NumPart; i++)
649 Sp.P[i].ID.set(id++);
650
651 Mem.myfree(numpart_list);
652}
653
659void sim::test_id_uniqueness(void)
660{
661 mpi_printf("INIT: Testing ID uniqueness...\n");
662
663 double t0 = Logs.second();
664
665 MyIDType *ids = (MyIDType *)Mem.mymalloc("ids", (Sp.NumPart + 1) * sizeof(MyIDType));
666 MyIDType *ids_first = (MyIDType *)Mem.mymalloc("ids_first", NTask * sizeof(MyIDType));
667 int *num_list = (int *)Mem.mymalloc("num_list", NTask * sizeof(int));
668
669 for(int i = 0; i < Sp.NumPart; i++)
670 ids[i] = Sp.P[i].ID.get();
671
673
674 for(int i = 1; i < Sp.NumPart; i++)
675 {
676 if(ids[i] == ids[i - 1])
677 Terminate("non-unique ID=%lld found on task=%d (i=%d Sp.NumPart=%d type=%d)\n", (long long)ids[i], ThisTask, i, Sp.NumPart,
678 Sp.P[i].getType());
679 }
680
681 MPI_Allgather(&ids[0], sizeof(MyIDType), MPI_BYTE, ids_first, sizeof(MyIDType), MPI_BYTE, Communicator);
682 MPI_Allgather(&Sp.NumPart, 1, MPI_INT, num_list, 1, MPI_INT, Communicator);
683
684 int next_non_empty_task = ThisTask + 1;
685
686 while(next_non_empty_task < NTask)
687 if(num_list[next_non_empty_task] == 0)
688 next_non_empty_task++;
689 else
690 break;
691
692 if(Sp.NumPart > 0 && next_non_empty_task < NTask)
693 {
694 if(ids[Sp.NumPart - 1] == ids_first[next_non_empty_task])
695 Terminate("non-unique ID=%lld found on task=%d\n", (long long)ids[Sp.NumPart - 1], ThisTask);
696 }
697
698 Mem.myfree(num_list);
699 Mem.myfree(ids_first);
700 Mem.myfree(ids);
701
702 double t1 = Logs.second();
703
704 mpi_printf("INIT: success. took=%g sec\n\n", Logs.timediff(t0, t1));
705}
global_data_all_processes All
Definition: main.cc:40
MyIDType get(void) const
Definition: idstorage.h:87
void set(MyIDType ID)
Definition: idstorage.h:96
MySignedIntPosType pos_to_signedintpos(T posdiff)
MyReal RegionLen
Definition: intposconvert.h:39
MyReal FacIntToCoord
Definition: intposconvert.h:37
double timediff(double t0, double t1)
Definition: logs.cc:488
void log_debug_md5(const char *msg)
double second(void)
Definition: logs.cc:471
void mpi_printf(const char *fmt,...)
Definition: setcomm.h:55
int ThisTask
Definition: setcomm.h:33
int NTask
Definition: setcomm.h:32
MPI_Comm Communicator
Definition: setcomm.h:31
domain< simparticles > Domain
Definition: simulation.h:58
gwalk GravTree
Definition: simulation.h:65
void init(int RestartSnapNum)
Prepares the loaded initial conditions for the run.
Definition: init.cc:51
sph NgbTree
Definition: simulation.h:68
void endrun(void)
This function aborts the simulations.
Definition: begrun.cc:395
simparticles Sp
Definition: simulation.h:56
long long TotNumPart
Definition: simparticles.h:46
long long TotNumGas
Definition: simparticles.h:47
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
static bool compare_IDs(const MyIDType &a, const MyIDType &b)
Definition: simparticles.h:72
void write_snapshot(int num, mysnaptype snap_type)
Save snapshot to disk.
Definition: snap_io.cc:559
void density(int *targetlist, int ntarget)
Definition: density.cc:434
void treeallocate(int max_partindex, partset *Pptr, domain< partset > *Dptr)
Definition: tree.cc:776
void treefree(void)
Definition: tree.cc:1232
node * get_nodep(int no)
Definition: tree.h:334
int treebuild(int ninsert, int *indexlist)
Definition: tree.cc:52
int * Father
Definition: tree.h:104
#define TIMEBINS
Definition: constants.h:332
#define NSOFTCLASSES_HYDRO
Definition: constants.h:321
#define HYDROGEN_MASSFRAC
Definition: constants.h:112
#define GAMMA_MINUS1
Definition: constants.h:110
#define PROTONMASS
Definition: constants.h:124
#define GAMMA
Definition: constants.h:99
#define BOLTZMANN
Definition: constants.h:120
#define NTYPES
Definition: constants.h:308
#define LIGHTCONE_MASSMAP_ALLOC_FAC
Definition: constants.h:66
#define STAR_TYPE
Definition: constants.h:351
#define TIMEBASE
Definition: constants.h:333
#define LIGHTCONE_ALLOC_FAC
Definition: constants.h:62
#define M_PI
Definition: constants.h:56
@ STANDARD
Definition: domain.h:23
@ MOST_BOUND_PARTICLE_SNAPHOT
Definition: dtypes.h:307
@ NORMAL_SNAPSHOT
Definition: dtypes.h:306
#define LONG_Y
Definition: dtypes.h:369
uint32_t MyIntPosType
Definition: dtypes.h:35
int32_t MySignedIntPosType
Definition: dtypes.h:36
@ RST_STARTFROMSNAP
Definition: dtypes.h:315
@ RST_POWERSPEC
Definition: dtypes.h:317
@ RST_CREATEICS
Definition: dtypes.h:319
@ RST_FOF
Definition: dtypes.h:316
@ RST_BEGIN
Definition: dtypes.h:313
@ RST_RESUME
Definition: dtypes.h:314
#define BITS_FOR_POSITIONS
Definition: dtypes.h:37
unsigned int MyIDType
Definition: dtypes.h:68
#define LONG_Z
Definition: dtypes.h:377
#define HALONR_MAX
Definition: idstorage.h:20
logs Logs
Definition: main.cc:43
#define Terminate(...)
Definition: macros.h:15
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
expr sqrt(half arg)
Definition: half.hpp:2777
STL namespace.
double mycxxsort_parallel(T *begin, T *end, Comp comp, MPI_Comm comm)
int NActiveParticles
Definition: timestep.h:18
int * ActiveParticleList
Definition: timestep.h:20
long long GlobalNActiveParticles
Definition: timestep.h:19
unsigned char level
Definition: tree.h:64
int father
Definition: tree.h:59
int SofteningClassOfPartType[NTYPES]
Definition: allvars.h:250
enum restart_options RestartFlag
Definition: allvars.h:68
double MassTable[NTYPES]
Definition: allvars.h:269
integertime Ti_Current
Definition: allvars.h:188
integertime Ti_begstep[TIMEBINS]
Definition: allvars.h:210
long long TotNumDirectForces
Definition: allvars.h:104
char SnapshotFileBase[MAXLEN_PATH]
Definition: allvars.h:272
void set_cosmo_factors_for_current_time(void)
Definition: allvars.cc:20
long long TotNumOfForces
Definition: allvars.h:103
char InitCondFile[MAXLEN_PATH]
Definition: allvars.h:272
MyDouble mass
Definition: gravtree.h:65
void setSofteningClass(unsigned char softclass)
MyDouble getMass(void)
void setMass(MyDouble mass)
MyFloat Vel[3]
Definition: particle_data.h:54
MyIDStorage ID
Definition: particle_data.h:70
unsigned char getType(void)
void setType(unsigned char type)
MyIntPosType IntPos[3]
Definition: particle_data.h:53
void subdivide_evenly(long long N, int pieces, int index_bin, long long *first, int *count)
Definition: system.cc:51