GADGET-4
starformation.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#ifdef STARFORMATION
15
16#include <assert.h>
17#include <math.h>
18#include <mpi.h>
19#include <stdio.h>
20#include <stdlib.h>
21#include <string.h>
22
23#include "../cooling_sfr/cooling.h"
24#include "../data/allvars.h"
25#include "../data/dtypes.h"
26#include "../data/mymalloc.h"
27#include "../logs/logs.h"
28#include "../logs/timer.h"
29#include "../system/system.h"
30#include "../time_integration/timestep.h"
31
40void coolsfr::sfr_create_star_particles(simparticles *Sp)
41{
42 TIMER_START(CPU_COOLING_SFR);
43
44 double dt, dtime;
45 MyDouble mass_of_star;
46 double sum_sm, total_sm, rate, sum_mass_stars, total_sum_mass_stars;
47 double p = 0, pall = 0, prob, p_decide;
48 double rate_in_msunperyear;
49 double totsfrrate;
50 double w = 0;
51
53
54 stars_spawned = stars_converted = 0;
55
56 sum_sm = sum_mass_stars = 0;
57
58 for(int i = 0; i < Sp->TimeBinsHydro.NActiveParticles; i++)
59 {
60 int target = Sp->TimeBinsHydro.ActiveParticleList[i];
61 if(Sp->P[target].getType() == 0)
62 {
63 if(Sp->P[target].getMass() == 0 && Sp->P[target].ID.get() == 0)
64 continue; /* skip cells that have been swallowed or eliminated */
65
66 dt = (Sp->P[target].getTimeBinHydro() ? (((integertime)1) << Sp->P[target].getTimeBinHydro()) : 0) * All.Timebase_interval;
67 /* the actual time-step */
68
69 dtime = All.cf_atime * dt / All.cf_atime_hubble_a;
70
71 mass_of_star = 0;
72 prob = 0;
73 p = 0;
74
75 if(Sp->SphP[target].Sfr > 0)
76 {
77 p = Sp->SphP[target].Sfr / ((All.UnitMass_in_g / SOLAR_MASS) / (All.UnitTime_in_s / SEC_PER_YEAR)) * dtime /
78 Sp->P[target].getMass();
79 pall = p;
80 sum_sm += Sp->P[target].getMass() * (1 - exp(-p));
81
83
84 Sp->SphP[target].Metallicity += w * METAL_YIELD * (1 - exp(-p));
85 Sp->SphP[target].MassMetallicity = Sp->SphP[target].Metallicity * Sp->P[target].getMass();
86 Sp->P[target].Metallicity = Sp->SphP[target].Metallicity;
87
88 mass_of_star = Sp->P[target].getMass();
89
90 prob = Sp->P[target].getMass() / mass_of_star * (1 - exp(-pall));
91 }
92
93 if(prob == 0)
94 continue;
95
96 if(prob < 0)
97 Terminate("prob < 0");
98
99 /* decide what process to consider (currently available: make a star or kick to wind) */
100 p_decide = get_random_number();
101
102 if(p_decide < p / pall) /* ok, a star formation is considered */
103 make_star(Sp, target, prob, mass_of_star, &sum_mass_stars);
104
105 if(Sp->SphP[target].Sfr > 0)
106 {
107 if(Sp->P[target].getType() == 0) /* to protect using a particle that has been turned into a star */
108 {
109 Sp->SphP[target].Metallicity += (1 - w) * METAL_YIELD * (1 - exp(-p));
110 Sp->SphP[target].MassMetallicity = Sp->SphP[target].Metallicity * Sp->P[target].getMass();
111 }
112 }
113 Sp->P[target].Metallicity = Sp->SphP[target].Metallicity;
114 }
115 } /* end of main loop over active gas particles */
116
117 MPI_Allreduce(&stars_spawned, &tot_stars_spawned, 1, MPI_INT, MPI_SUM, Communicator);
118 MPI_Allreduce(&stars_converted, &tot_stars_converted, 1, MPI_INT, MPI_SUM, Communicator);
119
120 if(tot_stars_spawned > 0 || tot_stars_converted > 0)
121 {
122 mpi_printf("SFR: spawned %d stars, converted %d gas particles into stars\n", tot_stars_spawned, tot_stars_converted);
123 }
124
125 tot_altogether_spawned = tot_stars_spawned;
126 altogether_spawned = stars_spawned;
127 if(tot_altogether_spawned)
128 {
129 /* need to assign new unique IDs to the spawned stars */
130
131 if(All.MaxID == 0) /* MaxID not calculated yet */
132 {
133 /* determine maximum ID */
134 MyIDType maxid = 0;
135 for(int i = 0; i < Sp->NumPart; i++)
136 if(Sp->P[i].ID.get() > maxid)
137 {
138 maxid = Sp->P[i].ID.get();
139 }
140
141 MyIDType *tmp = (MyIDType *)Mem.mymalloc("tmp", NTask * sizeof(MyIDType));
142
143 MPI_Allgather(&maxid, sizeof(MyIDType), MPI_BYTE, tmp, sizeof(MyIDType), MPI_BYTE, Communicator);
144
145 for(int i = 0; i < NTask; i++)
146 if(tmp[i] > maxid)
147 maxid = tmp[i];
148
149 All.MaxID = maxid;
150
151 Mem.myfree(tmp);
152 }
153
154 int *list = (int *)Mem.mymalloc("list", NTask * sizeof(int));
155
156 MPI_Allgather(&altogether_spawned, 1, MPI_INT, list, 1, MPI_INT, Communicator);
157
158 MyIDType newid = All.MaxID + 1;
159
160 for(int i = 0; i < ThisTask; i++)
161 newid += list[i];
162
163 Mem.myfree(list);
164
165 for(int i = 0; i < altogether_spawned; i++)
166 Sp->P[Sp->NumPart + i].ID.set(newid++);
167
168 All.MaxID += tot_altogether_spawned;
169 }
170
171 /* Note: New tree construction can be avoided because of `force_add_star_to_tree()' */
172 if(tot_stars_spawned > 0 || tot_stars_converted > 0)
173 {
174 Sp->TotNumPart += tot_stars_spawned;
175 Sp->TotNumGas -= tot_stars_converted;
176 Sp->NumPart += stars_spawned;
177 }
178
179 double sfrrate = 0;
180 for(int bin = 0; bin < TIMEBINS; bin++)
181 if(Sp->TimeBinsHydro.TimeBinCount[bin])
182 sfrrate += Sp->TimeBinSfr[bin];
183
184 MPI_Allreduce(&sfrrate, &totsfrrate, 1, MPI_DOUBLE, MPI_SUM, Communicator);
185
186 MPI_Reduce(&sum_sm, &total_sm, 1, MPI_DOUBLE, MPI_SUM, 0, Communicator);
187 MPI_Reduce(&sum_mass_stars, &total_sum_mass_stars, 1, MPI_DOUBLE, MPI_SUM, 0, Communicator);
188 if(ThisTask == 0)
189 {
190 if(All.TimeStep > 0)
191 rate = total_sm / (All.TimeStep / All.cf_atime_hubble_a);
192 else
193 rate = 0;
194
195 /* compute the cumulative mass of stars */
196 cum_mass_stars += total_sum_mass_stars;
197
198 /* convert to solar masses per yr */
199 rate_in_msunperyear = rate * (All.UnitMass_in_g / SOLAR_MASS) / (All.UnitTime_in_s / SEC_PER_YEAR);
200
201 fprintf(Logs.FdSfr, "%14e %14e %14e %14e %14e %14e\n", All.Time, total_sm, totsfrrate, rate_in_msunperyear, total_sum_mass_stars,
202 cum_mass_stars);
203 myflush(Logs.FdSfr);
204 }
205
206 TIMER_STOP(CPU_COOLING_SFR);
207}
208
219void coolsfr::convert_sph_particle_into_star(simparticles *Sp, int i, double birthtime)
220{
221 Sp->P[i].setType(STAR_TYPE);
222#if NSOFTCLASSES > 1
224#endif
225#ifdef INDIVIDUAL_GRAVITY_SOFTENING
226 if(((1 << Sp->P[i].getType()) & (INDIVIDUAL_GRAVITY_SOFTENING)))
227 Sp->P[i].setSofteningClass(Sp->get_softening_type_from_mass(Sp->P[i].getMass()));
228#endif
229
230 Sp->TimeBinSfr[Sp->P[i].getTimeBinHydro()] -= Sp->SphP[i].Sfr;
231
232 Sp->P[i].StellarAge = birthtime;
233
234 return;
235}
236
251void coolsfr::spawn_star_from_sph_particle(simparticles *Sp, int igas, double birthtime, int istar, MyDouble mass_of_star)
252{
253 Sp->P[istar] = Sp->P[igas];
254 Sp->P[istar].setType(STAR_TYPE);
255#if NSOFTCLASSES > 1
256 Sp->P[istar].setSofteningClass(All.SofteningClassOfPartType[Sp->P[istar].getType()]);
257#endif
258#ifdef INDIVIDUAL_GRAVITY_SOFTENING
259 if(((1 << Sp->P[istar].getType()) & (INDIVIDUAL_GRAVITY_SOFTENING)))
260 Sp->P[istar].setSofteningClass(Sp->get_softening_type_from_mass(Sp->P[istar].getMass()));
261#endif
262
264
265 Sp->TimeBinsGravity.timebin_add_particle(istar, igas, Sp->P[istar].TimeBinGrav, Sp->TimeBinSynchronized[Sp->P[istar].TimeBinGrav]);
266
267 Sp->P[istar].setMass(mass_of_star);
268
269 Sp->P[istar].StellarAge = birthtime;
270
271 /* now change the conserved quantities in the cell in proportion */
272 double fac = (Sp->P[igas].getMass() - Sp->P[istar].getMass()) / Sp->P[igas].getMass();
273
274 Sp->P[igas].setMass(fac * Sp->P[igas].getMass());
275
276 return;
277}
278
291void coolsfr::make_star(simparticles *Sp, int i, double prob, MyDouble mass_of_star, double *sum_mass_stars)
292{
293 if(mass_of_star > Sp->P[i].getMass())
294 Terminate("mass_of_star > P[i].Mass");
295
296 if(get_random_number() < prob)
297 {
298 if(mass_of_star == Sp->P[i].getMass())
299 {
300 /* here we turn the gas particle itself into a star particle */
301 stars_converted++;
302
303 *sum_mass_stars += Sp->P[i].getMass();
304
305 convert_sph_particle_into_star(Sp, i, All.Time);
306 }
307 else
308 {
309 /* in this case we spawn a new star particle, only reducing the mass in the cell by mass_of_star */
310 altogether_spawned = stars_spawned;
311 if(Sp->NumPart + altogether_spawned >= Sp->MaxPart)
312 Terminate("NumPart=%d spwawn %d particles no space left (Sp.MaxPart=%d)\n", Sp->NumPart, altogether_spawned, Sp->MaxPart);
313
314 int j = Sp->NumPart + altogether_spawned; /* index of new star */
315
316 spawn_star_from_sph_particle(Sp, i, All.Time, j, mass_of_star);
317
318 *sum_mass_stars += mass_of_star;
319 stars_spawned++;
320 }
321 }
322}
323
324#endif /* closes SFR */
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
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
int TimeBinSynchronized[TIMEBINS]
Definition: simparticles.h:203
#define TIMEBINS
Definition: constants.h:332
#define SOLAR_MASS
Definition: constants.h:119
#define METAL_YIELD
Definition: constants.h:114
#define SEC_PER_YEAR
Definition: constants.h:128
#define STAR_TYPE
Definition: constants.h:351
int integertime
Definition: constants.h:331
float MyDouble
Definition: dtypes.h:87
unsigned int MyIDType
Definition: dtypes.h:68
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:15
memory Mem
Definition: main.cc:44
expr exp(half arg)
Definition: half.hpp:2724
int NActiveParticles
Definition: timestep.h:18
int * ActiveParticleList
Definition: timestep.h:20
int TimeBinCount[TIMEBINS]
Definition: timestep.h:21
void timebin_add_particle(int i_new, int i_old, int timeBin, int addToListOfActiveParticles)
Definition: timestep.cc:599
int SofteningClassOfPartType[NTYPES]
Definition: allvars.h:250
void set_cosmo_factors_for_current_time(void)
Definition: allvars.cc:20
void setSofteningClass(unsigned char softclass)
MyDouble getMass(void)
void setMass(MyDouble mass)
unsigned char getTimeBinHydro(void)
MyIDStorage ID
Definition: particle_data.h:70
signed char TimeBinGrav
Definition: particle_data.h:71
unsigned char getType(void)
void setType(unsigned char type)
void myflush(FILE *fstream)
Definition: system.cc:125
double get_random_number(void)
Definition: system.cc:49