double INLINE_FUNC hubble_function(double a) { double hubble_a; #ifdef EXTERNALHUBBLE hubble_a = hubble_function_external(a); #else hubble_a = All.Omega0 / (a * a * a) + (1 - All.Omega0 - All.OmegaLambda) / (a * a) #ifdef DARKENERGY + DarkEnergy_a(a); #else + All.OmegaLambda; #endif hubble_a = All.Hubble * sqrt(hubble_a); #endif #ifdef TIMEDEPGRAV hubble_a *= dHfak(a); #endif return (hubble_a); }
/* This function determines the present equation of state * in case the HPM method is selected. It normally does this * by following the ionization history of two fiducial gas * elements at the mean density, and at 1.1 times the mean density. * From the pressure values, and effective EQS index is then computed. */ void hpm_find_eqs_selfconsistent(void) { double dt, dtime, hubble_a = 0, a3inv; double time_hubble_a, dmax1, dmax2; double u0, u1, meanWeight, temp; if(All.ComovingIntegrationOn) { /* Factors for comoving integration of hydro */ a3inv = 1 / (All.Time * All.Time * All.Time); #if defined(DEDM_HUBBLE) || defined(VDE) hubble_a = getH_a(All.Time); #else hubble_a = All.Omega0 / (All.Time * All.Time * All.Time) + (1 - All.Omega0 - All.OmegaLambda) / (All.Time * All.Time) #ifdef DARKENERGY + DarkEnergy_a(All.Time); #else + All.OmegaLambda; #endif #endif hubble_a = All.Hubble * sqrt(hubble_a); time_hubble_a = All.Time * hubble_a; } else a3inv = time_hubble_a = 1; dt = (P[0].Ti_endstep - P[0].Ti_begstep) * All.Timebase_interval; /* the time-step */ if(All.ComovingIntegrationOn) dtime = All.Time * dt / time_hubble_a; else dtime = dt; All.HPM_rho0 = All.OmegaBaryon * 3 * All.Hubble * All.Hubble / (8 * M_PI * All.G); All.HPM_rho1 = 1.1 * All.HPM_rho0; u0 = All.HPM_entr0 / GAMMA_MINUS1 * pow(All.HPM_rho0 * a3inv, GAMMA_MINUS1); u1 = All.HPM_entr1 / GAMMA_MINUS1 * pow(All.HPM_rho1 * a3inv, GAMMA_MINUS1); u0 = DoCooling(DMAX(All.MinEgySpec, u0), All.HPM_rho0 * a3inv, dtime, &All.HPM_ne0); u1 = DoCooling(DMAX(All.MinEgySpec, u1), All.HPM_rho1 * a3inv, dtime, &All.HPM_ne1); All.HPM_entr0 = u0 * GAMMA_MINUS1 / pow(All.HPM_rho0 * a3inv, GAMMA_MINUS1); All.HPM_entr1 = u1 * GAMMA_MINUS1 / pow(All.HPM_rho1 * a3inv, GAMMA_MINUS1); /* Note: All.HPM_P0 is a "comoving pressure", i.e. computed in gadget's * internal unit convention using the comoving density and the physical entropy. */ All.HPM_P0 = All.HPM_entr0 * pow(All.HPM_rho0, GAMMA); All.HPM_P1 = All.HPM_entr1 * pow(All.HPM_rho1, GAMMA); All.HPM_alpha = log(All.HPM_P1 / All.HPM_P0) / log(All.HPM_rho1 / All.HPM_rho0); /* compute the temperature corresponding to P0, for information purposes only */ meanWeight = 4.0 / (3 * HYDROGEN_MASSFRAC + 1 + 4 * HYDROGEN_MASSFRAC * All.HPM_ne0) * PROTONMASS; temp = meanWeight / BOLTZMANN * GAMMA_MINUS1 * u0 * All.UnitEnergy_in_cgs / All.UnitMass_in_g; if(ThisTask == 0) { printf("IGM-EQS: a=%g T0=%g P0=%g alpha=%g\n", All.Time, temp, All.HPM_P0, All.HPM_alpha); fflush(stdout); } }
void blackhole_accretion(void) { int i, j, k, n, ngrp, maxfill, source; int ntot, ntotleft, ndone, ndonetot; int *nbuffer, *noffset, *nsend_local, *nsend; int level, sendTask, recvTask; int nexport, place, num_blackholes; int Ntot_gas_swallowed, Ntot_BH_swallowed; int total_num_blackholes; double mdot, rho, bhvel, soundspeed, meddington, dt, mdot_in_msun_per_year; double mass_real, total_mass_real; double mass_holes, total_mass_holes, total_mdot; double fac; double mdoteddington, total_mdoteddington; MPI_Status status; #ifdef PERIODIC boxSize = All.BoxSize; boxHalf = 0.5 * All.BoxSize; #ifdef LONG_X boxHalf_X = boxHalf * LONG_X; boxSize_X = boxSize * LONG_X; #endif #ifdef LONG_Y boxHalf_Y = boxHalf * LONG_Y; boxSize_Y = boxSize * LONG_Y; #endif #ifdef LONG_Z boxHalf_Z = boxHalf * LONG_Z; boxSize_Z = boxSize * LONG_Z; #endif #endif if(ThisTask == 0) { printf("Beginning black-hole accretion\n"); fflush(stdout); } if(All.ComovingIntegrationOn) { ascale = All.Time; #if defined(DEDM_HUBBLE) || defined(VDE) hubble_a = getH_a(All.Time); #else hubble_a = All.Omega0 / (All.Time * All.Time * All.Time) + (1 - All.Omega0 - All.OmegaLambda) / (All.Time * All.Time) #ifdef DARKENERGY + DarkEnergy_a(All.Time); #else + All.OmegaLambda; #endif #endif hubble_a = All.Hubble * sqrt(hubble_a); } else hubble_a = ascale = 1; /* Let's first compute the Mdot values */ for(n = 0; n < NumPart; n++) if(P[n].Type == 5) if(P[n].Ti_endstep == All.Ti_Current) { mdot = 0; /* if no accretion model is enabled, we have mdot=0 */ rho = P[n].BH_Density; bhvel = sqrt(pow(P[n].Vel[0] - P[n].BH_SurroundingGasVel[0], 2) + pow(P[n].Vel[1] - P[n].BH_SurroundingGasVel[1], 2) + pow(P[n].Vel[2] - P[n].BH_SurroundingGasVel[2], 2)); if(All.ComovingIntegrationOn) { bhvel /= All.Time; rho /= pow(All.Time, 3); } soundspeed = sqrt(GAMMA * P[n].BH_Entropy * pow(rho, GAMMA_MINUS1)); /* Note: we take here a radiative efficiency of 0.1 for Eddington accretion */ meddington = (4 * M_PI * GRAVITY * C * PROTONMASS / (0.1 * C * C * THOMPSON)) * P[n].BH_Mass * All.UnitTime_in_s; #ifdef BONDI mdot = 4. * M_PI * All.BlackHoleAccretionFactor * All.G * All.G * P[n].BH_Mass * P[n].BH_Mass * rho / pow((pow(soundspeed, 2) + pow(bhvel, 2)), 1.5); #endif #ifdef HIGHDENS_ACCRETION if(rho >= 3 * All.PhysDensThresh) mdot = meddington; else mdot = 0; #endif #ifdef MODIFIEDBONDI mdot = All.BlackHoleAccretionFactor * meddington * (rho / All.BlackHoleRefDensity) / pow((pow(soundspeed, 2) + pow(bhvel, 2)), 1.5) * pow(All.BlackHoleRefSoundspeed, 3); #endif #ifdef ENFORCE_EDDINGTON_LIMIT if(mdot > All.BlackHoleEddingtonFactor * meddington) mdot = All.BlackHoleEddingtonFactor * meddington; #endif P[n].BH_Mdot = mdot; if(meddington > 0) P[n].BH_MdotEddington = mdot / meddington; /* this stores the accretion rate in units of the Eddington rate */ else P[n].BH_MdotEddington = 0; if(P[n].BH_Mass > 0) fprintf(FdBlackHolesDetails, "BH=%d %g %g %g %g %g\n", P[n].ID, All.Time, P[n].BH_Mass, mdot, rho, soundspeed); dt = (P[n].Ti_endstep - P[n].Ti_begstep) * All.Timebase_interval / hubble_a; #ifdef BH_DRAG /* add a drag force for the black-holes, accounting for the accretion */ if(P[n].BH_Mass > 0) { /* fac = P[n].BH_Mdot * dt / P[n].BH_Mass; */ fac = meddington * dt / P[n].BH_Mass; if(fac > 1) fac = 1; if(dt > 0) for(k = 0; k < 3; k++) P[n].GravAccel[k] += -ascale * ascale * fac / dt * (P[n].Vel[k] - P[n].BH_SurroundingGasVel[k]) / ascale; } #endif P[n].BH_Mass += P[n].BH_Mdot * dt; #ifdef BH_KINETICFEEDBACK if(mdot >= 0.99 * meddington) { P[n].ActiveTime += dt; P[n].ActiveEnergy += All.BlackHoleFeedbackFactor * 0.1 * mdot * dt * pow(C / All.UnitVelocity_in_cm_per_s, 2); } #endif } /* Now let's invoke the functions that stochasticall swallow gas * and deal with black hole mergers. */ if(ThisTask == 0) { printf("Start swallowing of gas particles and black holes\n"); fflush(stdout); } N_gas_swallowed = N_BH_swallowed = 0; for(n = 0, num_blackholes = 0; n < NumPart; n++) { if(P[n].Type == 5) if(P[n].Ti_endstep == All.Ti_Current) num_blackholes++; } MPI_Allreduce(&num_blackholes, &ntot, 1, MPI_INT, MPI_SUM, MPI_COMM_WORLD); noffset = malloc(sizeof(int) * NTask); /* offsets of bunches in common list */ nbuffer = malloc(sizeof(int) * NTask); nsend_local = malloc(sizeof(int) * NTask); nsend = malloc(sizeof(int) * NTask * NTask); /* to avoid race conditions in BH mergers, we have to do it twice */ /* for(mode = 1; mode <= 2; mode++) */ { i = 0; /* first particle for this task */ ntotleft = ntot; /* particles left for all tasks together */ while(ntotleft > 0) { for(j = 0; j < NTask; j++) nsend_local[j] = 0; /* do local particles and prepare export list */ for(nexport = 0, ndone = 0; i < NumPart && nexport < All.BunchSizeBlackhole - NTask; i++) if(P[i].Type == 5) if(P[i].Ti_endstep == All.Ti_Current) { ndone++; for(j = 0; j < NTask; j++) Exportflag[j] = 0; blackhole_evaluate(i, 0); for(j = 0; j < NTask; j++) { if(Exportflag[j]) { for(k = 0; k < 3; k++) BlackholeDataIn[nexport].Pos[k] = P[i].Pos[k]; BlackholeDataIn[nexport].Hsml = PPP[i].Hsml; BlackholeDataIn[nexport].Mass = P[i].Mass; BlackholeDataIn[nexport].BH_Mass = P[i].BH_Mass; BlackholeDataIn[nexport].Vel[k] = P[i].Vel[k]; #ifdef BH_KINETICFEEDBACK BlackholeDataIn[nexport].ActiveTime = P[i].ActiveTime; BlackholeDataIn[nexport].ActiveEnergy = P[i].ActiveEnergy; #endif BlackholeDataIn[nexport].Density = P[i].BH_Density; BlackholeDataIn[nexport].Mdot = P[i].BH_Mdot; BlackholeDataIn[nexport].Dt = (P[i].Ti_endstep - P[i].Ti_begstep) * All.Timebase_interval / hubble_a; BlackholeDataIn[nexport].ID = P[i].ID; BlackholeDataIn[nexport].Index = i; BlackholeDataIn[nexport].Task = j; nexport++; nsend_local[j]++; } } } qsort(BlackholeDataIn, nexport, sizeof(struct blackholedata_in), blackhole_compare_key); for(j = 1, noffset[0] = 0; j < NTask; j++) noffset[j] = noffset[j - 1] + nsend_local[j - 1]; MPI_Allgather(nsend_local, NTask, MPI_INT, nsend, NTask, MPI_INT, MPI_COMM_WORLD); /* now do the particles that need to be exported */ for(level = 1; level < (1 << PTask); level++) { for(j = 0; j < NTask; j++) nbuffer[j] = 0; for(ngrp = level; ngrp < (1 << PTask); ngrp++) { maxfill = 0; for(j = 0; j < NTask; j++) { if((j ^ ngrp) < NTask) if(maxfill < nbuffer[j] + nsend[(j ^ ngrp) * NTask + j]) maxfill = nbuffer[j] + nsend[(j ^ ngrp) * NTask + j]; } if(maxfill >= All.BunchSizeBlackhole) break; sendTask = ThisTask; recvTask = ThisTask ^ ngrp; if(recvTask < NTask) { if(nsend[ThisTask * NTask + recvTask] > 0 || nsend[recvTask * NTask + ThisTask] > 0) { /* get the particles */ MPI_Sendrecv(&BlackholeDataIn[noffset[recvTask]], nsend_local[recvTask] * sizeof(struct blackholedata_in), MPI_BYTE, recvTask, TAG_BH_A, &BlackholeDataGet[nbuffer[ThisTask]], nsend[recvTask * NTask + ThisTask] * sizeof(struct blackholedata_in), MPI_BYTE, recvTask, TAG_BH_A, MPI_COMM_WORLD, &status); } } for(j = 0; j < NTask; j++) if((j ^ ngrp) < NTask) nbuffer[j] += nsend[(j ^ ngrp) * NTask + j]; } /* now do the imported particles */ for(j = 0; j < nbuffer[ThisTask]; j++) blackhole_evaluate(j, 1); /* get the result */ for(j = 0; j < NTask; j++) nbuffer[j] = 0; for(ngrp = level; ngrp < (1 << PTask); ngrp++) { maxfill = 0; for(j = 0; j < NTask; j++) { if((j ^ ngrp) < NTask) if(maxfill < nbuffer[j] + nsend[(j ^ ngrp) * NTask + j]) maxfill = nbuffer[j] + nsend[(j ^ ngrp) * NTask + j]; } if(maxfill >= All.BunchSizeBlackhole) break; sendTask = ThisTask; recvTask = ThisTask ^ ngrp; if(recvTask < NTask) { if(nsend[ThisTask * NTask + recvTask] > 0 || nsend[recvTask * NTask + ThisTask] > 0) { /* send the results */ MPI_Sendrecv(&BlackholeDataResult[nbuffer[ThisTask]], nsend[recvTask * NTask + ThisTask] * sizeof(struct blackholedata_out), MPI_BYTE, recvTask, TAG_BH_B, &BlackholeDataPartialResult[noffset[recvTask]], nsend_local[recvTask] * sizeof(struct blackholedata_out), MPI_BYTE, recvTask, TAG_BH_B, MPI_COMM_WORLD, &status); /* add the result to the particles */ for(j = 0; j < nsend_local[recvTask]; j++) { source = j + noffset[recvTask]; place = BlackholeDataIn[source].Index; if(P[place].Mass + BlackholeDataPartialResult[source].Mass > 0) { for(k = 0; k < 3; k++) P[place].Vel[k] = (P[place].Vel[k] * P[place].Mass + BlackholeDataPartialResult[source].AccretedMomentum[k]) / (P[place].Mass + BlackholeDataPartialResult[source].Mass); } P[place].Mass += BlackholeDataPartialResult[source].Mass; P[place].BH_Mass += BlackholeDataPartialResult[source].BH_Mass; #ifdef REPOSITION_ON_POTMIN if(P[place].BH_MinPot > BlackholeDataPartialResult[source].BH_MinPot) { P[place].BH_MinPot = BlackholeDataPartialResult[source].BH_MinPot; for(k = 0; k < 3; k++) P[place].BH_MinPotPos[k] = BlackholeDataPartialResult[source].BH_MinPotPos[k]; } #endif } } } for(j = 0; j < NTask; j++) if((j ^ ngrp) < NTask) nbuffer[j] += nsend[(j ^ ngrp) * NTask + j]; } level = ngrp - 1; } MPI_Allreduce(&ndone, &ndonetot, 1, MPI_INT, MPI_SUM, MPI_COMM_WORLD); ntotleft -= ndonetot; } } free(nsend); free(nsend_local); free(nbuffer); free(noffset); MPI_Reduce(&N_gas_swallowed, &Ntot_gas_swallowed, 1, MPI_INT, MPI_SUM, 0, MPI_COMM_WORLD); MPI_Reduce(&N_BH_swallowed, &Ntot_BH_swallowed, 1, MPI_INT, MPI_SUM, 0, MPI_COMM_WORLD); if(ThisTask == 0) { printf("Accretion done: %d gas particles swallowed, %d BH particles swallowed\n", Ntot_gas_swallowed, Ntot_BH_swallowed); fflush(stdout); } #ifdef REPOSITION_ON_POTMIN for(n = 0; n < NumPart; n++) if(P[n].Ti_endstep == All.Ti_Current) if(P[n].Type == 5) for(k = 0; k < 3; k++) P[n].Pos[k] = P[n].BH_MinPotPos[k]; #endif for(n = 0, num_blackholes = 0, mdot = 0, mass_holes = 0, mass_real = 0, mdoteddington = 0; n < NumPart; n++) if(P[n].Type == 5) { #ifdef BH_KINETICFEEDBACK if(P[n].ActiveTime > All.BlackHoleActiveTime) { P[n].ActiveTime = 0; P[n].ActiveEnergy = 0; } #endif mass_holes += P[n].BH_Mass; mass_real += P[n].Mass; mdot += P[n].BH_Mdot; mdoteddington += P[n].BH_MdotEddington; num_blackholes++; } MPI_Reduce(&mass_holes, &total_mass_holes, 1, MPI_DOUBLE, MPI_SUM, 0, MPI_COMM_WORLD); MPI_Reduce(&mass_real, &total_mass_real, 1, MPI_DOUBLE, MPI_SUM, 0, MPI_COMM_WORLD); MPI_Reduce(&mdot, &total_mdot, 1, MPI_DOUBLE, MPI_SUM, 0, MPI_COMM_WORLD); MPI_Reduce(&mdoteddington, &total_mdoteddington, 1, MPI_DOUBLE, MPI_SUM, 0, MPI_COMM_WORLD); MPI_Reduce(&num_blackholes, &total_num_blackholes, 1, MPI_INT, MPI_SUM, 0, MPI_COMM_WORLD); if(ThisTask == 0) { /* convert to solar masses per yr */ mdot_in_msun_per_year = total_mdot * (All.UnitMass_in_g / SOLAR_MASS) / (All.UnitTime_in_s / SEC_PER_YEAR); fprintf(FdBlackHoles, "%g %d %g %g %g %g %g\n", All.Time, total_num_blackholes, total_mass_holes, total_mdot, mdot_in_msun_per_year, total_mass_real, total_mdoteddington); fflush(FdBlackHoles); } fflush(FdBlackHolesDetails); }