GADGET-4
sfr_eos.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 <math.h>
17 #include <mpi.h>
18 #include <stdio.h>
19 #include <stdlib.h>
20 #include <string.h>
21 
22 #include "../cooling_sfr/cooling.h"
23 #include "../data/allvars.h"
24 #include "../data/dtypes.h"
25 #include "../logs/logs.h"
26 #include "../system/system.h"
27 #include "../time_integration/timestep.h"
28 
38 void coolsfr::cooling_and_starformation(simparticles *Sp)
39 {
40  TIMER_START(CPU_COOLING_SFR);
41 
42  /* clear the SFR stored in the active timebins */
43  for(int bin = 0; bin < TIMEBINS; bin++)
44  if(Sp->TimeBinSynchronized[bin])
45  Sp->TimeBinSfr[bin] = 0;
46 
48 
49  gas_state gs = GasState;
50  do_cool_data dc = DoCoolData;
51 
52  for(int i = 0; i < Sp->TimeBinsHydro.NActiveParticles; i++)
53  {
54  int target = Sp->TimeBinsHydro.ActiveParticleList[i];
55  if(Sp->P[target].getType() == 0)
56  {
57  if(Sp->P[target].getMass() == 0 && Sp->P[target].ID.get() == 0)
58  continue; /* skip cells that have been swallowed or eliminated */
59 
60  double dens = Sp->SphP[target].Density;
61 
62  double dt =
63  (Sp->P[target].getTimeBinHydro() ? (((integertime)1) << Sp->P[target].getTimeBinHydro()) : 0) * All.Timebase_interval;
64  /* the actual time-step */
65 
66  double dtime = All.cf_atime * dt / All.cf_atime_hubble_a;
67 
68  /* check whether conditions for star formation are fulfilled.
69  *
70  * f=1 normal cooling
71  * f=0 star formation
72  */
73  int flag = 1; /* default is normal cooling */
74 
75  if(dens * All.cf_a3inv >= All.PhysDensThresh)
76  flag = 0;
77 
79  if(dens < All.OverDensThresh)
80  flag = 1;
81 
82  if(flag == 1) /* normal implicit isochoric cooling */
83  {
84  Sp->SphP[target].Sfr = 0;
85  cool_sph_particle(Sp, target, &gs, &dc);
86  }
87 
88  if(flag == 0) /* active star formation */
89  {
90  double tsfr = sqrt(All.PhysDensThresh / (dens * All.cf_a3inv)) * All.MaxSfrTimescale;
91 
92  double factorEVP = pow(dens * All.cf_a3inv / All.PhysDensThresh, -0.8) * All.FactorEVP;
93 
94  double egyhot = All.EgySpecSN / (1 + factorEVP) + All.EgySpecCold;
95 
96  double ne = Sp->SphP[target].Ne;
97 
98  double tcool = GetCoolingTime(egyhot, dens * All.cf_a3inv, &ne, &gs, &dc);
99  Sp->SphP[target].Ne = ne;
100 
101  double y = tsfr / tcool * egyhot / (All.FactorSN * All.EgySpecSN - (1 - All.FactorSN) * All.EgySpecCold);
102 
103  double x = 1 + 1 / (2 * y) - sqrt(1 / y + 1 / (4 * y * y));
104 
105  double egyeff = egyhot * (1 - x) + All.EgySpecCold * x;
106 
107  double cloudmass = x * Sp->P[target].getMass();
108  double utherm = Sp->get_utherm_from_entropy(target);
109 
110  if(tsfr < dtime)
111  tsfr = dtime;
112 
113  if(dt > 0)
114  {
115  if(Sp->P[target].getTimeBinHydro()) /* upon start-up, we need to protect against dt==0 */
116  {
117  double trelax = tsfr * (1 - x) / x / (All.FactorSN * (1 + factorEVP));
118  double egycurrent = utherm;
119 
120  double unew;
121  if(egycurrent > egyeff)
122  {
123  double dtcool = dtime;
124  ne = Sp->SphP[target].Ne;
125 
126  unew = DoCooling(egycurrent, dens * All.cf_a3inv, dtcool, &ne, &gs, &dc);
127 
128  if(unew < egyeff)
129  {
130  unew = egyeff;
131  }
132  }
133  else
134  unew = (egyeff + (egycurrent - egyeff) * exp(-dtime / trelax));
135 
136  double du = unew - utherm;
137  if(unew < All.MinEgySpec)
138  du = All.MinEgySpec - utherm;
139 
140  utherm += du;
141  Sp->set_entropy_from_utherm(utherm, target);
142  Sp->SphP[target].DtEntropy = 0.0;
143 
144 #ifdef OUTPUT_COOLHEAT
145  if(dtime > 0)
146  Sp->SphP[target].CoolHeat = du * Sp->P[target].getMass() / dtime;
147 #endif
148  Sp->SphP[target].set_thermodynamic_variables();
149  }
150  }
151 
152  if(utherm > 1.01 * egyeff)
153  Sp->SphP[target].Sfr = 0;
154  else
155  {
156  /* note that we convert the star formation rate to solar masses per year */
157  Sp->SphP[target].Sfr =
158  (1 - All.FactorSN) * cloudmass / tsfr * (All.UnitMass_in_g / SOLAR_MASS) / (All.UnitTime_in_s / SEC_PER_YEAR);
159  }
160 
161  Sp->TimeBinSfr[Sp->P[target].getTimeBinHydro()] += Sp->SphP[target].Sfr;
162  }
163  }
164  } /* end of main loop over active particles */
165 
166  TIMER_STOP(CPU_COOLING_SFR);
167 }
168 
175 void coolsfr::init_clouds(void)
176 {
177  gas_state gs = GasState;
178  do_cool_data dc = DoCoolData;
179 
180  if(All.PhysDensThresh == 0)
181  {
182  double A0 = All.FactorEVP;
183 
184  double egyhot = All.EgySpecSN / A0;
185 
186  double meanweight = 4 / (8 - 5 * (1 - HYDROGEN_MASSFRAC)); /* note: assuming FULL ionization */
187 
188  double u4 = 1 / meanweight * (1.0 / GAMMA_MINUS1) * (BOLTZMANN / PROTONMASS) * 1.0e4;
190 
191  double dens;
193  dens = 1.0e6 * 3 * All.Hubble * All.Hubble / (8 * M_PI * All.G);
194  else
195  dens = 1.0e6 * (1.0e-29 / All.UnitDensity_in_cgs);
196 
198  {
199  All.Time = 1.0; /* to be guaranteed to get z=0 rate */
201  IonizeParams();
202  }
203 
204  double ne = 1.0;
205  SetZeroIonization();
206 
207  double tcool = GetCoolingTime(egyhot, dens, &ne, &gs, &dc);
208 
209  double coolrate = egyhot / tcool / dens;
210 
211  double x = (egyhot - u4) / (egyhot - All.EgySpecCold);
212 
213  All.PhysDensThresh =
214  x / pow(1 - x, 2) * (All.FactorSN * All.EgySpecSN - (1 - All.FactorSN) * All.EgySpecCold) / (All.MaxSfrTimescale * coolrate);
215 
216  mpi_printf(
217  "\nA0= %g \nComputed: PhysDensThresh= %g (int units) %g h^2 cm^-3\nEXPECTED FRACTION OF COLD GAS AT THRESHOLD = "
218  "%g\n\ntcool=%g dens=%g egyhot=%g\n",
219  A0, All.PhysDensThresh, All.PhysDensThresh / (PROTONMASS / HYDROGEN_MASSFRAC / All.UnitDensity_in_cgs), x, tcool, dens,
220  egyhot);
221 
222  dens = All.PhysDensThresh * 10;
223 
224  double neff;
225  do
226  {
227  double tsfr = sqrt(All.PhysDensThresh / (dens)) * All.MaxSfrTimescale;
228  double factorEVP = pow(dens / All.PhysDensThresh, -0.8) * All.FactorEVP;
229  egyhot = All.EgySpecSN / (1 + factorEVP) + All.EgySpecCold;
230 
231  ne = 0.5;
232 
233  tcool = GetCoolingTime(egyhot, dens, &ne, &gs, &dc);
234 
235  double y = tsfr / tcool * egyhot / (All.FactorSN * All.EgySpecSN - (1 - All.FactorSN) * All.EgySpecCold);
236  x = 1 + 1 / (2 * y) - sqrt(1 / y + 1 / (4 * y * y));
237  double egyeff = egyhot * (1 - x) + All.EgySpecCold * x;
238 
239  double peff = GAMMA_MINUS1 * dens * egyeff;
240 
241  double fac = 1 / (log(dens * 1.025) - log(dens));
242  dens *= 1.025;
243 
244  neff = -log(peff) * fac;
245 
246  tsfr = sqrt(All.PhysDensThresh / (dens)) * All.MaxSfrTimescale;
247  factorEVP = pow(dens / All.PhysDensThresh, -0.8) * All.FactorEVP;
248  egyhot = All.EgySpecSN / (1 + factorEVP) + All.EgySpecCold;
249 
250  ne = 0.5;
251 
252  tcool = GetCoolingTime(egyhot, dens, &ne, &gs, &dc);
253 
254  y = tsfr / tcool * egyhot / (All.FactorSN * All.EgySpecSN - (1 - All.FactorSN) * All.EgySpecCold);
255  x = 1 + 1 / (2 * y) - sqrt(1 / y + 1 / (4 * y * y));
256  egyeff = egyhot * (1 - x) + All.EgySpecCold * x;
257 
258  peff = GAMMA_MINUS1 * dens * egyeff;
259 
260  neff += log(peff) * fac;
261  }
262  while(neff > 4.0 / 3);
263 
264  double thresholdStarburst = dens;
265 
266  mpi_printf("Run-away sets in for dens=%g\nDynamic range for quiescent star formation= %g\n", thresholdStarburst,
267  thresholdStarburst / All.PhysDensThresh);
268 
269  integrate_sfr();
270 
271  if(ThisTask == 0)
272  {
273  double sigma = 10.0 / All.Hubble * 1.0e-10 / pow(1.0e-3, 2);
274 
275  printf("Isotherm sheet central density: %g z0=%g\n", M_PI * All.G * sigma * sigma / (2 * GAMMA_MINUS1) / u4,
276  GAMMA_MINUS1 * u4 / (2 * M_PI * All.G * sigma));
277  myflush(stdout);
278  }
279 
281  {
282  All.Time = All.TimeBegin;
284  IonizeParams();
285  }
286  }
287 }
288 
300 void coolsfr::integrate_sfr(void)
301 {
302  double meanweight = 4 / (8 - 5 * (1 - HYDROGEN_MASSFRAC)); /* note: assuming FULL ionization */
303  double u4 = 1 / meanweight * (1.0 / GAMMA_MINUS1) * (BOLTZMANN / PROTONMASS) * 1.0e4;
305  gas_state gs = GasState;
306  do_cool_data dc = DoCoolData;
307 
309  {
310  All.Time = 1.0; /* to be guaranteed to get z=0 rate */
312  IonizeParams();
313  }
314 
315  FILE *fd = (WriteMiscFiles && (ThisTask == 0)) ? fopen("eos.txt", "w") : NULL;
316 
317  for(double rho = All.PhysDensThresh; rho <= 1000 * All.PhysDensThresh; rho *= 1.1)
318  {
319  double tsfr = sqrt(All.PhysDensThresh / rho) * All.MaxSfrTimescale;
320 
321  double factorEVP = pow(rho / All.PhysDensThresh, -0.8) * All.FactorEVP;
322 
323  double egyhot = All.EgySpecSN / (1 + factorEVP) + All.EgySpecCold;
324 
325  double ne = 1.0;
326 
327  double tcool = GetCoolingTime(egyhot, rho, &ne, &gs, &dc);
328 
329  double y = tsfr / tcool * egyhot / (All.FactorSN * All.EgySpecSN - (1 - All.FactorSN) * All.EgySpecCold);
330  double x = 1 + 1 / (2 * y) - sqrt(1 / y + 1 / (4 * y * y));
331 
332  double egyeff = egyhot * (1 - x) + All.EgySpecCold * x;
333 
334  double P0 = GAMMA_MINUS1 * rho * egyeff;
335 
336  if(WriteMiscFiles && (ThisTask == 0))
337  {
338  fprintf(fd, "%g %g\n", rho, P0);
339  }
340  }
341 
342  if(WriteMiscFiles && (ThisTask == 0))
343  {
344  fclose(fd);
345  fd = fopen("sfrrate.txt", "w");
346  }
347 
348  for(double rho0 = All.PhysDensThresh; rho0 <= 10000 * All.PhysDensThresh; rho0 *= 1.02)
349  {
350  double rho = rho0;
351  double q = 0;
352  double dz = 0.001;
353 
354  double sigma = 0, sigmasfr = 0, sigma_u4 = 0, x = 0;
355 
356  while(rho > 0.0001 * rho0)
357  {
358  double tsfr, P0, gam;
359  if(rho > All.PhysDensThresh)
360  {
361  tsfr = sqrt(All.PhysDensThresh / rho) * All.MaxSfrTimescale;
362 
363  double factorEVP = pow(rho / All.PhysDensThresh, -0.8) * All.FactorEVP;
364 
365  double egyhot = All.EgySpecSN / (1 + factorEVP) + All.EgySpecCold;
366 
367  double ne = 1.0;
368 
369  double tcool = GetCoolingTime(egyhot, rho, &ne, &gs, &dc);
370 
371  double y = tsfr / tcool * egyhot / (All.FactorSN * All.EgySpecSN - (1 - All.FactorSN) * All.EgySpecCold);
372  x = 1 + 1 / (2 * y) - sqrt(1 / y + 1 / (4 * y * y));
373 
374  double egyeff = egyhot * (1 - x) + All.EgySpecCold * x;
375 
376  P0 = GAMMA_MINUS1 * rho * egyeff;
377  double P1 = P0;
378 
379  double rho2 = 1.1 * rho;
380  double tsfr2 = sqrt(All.PhysDensThresh / rho2) * All.MaxSfrTimescale;
381  double factorEVP2 = pow(rho2 / All.PhysDensThresh, -0.8) * All.FactorEVP;
382  double egyhot2 = All.EgySpecSN / (1 + factorEVP2) + All.EgySpecCold;
383 
384  double tcool2 = GetCoolingTime(egyhot2, rho2, &ne, &gs, &dc);
385  double y2 = tsfr2 / tcool2 * egyhot2 / (All.FactorSN * All.EgySpecSN - (1 - All.FactorSN) * All.EgySpecCold);
386  double x2 = 1 + 1 / (2 * y2) - sqrt(1 / y2 + 1 / (4 * y2 * y2));
387  double egyeff2 = egyhot2 * (1 - x2) + All.EgySpecCold * x2;
388  double P2 = GAMMA_MINUS1 * rho2 * egyeff2;
389 
390  gam = log(P2 / P1) / log(rho2 / rho);
391  }
392  else
393  {
394  tsfr = 0;
395 
396  P0 = GAMMA_MINUS1 * rho * u4;
397  gam = 1.0;
398 
399  sigma_u4 += rho * dz;
400  }
401 
402  double drho = q;
403  double dq = -(gam - 2) / rho * q * q - 4 * M_PI * All.G / (gam * P0) * rho * rho * rho;
404 
405  sigma += rho * dz;
406  if(tsfr > 0)
407  {
408  sigmasfr += (1 - All.FactorSN) * rho * x / tsfr * dz;
409  }
410 
411  rho += drho * dz;
412  q += dq * dz;
413  }
414 
415  sigma *= 2; /* to include the other side */
416  sigmasfr *= 2;
417  sigma_u4 *= 2;
418 
422 
423  if(WriteMiscFiles && (ThisTask == 0))
424  {
425  fprintf(fd, "%g %g %g %g\n", rho0, sigma, sigmasfr, sigma_u4);
426  }
427  }
428 
430  {
431  All.Time = All.TimeBegin;
433  IonizeParams();
434  }
435 
436  if(WriteMiscFiles && (ThisTask == 0))
437  fclose(fd);
438 }
439 
442 void coolsfr::set_units_sfr(void)
443 {
444  All.OverDensThresh = All.CritOverDensity * All.OmegaBaryon * 3 * All.Hubble * All.Hubble / (8 * M_PI * All.G);
445 
446  All.PhysDensThresh = All.CritPhysDensity * PROTONMASS / HYDROGEN_MASSFRAC / All.UnitDensity_in_cgs;
447 
448  double meanweight = 4 / (1 + 3 * HYDROGEN_MASSFRAC); /* note: assuming NEUTRAL GAS */
449 
450  All.EgySpecCold = 1 / meanweight * (1.0 / GAMMA_MINUS1) * (BOLTZMANN / PROTONMASS) * All.TempClouds;
451  All.EgySpecCold *= All.UnitMass_in_g / All.UnitEnergy_in_cgs;
452 
453  meanweight = 4 / (8 - 5 * (1 - HYDROGEN_MASSFRAC)); /* note: assuming FULL ionization */
454 
455  All.EgySpecSN = 1 / meanweight * (1.0 / GAMMA_MINUS1) * (BOLTZMANN / PROTONMASS) * All.TempSupernova;
456  All.EgySpecSN *= All.UnitMass_in_g / All.UnitEnergy_in_cgs;
457 }
458 #endif
global_data_all_processes All
Definition: main.cc:40
MyIDType get(void) const
Definition: idstorage.h:87
double get_utherm_from_entropy(int i)
Definition: simparticles.h:238
void set_entropy_from_utherm(double utherm, int i)
Definition: simparticles.h:249
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 HYDROGEN_MASSFRAC
Definition: constants.h:112
#define GAMMA_MINUS1
Definition: constants.h:110
#define PROTONMASS
Definition: constants.h:124
#define SEC_PER_YEAR
Definition: constants.h:128
#define BOLTZMANN
Definition: constants.h:120
int integertime
Definition: constants.h:331
#define PARSEC
Definition: constants.h:123
#define M_PI
Definition: constants.h:56
#define TIMER_START(counter)
Starts the timer counter.
Definition: logs.h:197
#define TIMER_STOP(counter)
Stops the timer counter.
Definition: logs.h:220
expr exp(half arg)
Definition: half.hpp:2724
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
int NActiveParticles
Definition: timestep.h:18
int * ActiveParticleList
Definition: timestep.h:20
void set_cosmo_factors_for_current_time(void)
Definition: allvars.cc:20
MyDouble getMass(void)
unsigned char getTimeBinHydro(void)
MyIDStorage ID
Definition: particle_data.h:70
unsigned char getType(void)
void set_thermodynamic_variables(void)
void myflush(FILE *fstream)
Definition: system.cc:125