GADGET-4
cooling.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 COOLING
15 
16 #include <gsl/gsl_math.h>
17 #include <math.h>
18 #include <mpi.h>
19 #include <stdio.h>
20 #include <stdlib.h>
21 #include <string.h>
22 #include <algorithm>
23 
24 #include "../cooling_sfr/cooling.h"
25 #include "../data/allvars.h"
26 #include "../data/dtypes.h"
27 #include "../data/mymalloc.h"
28 #include "../logs/logs.h"
29 #include "../logs/timer.h"
30 #include "../system/system.h"
31 #include "../time_integration/timestep.h"
32 
46 double coolsfr::DoCooling(double u_old, double rho, double dt, double *ne_guess, gas_state *gs, do_cool_data *DoCool)
47 {
48  DoCool->u_old_input = u_old;
49  DoCool->rho_input = rho;
50  DoCool->dt_input = dt;
51  DoCool->ne_guess_input = *ne_guess;
52 
53  if(!gsl_finite(u_old))
54  Terminate("invalid input: u_old=%g\n", u_old);
55 
56  if(u_old < 0 || rho < 0)
57  Terminate("invalid input: u_old=%g rho=%g dt=%g All.MinEgySpec=%g\n", u_old, rho, dt, All.MinEgySpec);
58 
59  rho *= All.UnitDensity_in_cgs * All.HubbleParam * All.HubbleParam; /* convert to physical cgs units */
62 
63  gs->nHcgs = gs->XH * rho / PROTONMASS; /* hydrogen number dens in cgs units */
64  double ratefact = gs->nHcgs * gs->nHcgs / rho;
65 
66  double u = u_old;
67  double u_lower = u;
68  double u_upper = u;
69 
70  double LambdaNet = CoolingRateFromU(u, rho, ne_guess, gs, DoCool);
71 
72  /* bracketing */
73 
74  if(u - u_old - ratefact * LambdaNet * dt < 0) /* heating */
75  {
76  u_upper *= sqrt(1.1);
77  u_lower /= sqrt(1.1);
78  while(u_upper - u_old - ratefact * CoolingRateFromU(u_upper, rho, ne_guess, gs, DoCool) * dt < 0)
79  {
80  u_upper *= 1.1;
81  u_lower *= 1.1;
82  }
83  }
84 
85  if(u - u_old - ratefact * LambdaNet * dt > 0)
86  {
87  u_lower /= sqrt(1.1);
88  u_upper *= sqrt(1.1);
89  while(u_lower - u_old - ratefact * CoolingRateFromU(u_lower, rho, ne_guess, gs, DoCool) * dt > 0)
90  {
91  u_upper /= 1.1;
92  u_lower /= 1.1;
93  }
94  }
95 
96  int iter = 0;
97  double du;
98  do
99  {
100  u = 0.5 * (u_lower + u_upper);
101 
102  LambdaNet = CoolingRateFromU(u, rho, ne_guess, gs, DoCool);
103 
104  if(u - u_old - ratefact * LambdaNet * dt > 0)
105  {
106  u_upper = u;
107  }
108  else
109  {
110  u_lower = u;
111  }
112 
113  du = u_upper - u_lower;
114 
115  iter++;
116 
117  if(iter >= (MAXITER - 10))
118  printf("u= %g\n", u);
119  }
120  while(fabs(du / u) > 1.0e-6 && iter < MAXITER);
121 
122  if(iter >= MAXITER)
123  Terminate(
124  "failed to converge in DoCooling(): DoCool->u_old_input=%g\nDoCool->rho_input= %g\nDoCool->dt_input= "
125  "%g\nDoCool->ne_guess_input= %g\n",
126  DoCool->u_old_input, DoCool->rho_input, DoCool->dt_input, DoCool->ne_guess_input);
127 
128  u *= All.UnitDensity_in_cgs / All.UnitPressure_in_cgs; /* to internal units */
129 
130  return u;
131 }
132 
141 double coolsfr::GetCoolingTime(double u_old, double rho, double *ne_guess, gas_state *gs, do_cool_data *DoCool)
142 {
143  DoCool->u_old_input = u_old;
144  DoCool->rho_input = rho;
145  DoCool->ne_guess_input = *ne_guess;
146 
147  rho *= All.UnitDensity_in_cgs * All.HubbleParam * All.HubbleParam; /* convert to physical cgs units */
149 
150  gs->nHcgs = gs->XH * rho / PROTONMASS; /* hydrogen number dens in cgs units */
151  double ratefact = gs->nHcgs * gs->nHcgs / rho;
152 
153  double u = u_old;
154 
155  double LambdaNet = CoolingRateFromU(u, rho, ne_guess, gs, DoCool);
156 
157  /* bracketing */
158 
159  if(LambdaNet >= 0) /* ups, we have actually heating due to UV background */
160  return 0;
161 
162  double coolingtime = u_old / (-ratefact * LambdaNet);
163 
164  coolingtime *= All.HubbleParam / All.UnitTime_in_s;
165 
166  return coolingtime;
167 }
168 
180 double coolsfr::convert_u_to_temp(double u, double rho, double *ne_guess, gas_state *gs, const do_cool_data *DoCool)
181 {
182  double u_input = u;
183  double rho_input = rho;
184  double ne_input = *ne_guess;
185 
186  double mu = (1 + 4 * gs->yhelium) / (1 + gs->yhelium + *ne_guess);
187  double temp = GAMMA_MINUS1 / BOLTZMANN * u * PROTONMASS * mu;
188 
189  double max = 0;
190  int iter = 0;
191  double temp_old;
192  do
193  {
194  double ne_old = *ne_guess;
195 
196  find_abundances_and_rates(log10(temp), rho, ne_guess, gs, DoCool);
197  temp_old = temp;
198 
199  mu = (1 + 4 * gs->yhelium) / (1 + gs->yhelium + *ne_guess);
200 
201  double temp_new = GAMMA_MINUS1 / BOLTZMANN * u * PROTONMASS * mu;
202 
203  max = std::max<double>(max, temp_new / (1 + gs->yhelium + *ne_guess) * fabs((*ne_guess - ne_old) / (temp_new - temp_old + 1.0)));
204 
205  temp = temp_old + (temp_new - temp_old) / (1 + max);
206  iter++;
207 
208  if(iter > (MAXITER - 10))
209  printf("-> temp= %g ne=%g\n", temp, *ne_guess);
210  }
211  while(fabs(temp - temp_old) > 1.0e-3 * temp && iter < MAXITER);
212 
213  if(iter >= MAXITER)
214  {
215  printf("failed to converge in convert_u_to_temp()\n");
216  printf("u_input= %g\nrho_input=%g\n ne_input=%g\n", u_input, rho_input, ne_input);
217  printf("DoCool->u_old_input=%g\nDoCool->rho_input= %g\nDoCool->dt_input= %g\nDoCool->ne_guess_input= %g\n", DoCool->u_old_input,
218  DoCool->rho_input, DoCool->dt_input, DoCool->ne_guess_input);
219  Terminate("convergence failure");
220  }
221 
222  return temp;
223 }
224 
233 void coolsfr::find_abundances_and_rates(double logT, double rho, double *ne_guess, gas_state *gs, const do_cool_data *DoCool)
234 {
235  double logT_input = logT;
236  double rho_input = rho;
237  double ne_input = *ne_guess;
238 
239  if(!gsl_finite(logT))
240  Terminate("logT=%g\n", logT);
241 
242  if(logT <= Tmin) /* everything neutral */
243  {
244  gs->nH0 = 1.0;
245  gs->nHe0 = gs->yhelium;
246  gs->nHp = 0;
247  gs->nHep = 0;
248  gs->nHepp = 0;
249  gs->ne = 0;
250  *ne_guess = 0;
251  return;
252  }
253 
254  if(logT >= Tmax) /* everything is ionized */
255  {
256  gs->nH0 = 0;
257  gs->nHe0 = 0;
258  gs->nHp = 1.0;
259  gs->nHep = 0;
260  gs->nHepp = gs->yhelium;
261  gs->ne = gs->nHp + 2.0 * gs->nHepp;
262  *ne_guess = gs->ne; /* note: in units of the hydrogen number density */
263  return;
264  }
265 
266  double t = (logT - Tmin) / deltaT;
267  int j = (int)t;
268  double fhi = t - j;
269  double flow = 1 - fhi;
270 
271  if(*ne_guess == 0)
272  *ne_guess = 1.0;
273 
274  gs->nHcgs = gs->XH * rho / PROTONMASS; /* hydrogen number dens in cgs units */
275 
276  gs->ne = *ne_guess;
277  double neold = gs->ne;
278  int niter = 0;
279  gs->necgs = gs->ne * gs->nHcgs;
280 
281  /* evaluate number densities iteratively (cf KWH eqns 33-38) in units of nH */
282  do
283  {
284  niter++;
285 
286  gs->aHp = flow * RateT[j].AlphaHp + fhi * RateT[j + 1].AlphaHp;
287  gs->aHep = flow * RateT[j].AlphaHep + fhi * RateT[j + 1].AlphaHep;
288  gs->aHepp = flow * RateT[j].AlphaHepp + fhi * RateT[j + 1].AlphaHepp;
289  gs->ad = flow * RateT[j].Alphad + fhi * RateT[j + 1].Alphad;
290  gs->geH0 = flow * RateT[j].GammaeH0 + fhi * RateT[j + 1].GammaeH0;
291  gs->geHe0 = flow * RateT[j].GammaeHe0 + fhi * RateT[j + 1].GammaeHe0;
292  gs->geHep = flow * RateT[j].GammaeHep + fhi * RateT[j + 1].GammaeHep;
293 
294  if(gs->necgs <= 1.e-25 || pc.J_UV == 0)
295  {
296  gs->gJH0ne = gs->gJHe0ne = gs->gJHepne = 0;
297  }
298  else
299  {
300  gs->gJH0ne = pc.gJH0 / gs->necgs;
301  gs->gJHe0ne = pc.gJHe0 / gs->necgs;
302  gs->gJHepne = pc.gJHep / gs->necgs;
303  }
304 
305  gs->nH0 = gs->aHp / (gs->aHp + gs->geH0 + gs->gJH0ne); /* eqn (33) */
306  gs->nHp = 1.0 - gs->nH0; /* eqn (34) */
307 
308  if((gs->gJHe0ne + gs->geHe0) <= SMALLNUM) /* no ionization at all */
309  {
310  gs->nHep = 0.0;
311  gs->nHepp = 0.0;
312  gs->nHe0 = gs->yhelium;
313  }
314  else
315  {
316  gs->nHep = gs->yhelium /
317  (1.0 + (gs->aHep + gs->ad) / (gs->geHe0 + gs->gJHe0ne) + (gs->geHep + gs->gJHepne) / gs->aHepp); /* eqn (35) */
318  gs->nHe0 = gs->nHep * (gs->aHep + gs->ad) / (gs->geHe0 + gs->gJHe0ne); /* eqn (36) */
319  gs->nHepp = gs->nHep * (gs->geHep + gs->gJHepne) / gs->aHepp; /* eqn (37) */
320  }
321 
322  neold = gs->ne;
323 
324  gs->ne = gs->nHp + gs->nHep + 2 * gs->nHepp; /* eqn (38) */
325  gs->necgs = gs->ne * gs->nHcgs;
326 
327  if(pc.J_UV == 0)
328  break;
329 
330  double nenew = 0.5 * (gs->ne + neold);
331  gs->ne = nenew;
332  gs->necgs = gs->ne * gs->nHcgs;
333 
334  if(fabs(gs->ne - neold) < 1.0e-4)
335  break;
336 
337  if(niter > (MAXITER - 10))
338  printf("ne= %g niter=%d\n", gs->ne, niter);
339  }
340  while(niter < MAXITER);
341 
342  if(niter >= MAXITER)
343  Terminate(
344  "no convergence reached in find_abundances_and_rates(): logT_input= %g rho_input= %g ne_input= %g "
345  "DoCool->u_old_input=%g\nDoCool->rho_input= %g\nDoCool->dt_input= %g\nDoCool->ne_guess_input= %g\n",
346  logT_input, rho_input, ne_input, DoCool->u_old_input, DoCool->rho_input, DoCool->dt_input, DoCool->ne_guess_input);
347 
348  gs->bH0 = flow * RateT[j].BetaH0 + fhi * RateT[j + 1].BetaH0;
349  gs->bHep = flow * RateT[j].BetaHep + fhi * RateT[j + 1].BetaHep;
350  gs->bff = flow * RateT[j].Betaff + fhi * RateT[j + 1].Betaff;
351 
352  *ne_guess = gs->ne;
353 }
354 
365 double coolsfr::CoolingRateFromU(double u, double rho, double *ne_guess, gas_state *gs, const do_cool_data *DoCool)
366 {
367  double temp = convert_u_to_temp(u, rho, ne_guess, gs, DoCool);
368 
369  return CoolingRate(log10(temp), rho, ne_guess, gs, DoCool);
370 }
371 
382 double coolsfr::AbundanceRatios(double u, double rho, double *ne_guess, double *nH0_pointer, double *nHeII_pointer)
383 {
384  gas_state gs = GasState;
385  do_cool_data DoCool = DoCoolData;
386  DoCool.u_old_input = u;
387  DoCool.rho_input = rho;
388  DoCool.ne_guess_input = *ne_guess;
389 
390  rho *= All.UnitDensity_in_cgs * All.HubbleParam * All.HubbleParam; /* convert to physical cgs units */
392 
393  double temp = convert_u_to_temp(u, rho, ne_guess, &gs, &DoCool);
394 
395  *nH0_pointer = gs.nH0;
396  *nHeII_pointer = gs.nHep;
397 
398  return temp;
399 }
400 
408 double coolsfr::CoolingRate(double logT, double rho, double *nelec, gas_state *gs, const do_cool_data *DoCool)
409 {
410  double Lambda, Heat;
411 
412  if(logT <= Tmin)
413  logT = Tmin + 0.5 * deltaT; /* floor at Tmin */
414 
415  gs->nHcgs = gs->XH * rho / PROTONMASS; /* hydrogen number dens in cgs units */
416 
417  if(logT < Tmax)
418  {
419  find_abundances_and_rates(logT, rho, nelec, gs, DoCool);
420 
421  /* Compute cooling and heating rate (cf KWH Table 1) in units of nH**2 */
422  double T = pow(10.0, logT);
423 
424  double LambdaExcH0 = gs->bH0 * gs->ne * gs->nH0;
425  double LambdaExcHep = gs->bHep * gs->ne * gs->nHep;
426  double LambdaExc = LambdaExcH0 + LambdaExcHep; /* excitation */
427  double LambdaIonH0 = 2.18e-11 * gs->geH0 * gs->ne * gs->nH0;
428  double LambdaIonHe0 = 3.94e-11 * gs->geHe0 * gs->ne * gs->nHe0;
429  double LambdaIonHep = 8.72e-11 * gs->geHep * gs->ne * gs->nHep;
430  double LambdaIon = LambdaIonH0 + LambdaIonHe0 + LambdaIonHep; /* ionization */
431  double LambdaRecHp = 1.036e-16 * T * gs->ne * (gs->aHp * gs->nHp);
432  double LambdaRecHep = 1.036e-16 * T * gs->ne * (gs->aHep * gs->nHep);
433  double LambdaRecHepp = 1.036e-16 * T * gs->ne * (gs->aHepp * gs->nHepp);
434  double LambdaRecHepd = 6.526e-11 * gs->ad * gs->ne * gs->nHep;
435  double LambdaRec = LambdaRecHp + LambdaRecHep + LambdaRecHepp + LambdaRecHepd;
436  double LambdaFF = gs->bff * (gs->nHp + gs->nHep + 4 * gs->nHepp) * gs->ne;
437  Lambda = LambdaExc + LambdaIon + LambdaRec + LambdaFF;
438 
440  {
441  double redshift = 1 / All.Time - 1;
442  double LambdaCmptn = 5.65e-36 * gs->ne * (T - 2.73 * (1. + redshift)) * pow(1. + redshift, 4.) / gs->nHcgs;
443 
444  Lambda += LambdaCmptn;
445  }
446 
447  Heat = 0;
448  if(pc.J_UV != 0)
449  Heat += (gs->nH0 * pc.epsH0 + gs->nHe0 * pc.epsHe0 + gs->nHep * pc.epsHep) / gs->nHcgs;
450  }
451  else /* here we're outside of tabulated rates, T>Tmax K */
452  {
453  /* at high T (fully ionized); only free-free and Compton cooling are present. Assumes no heating. */
454 
455  Heat = 0;
456 
457  /* very hot: H and He both fully ionized */
458  gs->nHp = 1.0;
459  gs->nHep = 0;
460  gs->nHepp = gs->yhelium;
461  gs->ne = gs->nHp + 2.0 * gs->nHepp;
462  *nelec = gs->ne; /* note: in units of the hydrogen number density */
463 
464  double T = pow(10.0, logT);
465  double LambdaFF = 1.42e-27 * sqrt(T) * (1.1 + 0.34 * exp(-(5.5 - logT) * (5.5 - logT) / 3)) * (gs->nHp + 4 * gs->nHepp) * gs->ne;
466  double LambdaCmptn;
468  {
469  double redshift = 1 / All.Time - 1;
470  /* add inverse Compton cooling off the microwave background */
471  LambdaCmptn = 5.65e-36 * gs->ne * (T - 2.73 * (1. + redshift)) * pow(1. + redshift, 4.) / gs->nHcgs;
472  }
473  else
474  LambdaCmptn = 0;
475 
476  Lambda = LambdaFF + LambdaCmptn;
477  }
478 
479  return (Heat - Lambda);
480 }
481 
486 void coolsfr::MakeRateTable(void)
487 {
488  GasState.yhelium = (1 - GasState.XH) / (4 * GasState.XH);
489  GasState.mhboltz = PROTONMASS / BOLTZMANN;
490 
491  deltaT = (Tmax - Tmin) / NCOOLTAB;
492  GasState.ethmin = pow(10.0, Tmin) * (1. + GasState.yhelium) / ((1. + 4. * GasState.yhelium) * GasState.mhboltz * GAMMA_MINUS1);
493  /* minimum internal energy for neutral gas */
494 
495  for(int i = 0; i <= NCOOLTAB; i++)
496  {
497  RateT[i].BetaH0 = RateT[i].BetaHep = RateT[i].Betaff = RateT[i].AlphaHp = RateT[i].AlphaHep = RateT[i].AlphaHepp =
498  RateT[i].Alphad = RateT[i].GammaeH0 = RateT[i].GammaeHe0 = RateT[i].GammaeHep = 0;
499 
500  double T = pow(10.0, Tmin + deltaT * i);
501  double Tfact = 1.0 / (1 + sqrt(T / 1.0e5));
502 
503  /* collisional excitation */
504  /* Cen 1992 */
505  if(118348 / T < 70)
506  RateT[i].BetaH0 = 7.5e-19 * exp(-118348 / T) * Tfact;
507  if(473638 / T < 70)
508  RateT[i].BetaHep = 5.54e-17 * pow(T, -0.397) * exp(-473638 / T) * Tfact;
509 
510  /* free-free */
511  RateT[i].Betaff = 1.43e-27 * sqrt(T) * (1.1 + 0.34 * exp(-(5.5 - log10(T)) * (5.5 - log10(T)) / 3));
512 
513  /* recombination */
514 
515  /* Cen 1992 */
516  /* Hydrogen II */
517  RateT[i].AlphaHp = 8.4e-11 * pow(T / 1000, -0.2) / (1. + pow(T / 1.0e6, 0.7)) / sqrt(T);
518  /* Helium II */
519  RateT[i].AlphaHep = 1.5e-10 * pow(T, -0.6353);
520  /* Helium III */
521  RateT[i].AlphaHepp = 4. * RateT[i].AlphaHp;
522  /* Cen 1992 */
523  /* dielectric recombination */
524  if(470000 / T < 70)
525  RateT[i].Alphad = 1.9e-3 * pow(T, -1.5) * exp(-470000 / T) * (1. + 0.3 * exp(-94000 / T));
526 
527  /* collisional ionization */
528  /* Cen 1992 */
529  /* Hydrogen */
530  if(157809.1 / T < 70)
531  RateT[i].GammaeH0 = 5.85e-11 * sqrt(T) * exp(-157809.1 / T) * Tfact;
532  /* Helium */
533  if(285335.4 / T < 70)
534  RateT[i].GammaeHe0 = 2.38e-11 * sqrt(T) * exp(-285335.4 / T) * Tfact;
535  /* Hellium II */
536  if(631515.0 / T < 70)
537  RateT[i].GammaeHep = 5.68e-12 * sqrt(T) * exp(-631515.0 / T) * Tfact;
538  }
539 }
540 
545 void coolsfr::ReadIonizeParams(char *fname)
546 {
547  NheattabUVB = 0;
548  int i, iter;
549  for(iter = 0, i = 0; iter < 2; iter++)
550  {
551  FILE *fdcool;
552  if(!(fdcool = fopen(fname, "r")))
553  Terminate(" Cannot read ionization table in file `%s'\n", fname);
554  if(iter == 0)
555  while(fscanf(fdcool, "%*g %*g %*g %*g %*g %*g %*g") != EOF)
556  NheattabUVB++;
557  if(iter == 1)
558  while(fscanf(fdcool, "%g %g %g %g %g %g %g", &PhotoTUVB[i].variable, &PhotoTUVB[i].gH0, &PhotoTUVB[i].gHe, &PhotoTUVB[i].gHep,
559  &PhotoTUVB[i].eH0, &PhotoTUVB[i].eHe, &PhotoTUVB[i].eHep) != EOF)
560  i++;
561  fclose(fdcool);
562 
563  if(iter == 0)
564  {
565  PhotoTUVB = (photo_table *)Mem.mymalloc("PhotoT", NheattabUVB * sizeof(photo_table));
566  mpi_printf("COOLING: read ionization table with %d entries in file `%s'.\n", NheattabUVB, fname);
567  }
568  }
569  /* ignore zeros at end of treecool file */
570  for(i = 0; i < NheattabUVB; ++i)
571  if(PhotoTUVB[i].gH0 == 0.0)
572  break;
573 
574  NheattabUVB = i;
575  mpi_printf("COOLING: using %d ionization table entries from file `%s'.\n", NheattabUVB, fname);
576 
577  if(NheattabUVB < 1)
578  Terminate("The length of the cooling table has to have at least one entry");
579 }
580 
583 void coolsfr::IonizeParamsUVB(void)
584 {
586  {
587  SetZeroIonization();
588  return;
589  }
590 
591  if(NheattabUVB == 1)
592  {
593  /* treat the one value given as constant with redshift */
594  pc.J_UV = 1;
595  pc.gJH0 = PhotoTUVB[0].gH0;
596  pc.gJHe0 = PhotoTUVB[0].gHe;
597  pc.gJHep = PhotoTUVB[0].gHep;
598  pc.epsH0 = PhotoTUVB[0].eH0;
599  pc.epsHe0 = PhotoTUVB[0].eHe;
600  pc.epsHep = PhotoTUVB[0].eHep;
601  }
602  else
603  {
604  double redshift = 1 / All.Time - 1;
605  double logz = log10(redshift + 1.0);
606  int ilow = 0;
607  for(int i = 0; i < NheattabUVB; i++)
608  {
609  if(PhotoTUVB[i].variable < logz)
610  ilow = i;
611  else
612  break;
613  }
614 
615  if(logz > PhotoTUVB[NheattabUVB - 1].variable || ilow >= NheattabUVB - 1)
616  {
617  SetZeroIonization();
618  }
619  else
620  {
621  double dzlow = logz - PhotoTUVB[ilow].variable;
622  double dzhi = PhotoTUVB[ilow + 1].variable - logz;
623 
624  if(PhotoTUVB[ilow].gH0 == 0 || PhotoTUVB[ilow + 1].gH0 == 0)
625  {
626  SetZeroIonization();
627  }
628  else
629  {
630  pc.J_UV = 1;
631  pc.gJH0 = pow(10., (dzhi * log10(PhotoTUVB[ilow].gH0) + dzlow * log10(PhotoTUVB[ilow + 1].gH0)) / (dzlow + dzhi));
632  pc.gJHe0 = pow(10., (dzhi * log10(PhotoTUVB[ilow].gHe) + dzlow * log10(PhotoTUVB[ilow + 1].gHe)) / (dzlow + dzhi));
633  pc.gJHep = pow(10., (dzhi * log10(PhotoTUVB[ilow].gHep) + dzlow * log10(PhotoTUVB[ilow + 1].gHep)) / (dzlow + dzhi));
634  pc.epsH0 = pow(10., (dzhi * log10(PhotoTUVB[ilow].eH0) + dzlow * log10(PhotoTUVB[ilow + 1].eH0)) / (dzlow + dzhi));
635  pc.epsHe0 = pow(10., (dzhi * log10(PhotoTUVB[ilow].eHe) + dzlow * log10(PhotoTUVB[ilow + 1].eHe)) / (dzlow + dzhi));
636  pc.epsHep = pow(10., (dzhi * log10(PhotoTUVB[ilow].eHep) + dzlow * log10(PhotoTUVB[ilow + 1].eHep)) / (dzlow + dzhi));
637  }
638  }
639  }
640 }
641 
644 void coolsfr::SetZeroIonization(void) { memset(&pc, 0, sizeof(photo_current)); }
645 
648 void coolsfr::IonizeParams(void) { IonizeParamsUVB(); }
649 
656 void coolsfr::InitCool(void)
657 {
658  /* set default hydrogen mass fraction */
659  GasState.XH = HYDROGEN_MASSFRAC;
660 
661  /* zero photo-ionization/heating rates */
662  SetZeroIonization();
663 
664  /* allocate and construct rate table */
665  RateT = (rate_table *)Mem.mymalloc("RateT", (NCOOLTAB + 1) * sizeof(rate_table));
666  ;
667  MakeRateTable();
668 
669  /* read photo tables */
670  ReadIonizeParams(All.TreecoolFile);
671 
672  All.Time = All.TimeBegin;
674 
675  IonizeParams();
676 }
677 
681 void coolsfr::cooling_only(simparticles *Sp) /* normal cooling routine when star formation is disabled */
682 {
683  TIMER_START(CPU_COOLING_SFR);
685 
686  gas_state gs = GasState;
687  do_cool_data DoCool = DoCoolData;
688 
689  for(int i = 0; i < Sp->TimeBinsHydro.NActiveParticles; i++)
690  {
691  int target = Sp->TimeBinsHydro.ActiveParticleList[i];
692  if(Sp->P[target].getType() == 0)
693  {
694  if(Sp->P[target].getMass() == 0 && Sp->P[target].ID.get() == 0)
695  continue; /* skip particles that have been swallowed or eliminated */
696 
697  cool_sph_particle(Sp, target, &gs, &DoCool);
698  }
699  }
700  TIMER_STOP(CPU_COOLING_SFR);
701 }
702 
711 void coolsfr::cool_sph_particle(simparticles *Sp, int i, gas_state *gs, do_cool_data *DoCool)
712 {
713  double dens = Sp->SphP[i].Density;
714 
715  double dt = (Sp->P[i].getTimeBinHydro() ? (((integertime)1) << Sp->P[i].getTimeBinHydro()) : 0) * All.Timebase_interval;
716 
717  double dtime = All.cf_atime * dt / All.cf_atime_hubble_a;
718 
719  double utherm = Sp->get_utherm_from_entropy(i);
720 
721  double ne = Sp->SphP[i].Ne; /* electron abundance (gives ionization state and mean molecular weight) */
722  double unew = DoCooling(std::max<double>(All.MinEgySpec, utherm), dens * All.cf_a3inv, dtime, &ne, gs, DoCool);
723  Sp->SphP[i].Ne = ne;
724 
725  if(unew < 0)
726  Terminate("invalid temperature: i=%d unew=%g\n", i, unew);
727 
728  double du = unew - utherm;
729 
730  if(unew < All.MinEgySpec)
731  du = All.MinEgySpec - utherm;
732 
733  utherm += du;
734 
735 #ifdef OUTPUT_COOLHEAT
736  if(dtime > 0)
737  Sp->SphP[i].CoolHeat = du * Sp->P[i].getMass() / dtime;
738 #endif
739 
740  Sp->set_entropy_from_utherm(utherm, i);
742 }
743 
744 #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
#define HYDROGEN_MASSFRAC
Definition: constants.h:112
#define GAMMA_MINUS1
Definition: constants.h:110
#define PROTONMASS
Definition: constants.h:124
#define SMALLNUM
Definition: constants.h:83
#define MAXITER
Definition: constants.h:305
#define BOLTZMANN
Definition: constants.h:120
int integertime
Definition: constants.h:331
#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
memory Mem
Definition: main.cc:44
expr exp(half arg)
Definition: half.hpp:2724
half fabs(half arg)
Definition: half.hpp:2627
expr log10(half arg)
Definition: half.hpp:2752
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)