12#include "gadgetconfig.h"
16#include <gsl/gsl_math.h>
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"
46double coolsfr::DoCooling(
double u_old,
double rho,
double dt,
double *ne_guess, gas_state *gs, do_cool_data *DoCool)
48 DoCool->u_old_input = u_old;
49 DoCool->rho_input = rho;
50 DoCool->dt_input = dt;
51 DoCool->ne_guess_input = *ne_guess;
53 if(!gsl_finite(u_old))
54 Terminate(
"invalid input: u_old=%g\n", u_old);
56 if(u_old < 0 || rho < 0)
64 double ratefact = gs->nHcgs * gs->nHcgs / rho;
70 double LambdaNet = CoolingRateFromU(u, rho, ne_guess, gs, DoCool);
74 if(u - u_old - ratefact * LambdaNet * dt < 0)
78 while(u_upper - u_old - ratefact * CoolingRateFromU(u_upper, rho, ne_guess, gs, DoCool) * dt < 0)
85 if(u - u_old - ratefact * LambdaNet * dt > 0)
89 while(u_lower - u_old - ratefact * CoolingRateFromU(u_lower, rho, ne_guess, gs, DoCool) * dt > 0)
100 u = 0.5 * (u_lower + u_upper);
102 LambdaNet = CoolingRateFromU(u, rho, ne_guess, gs, DoCool);
104 if(u - u_old - ratefact * LambdaNet * dt > 0)
113 du = u_upper - u_lower;
118 printf(
"u= %g\n", u);
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);
141double coolsfr::GetCoolingTime(
double u_old,
double rho,
double *ne_guess, gas_state *gs, do_cool_data *DoCool)
143 DoCool->u_old_input = u_old;
144 DoCool->rho_input = rho;
145 DoCool->ne_guess_input = *ne_guess;
151 double ratefact = gs->nHcgs * gs->nHcgs / rho;
155 double LambdaNet = CoolingRateFromU(u, rho, ne_guess, gs, DoCool);
162 double coolingtime = u_old / (-ratefact * LambdaNet);
180double coolsfr::convert_u_to_temp(
double u,
double rho,
double *ne_guess, gas_state *gs,
const do_cool_data *DoCool)
183 double rho_input = rho;
184 double ne_input = *ne_guess;
186 double mu = (1 + 4 * gs->yhelium) / (1 + gs->yhelium + *ne_guess);
194 double ne_old = *ne_guess;
196 find_abundances_and_rates(
log10(temp), rho, ne_guess, gs, DoCool);
199 mu = (1 + 4 * gs->yhelium) / (1 + gs->yhelium + *ne_guess);
203 max = std::max<double>(max, temp_new / (1 + gs->yhelium + *ne_guess) *
fabs((*ne_guess - ne_old) / (temp_new - temp_old + 1.0)));
205 temp = temp_old + (temp_new - temp_old) / (1 + max);
209 printf(
"-> temp= %g ne=%g\n", temp, *ne_guess);
211 while(
fabs(temp - temp_old) > 1.0e-3 * temp && iter <
MAXITER);
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);
233void coolsfr::find_abundances_and_rates(
double logT,
double rho,
double *ne_guess, gas_state *gs,
const do_cool_data *DoCool)
235 double logT_input = logT;
236 double rho_input = rho;
237 double ne_input = *ne_guess;
239 if(!gsl_finite(logT))
245 gs->nHe0 = gs->yhelium;
260 gs->nHepp = gs->yhelium;
261 gs->ne = gs->nHp + 2.0 * gs->nHepp;
266 double t = (logT - Tmin) / deltaT;
269 double flow = 1 - fhi;
277 double neold = gs->ne;
279 gs->necgs = gs->ne * gs->nHcgs;
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;
294 if(gs->necgs <= 1.e-25 || pc.J_UV == 0)
296 gs->gJH0ne = gs->gJHe0ne = gs->gJHepne = 0;
300 gs->gJH0ne = pc.gJH0 / gs->necgs;
301 gs->gJHe0ne = pc.gJHe0 / gs->necgs;
302 gs->gJHepne = pc.gJHep / gs->necgs;
305 gs->nH0 = gs->aHp / (gs->aHp + gs->geH0 + gs->gJH0ne);
306 gs->nHp = 1.0 - gs->nH0;
308 if((gs->gJHe0ne + gs->geHe0) <=
SMALLNUM)
312 gs->nHe0 = gs->yhelium;
316 gs->nHep = gs->yhelium /
317 (1.0 + (gs->aHep + gs->ad) / (gs->geHe0 + gs->gJHe0ne) + (gs->geHep + gs->gJHepne) / gs->aHepp);
318 gs->nHe0 = gs->nHep * (gs->aHep + gs->ad) / (gs->geHe0 + gs->gJHe0ne);
319 gs->nHepp = gs->nHep * (gs->geHep + gs->gJHepne) / gs->aHepp;
324 gs->ne = gs->nHp + gs->nHep + 2 * gs->nHepp;
325 gs->necgs = gs->ne * gs->nHcgs;
330 double nenew = 0.5 * (gs->ne + neold);
332 gs->necgs = gs->ne * gs->nHcgs;
334 if(
fabs(gs->ne - neold) < 1.0e-4)
338 printf(
"ne= %g niter=%d\n", gs->ne, niter);
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);
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;
365double coolsfr::CoolingRateFromU(
double u,
double rho,
double *ne_guess, gas_state *gs,
const do_cool_data *DoCool)
367 double temp = convert_u_to_temp(u, rho, ne_guess, gs, DoCool);
369 return CoolingRate(
log10(temp), rho, ne_guess, gs, DoCool);
382double coolsfr::AbundanceRatios(
double u,
double rho,
double *ne_guess,
double *nH0_pointer,
double *nHeII_pointer)
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;
393 double temp = convert_u_to_temp(u, rho, ne_guess, &gs, &DoCool);
395 *nH0_pointer = gs.nH0;
396 *nHeII_pointer = gs.nHep;
408double coolsfr::CoolingRate(
double logT,
double rho,
double *nelec, gas_state *gs,
const do_cool_data *DoCool)
413 logT = Tmin + 0.5 * deltaT;
419 find_abundances_and_rates(logT, rho, nelec, gs, DoCool);
422 double T =
pow(10.0, logT);
424 double LambdaExcH0 = gs->bH0 * gs->ne * gs->nH0;
425 double LambdaExcHep = gs->bHep * gs->ne * gs->nHep;
426 double LambdaExc = LambdaExcH0 + LambdaExcHep;
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;
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;
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;
444 Lambda += LambdaCmptn;
449 Heat += (gs->nH0 * pc.epsH0 + gs->nHe0 * pc.epsHe0 + gs->nHep * pc.epsHep) / gs->nHcgs;
460 gs->nHepp = gs->yhelium;
461 gs->ne = gs->nHp + 2.0 * gs->nHepp;
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;
469 double redshift = 1 /
All.
Time - 1;
471 LambdaCmptn = 5.65e-36 * gs->ne * (T - 2.73 * (1. + redshift)) *
pow(1. + redshift, 4.) / gs->nHcgs;
476 Lambda = LambdaFF + LambdaCmptn;
479 return (Heat - Lambda);
486void coolsfr::MakeRateTable(
void)
488 GasState.yhelium = (1 - GasState.XH) / (4 * GasState.XH);
491 deltaT = (Tmax - Tmin) / NCOOLTAB;
492 GasState.ethmin =
pow(10.0, Tmin) * (1. + GasState.yhelium) / ((1. + 4. * GasState.yhelium) * GasState.mhboltz *
GAMMA_MINUS1);
495 for(
int i = 0; i <= NCOOLTAB; i++)
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;
500 double T =
pow(10.0, Tmin + deltaT * i);
501 double Tfact = 1.0 / (1 +
sqrt(T / 1.0e5));
506 RateT[i].BetaH0 = 7.5e-19 *
exp(-118348 / T) * Tfact;
508 RateT[i].BetaHep = 5.54e-17 *
pow(T, -0.397) *
exp(-473638 / T) * Tfact;
511 RateT[i].Betaff = 1.43e-27 *
sqrt(T) * (1.1 + 0.34 *
exp(-(5.5 -
log10(T)) * (5.5 -
log10(T)) / 3));
517 RateT[i].AlphaHp = 8.4e-11 *
pow(T / 1000, -0.2) / (1. +
pow(T / 1.0e6, 0.7)) /
sqrt(T);
519 RateT[i].AlphaHep = 1.5e-10 *
pow(T, -0.6353);
521 RateT[i].AlphaHepp = 4. * RateT[i].AlphaHp;
525 RateT[i].Alphad = 1.9e-3 *
pow(T, -1.5) *
exp(-470000 / T) * (1. + 0.3 *
exp(-94000 / T));
530 if(157809.1 / T < 70)
531 RateT[i].GammaeH0 = 5.85e-11 *
sqrt(T) *
exp(-157809.1 / T) * Tfact;
533 if(285335.4 / T < 70)
534 RateT[i].GammaeHe0 = 2.38e-11 *
sqrt(T) *
exp(-285335.4 / T) * Tfact;
536 if(631515.0 / T < 70)
537 RateT[i].GammaeHep = 5.68e-12 *
sqrt(T) *
exp(-631515.0 / T) * Tfact;
545void coolsfr::ReadIonizeParams(
char *fname)
549 for(iter = 0, i = 0; iter < 2; iter++)
552 if(!(fdcool = fopen(fname,
"r")))
553 Terminate(
" Cannot read ionization table in file `%s'\n", fname);
555 while(fscanf(fdcool,
"%*g %*g %*g %*g %*g %*g %*g") != EOF)
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)
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);
570 for(i = 0; i < NheattabUVB; ++i)
571 if(PhotoTUVB[i].gH0 == 0.0)
575 mpi_printf(
"COOLING: using %d ionization table entries from file `%s'.\n", NheattabUVB, fname);
578 Terminate(
"The length of the cooling table has to have at least one entry");
583void coolsfr::IonizeParamsUVB(
void)
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;
604 double redshift = 1 /
All.
Time - 1;
605 double logz =
log10(redshift + 1.0);
607 for(
int i = 0; i < NheattabUVB; i++)
609 if(PhotoTUVB[i].variable < logz)
615 if(logz > PhotoTUVB[NheattabUVB - 1].variable || ilow >= NheattabUVB - 1)
621 double dzlow = logz - PhotoTUVB[ilow].variable;
622 double dzhi = PhotoTUVB[ilow + 1].variable - logz;
624 if(PhotoTUVB[ilow].gH0 == 0 || PhotoTUVB[ilow + 1].gH0 == 0)
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));
644void coolsfr::SetZeroIonization(
void) { memset(&pc, 0,
sizeof(photo_current)); }
648void coolsfr::IonizeParams(
void) { IonizeParamsUVB(); }
656void coolsfr::InitCool(
void)
665 RateT = (rate_table *)
Mem.mymalloc(
"RateT", (NCOOLTAB + 1) *
sizeof(rate_table));
670 ReadIonizeParams(
All.TreecoolFile);
686 gas_state gs = GasState;
687 do_cool_data DoCool = DoCoolData;
697 cool_sph_particle(Sp, target, &gs, &DoCool);
711void coolsfr::cool_sph_particle(
simparticles *Sp,
int i, gas_state *gs, do_cool_data *DoCool)
721 double ne = Sp->
SphP[i].Ne;
726 Terminate(
"invalid temperature: i=%d unew=%g\n", i, unew);
728 double du = unew - utherm;
735#ifdef OUTPUT_COOLHEAT
global_data_all_processes All
double get_utherm_from_entropy(int i)
void set_entropy_from_utherm(double utherm, int i)
TimeBinData TimeBinsHydro
#define HYDROGEN_MASSFRAC
#define TIMER_START(counter)
Starts the timer counter.
#define TIMER_STOP(counter)
Stops the timer counter.
expr pow(half base, half exp)
double UnitPressure_in_cgs
double UnitDensity_in_cgs
int ComovingIntegrationOn
void set_cosmo_factors_for_current_time(void)
unsigned char getTimeBinHydro(void)
unsigned char getType(void)
void set_thermodynamic_variables(void)