MP-Gadget  5.0.1.dev1-76bc7d4726-dirty
Enumerations | Functions
blackhole.h File Reference
#include "utils/paramset.h"
#include "timestep.h"
#include "forcetree.h"
#include "physconst.h"
Include dependency graph for blackhole.h:
This graph shows which files directly or indirectly include this file:

Go to the source code of this file.

Enumerations

enum  BlackHoleFeedbackMethod {
  BH_FEEDBACK_TOPHAT = 0x2 , BH_FEEDBACK_SPLINE = 0x4 , BH_FEEDBACK_MASS = 0x8 , BH_FEEDBACK_VOLUME = 0x10 ,
  BH_FEEDBACK_OPTTHIN = 0x20
}
 

Functions

void set_blackhole_params (ParameterSet *ps)
 
void blackhole (const ActiveParticles *act, double atime, Cosmology *CP, ForceTree *tree, const struct UnitSystem units, FILE *FdBlackHoles, FILE *FdBlackholeDetails)
 
void blackhole_make_one (int index, const double atime)
 
int BHGetRepositionEnabled (void)
 

Enumeration Type Documentation

◆ BlackHoleFeedbackMethod

Enumerator
BH_FEEDBACK_TOPHAT 
BH_FEEDBACK_SPLINE 
BH_FEEDBACK_MASS 
BH_FEEDBACK_VOLUME 
BH_FEEDBACK_OPTTHIN 

Definition at line 8 of file blackhole.h.

8  {
9  BH_FEEDBACK_TOPHAT = 0x2,
10  BH_FEEDBACK_SPLINE = 0x4,
11  BH_FEEDBACK_MASS = 0x8,
12  BH_FEEDBACK_VOLUME = 0x10,
13  BH_FEEDBACK_OPTTHIN = 0x20,
14 };
@ BH_FEEDBACK_OPTTHIN
Definition: blackhole.h:13
@ BH_FEEDBACK_MASS
Definition: blackhole.h:11
@ BH_FEEDBACK_SPLINE
Definition: blackhole.h:10
@ BH_FEEDBACK_VOLUME
Definition: blackhole.h:12
@ BH_FEEDBACK_TOPHAT
Definition: blackhole.h:9

Function Documentation

◆ BHGetRepositionEnabled()

int BHGetRepositionEnabled ( void  )

Definition at line 60 of file blackhole.c.

61 {
63 }
struct BlackholeParams blackhole_params
int BlackHoleRepositionEnabled
Definition: blackhole.c:32

Referenced by real_drift_particle().

Here is the caller graph for this function:

◆ blackhole()

void blackhole ( const ActiveParticles act,
double  atime,
Cosmology CP,
ForceTree tree,
const struct UnitSystem  units,
FILE *  FdBlackHoles,
FILE *  FdBlackholeDetails 
)

Definition at line 516 of file blackhole.c.

517 {
518  /* Do nothing if no black holes*/
519  int64_t totbh;
520  MPI_Allreduce(&SlotsManager->info[5].size, &totbh, 1, MPI_INT64, MPI_SUM, MPI_COMM_WORLD);
521  if(totbh == 0)
522  return;
523  int i;
524 
525  walltime_measure("/Misc");
526  struct BHPriv priv[1] = {0};
527  priv->units = units;
528 
529  /*************************************************************************/
530  TreeWalk tw_dynfric[1] = {{0}};
531  tw_dynfric->ev_label = "BH_DYNFRIC";
533  tw_dynfric->ngbiter_type_elsize = sizeof(TreeWalkNgbIterBHDynfric);
535  tw_dynfric->haswork = blackhole_dynfric_haswork;
539  tw_dynfric->query_type_elsize = sizeof(TreeWalkQueryBHDynfric);
540  tw_dynfric->result_type_elsize = sizeof(TreeWalkResultBHDynfric);
541  tw_dynfric->tree = tree;
542  tw_dynfric->priv = priv;
543 
544  /*************************************************************************/
545  TreeWalk tw_accretion[1] = {{0}};
546 
547  tw_accretion->ev_label = "BH_ACCRETION";
549  tw_accretion->ngbiter_type_elsize = sizeof(TreeWalkNgbIterBHAccretion);
551  tw_accretion->haswork = blackhole_dynfric_haswork;
556  tw_accretion->query_type_elsize = sizeof(TreeWalkQueryBHAccretion);
557  tw_accretion->result_type_elsize = sizeof(TreeWalkResultBHAccretion);
558  tw_accretion->tree = tree;
559  tw_accretion->priv = priv;
560 
561  /*************************************************************************/
562 
563  TreeWalk tw_feedback[1] = {{0}};
564  tw_feedback->ev_label = "BH_FEEDBACK";
566  tw_feedback->ngbiter_type_elsize = sizeof(TreeWalkNgbIterBHFeedback);
568  tw_feedback->haswork = blackhole_feedback_haswork;
572  tw_feedback->query_type_elsize = sizeof(TreeWalkQueryBHFeedback);
573  tw_feedback->result_type_elsize = sizeof(TreeWalkResultBHFeedback);
574  tw_feedback->tree = tree;
575  tw_feedback->priv = priv;
576  tw_feedback->repeatdisallowed = 1;
577 
578  priv->atime = atime;
579  priv->a3inv = 1./(atime * atime * atime);
580  priv->hubble = hubble_function(CP, atime);
581  priv->CP = CP;
582 
583  /* Build the queue once, since it is really 'all black holes' and similar for all treewalks*/
584  treewalk_build_queue(tw_dynfric, act->ActiveParticle, act->NumActiveParticle, 0);
585  /* Now we have a BH queue and we can re-use it*/
586  int * ActiveBlackHoles = tw_dynfric->WorkSet;
587  int64_t NumActiveBlackHoles = tw_dynfric->WorkSetSize;
588  /* If this queue is empty, nothing to do.*/
589  MPI_Allreduce(&NumActiveBlackHoles, &totbh, 1, MPI_INT64, MPI_SUM, MPI_COMM_WORLD);
590  if(totbh == 0) {
591  myfree(ActiveBlackHoles);
592  return;
593  }
594 
595  /* We can re-use the current queue for these treewalks*/
596  tw_accretion->haswork = NULL;
597  tw_dynfric->haswork = NULL;
598 
599  /*************************************************************************/
600  /* Dynamical Friction Treewalk */
601 
602  /* Environment variables for DF */
603  priv->BH_SurroundingRmsVel = (MyFloat *) mymalloc("BH_SurroundingRmsVel", SlotsManager->info[5].size * sizeof(priv->BH_SurroundingRmsVel));
604  priv->BH_SurroundingVel = (MyFloat (*) [3]) mymalloc("BH_SurroundingVel", 3* SlotsManager->info[5].size * sizeof(priv->BH_SurroundingVel[0]));
605  priv->BH_SurroundingParticles = (int *)mymalloc("BH_SurroundingParticles", SlotsManager->info[5].size * sizeof(priv->BH_SurroundingParticles));
606  priv->BH_SurroundingDensity = (MyFloat *) mymalloc("BH_SurroundingDensity", SlotsManager->info[5].size * sizeof(priv->BH_SurroundingDensity));
607  /* guard treewalk */
609  treewalk_run(tw_dynfric, ActiveBlackHoles, NumActiveBlackHoles);
610 
611  /*************************************************************************/
612 
613  walltime_measure("/BH/DynFric");
614 
615  /* Let's determine which particles may be swallowed and calculate total feedback weights */
616  priv->SPH_SwallowID = (MyIDType *) mymalloc("SPH_SwallowID", SlotsManager->info[0].size * sizeof(MyIDType));
617  memset(priv->SPH_SwallowID, 0, SlotsManager->info[0].size * sizeof(MyIDType));
618 
619  /* Computed in accretion, used in feedback*/
620  priv->BH_FeedbackWeightSum = (MyFloat *) mymalloc("BH_FeedbackWeightSum", SlotsManager->info[5].size * sizeof(MyFloat));
621 
622  /* These are initialized in preprocess and used to reposition the BH in postprocess*/
623  priv->MinPot = (MyFloat *) mymalloc("BH_MinPot", SlotsManager->info[5].size * sizeof(MyFloat));
624 
625  /* Local to this treewalk*/
626  priv->BH_Entropy = (MyFloat *) mymalloc("BH_Entropy", SlotsManager->info[5].size * sizeof(MyFloat));
627  priv->BH_SurroundingGasVel = (MyFloat (*) [3]) mymalloc("BH_SurroundVel", 3* SlotsManager->info[5].size * sizeof(priv->BH_SurroundingGasVel[0]));
628 
629  /* For AGN kinetic feedback */
630  priv->NumDM = mymalloc("NumDM", SlotsManager->info[5].size * sizeof(MyFloat));
631  priv->V2sumDM = mymalloc("V2sumDM", SlotsManager->info[5].size * sizeof(MyFloat));
632  priv->V1sumDM = (MyFloat (*) [3]) mymalloc("V1sumDM", 3* SlotsManager->info[5].size * sizeof(priv->V1sumDM[0]));
633  priv->MgasEnc = mymalloc("MgasEnc", SlotsManager->info[5].size * sizeof(MyFloat));
634  /* mark the state of AGN kinetic feedback */
635  priv->KEflag = mymalloc("KEflag", SlotsManager->info[5].size * sizeof(int));
636 
637  /* This allocates memory*/
638  treewalk_run(tw_accretion, ActiveBlackHoles, NumActiveBlackHoles);
639 
640  /*************************************************************************/
641 
642  walltime_measure("/BH/Accretion");
643  MPIU_Barrier(MPI_COMM_WORLD);
644 
645  /* Now do the swallowing of particles and dump feedback energy */
646 
647  /* Ionization counters*/
648  priv[0].N_sph_swallowed = ta_malloc("n_sph_swallowed", int64_t, omp_get_max_threads());
649  priv[0].N_BH_swallowed = ta_malloc("n_BH_swallowed", int64_t, omp_get_max_threads());
650  memset(priv[0].N_sph_swallowed, 0, sizeof(int64_t) * omp_get_max_threads());
651  memset(priv[0].N_BH_swallowed, 0, sizeof(int64_t) * omp_get_max_threads());
652 
653  priv->BH_accreted_Mass = (MyFloat *) mymalloc("BH_accretedmass", SlotsManager->info[5].size * sizeof(MyFloat));
654  priv->BH_accreted_BHMass = (MyFloat *) mymalloc("BH_accreted_BHMass", SlotsManager->info[5].size * sizeof(MyFloat));
655  priv->BH_accreted_Mtrack = (MyFloat *) mymalloc("BH_accreted_Mtrack", SlotsManager->info[5].size * sizeof(MyFloat));
656  priv->BH_accreted_momentum = (MyFloat (*) [3]) mymalloc("BH_accretemom", 3* SlotsManager->info[5].size * sizeof(priv->BH_accreted_momentum[0]));
657 
658  treewalk_run(tw_feedback, ActiveBlackHoles, NumActiveBlackHoles);
659 
660  /*************************************************************************/
661  walltime_measure("/BH/Feedback");
662 
663  if(FdBlackholeDetails){
664  collect_BH_info(ActiveBlackHoles, NumActiveBlackHoles, priv, FdBlackholeDetails);
665  }
666 
668  myfree(priv->BH_accreted_Mtrack);
669  myfree(priv->BH_accreted_BHMass);
670  myfree(priv->BH_accreted_Mass);
671 
672  /*****************************************************************/
673  myfree(priv->KEflag);
674  myfree(priv->MgasEnc);
675  myfree(priv->V1sumDM);
676  myfree(priv->V2sumDM);
677  myfree(priv->NumDM);
678 
680  myfree(priv->BH_Entropy);
681  myfree(priv->MinPot);
682 
684  myfree(priv->SPH_SwallowID);
685 
686  /*****************************************************************/
689  myfree(priv->BH_SurroundingVel);
691 
692  /*****************************************************************/
693  myfree(ActiveBlackHoles);
694 
695  int64_t Ntot_gas_swallowed, Ntot_BH_swallowed;
696  int64_t N_sph_swallowed = 0, N_BH_swallowed = 0;
697  for(i = 0; i < omp_get_max_threads(); i++) {
698  N_sph_swallowed += priv[0].N_sph_swallowed[i];
699  N_BH_swallowed += priv[0].N_BH_swallowed[i];
700  }
701  ta_free(priv[0].N_BH_swallowed);
702  ta_free(priv[0].N_sph_swallowed);
703 
704  MPI_Reduce(&N_sph_swallowed, &Ntot_gas_swallowed, 1, MPI_INT64, MPI_SUM, 0, MPI_COMM_WORLD);
705  MPI_Reduce(&N_BH_swallowed, &Ntot_BH_swallowed, 1, MPI_INT64, MPI_SUM, 0, MPI_COMM_WORLD);
706 
707  MPIU_Barrier(MPI_COMM_WORLD);
708  message(0, "Accretion done: %d gas particles swallowed, %d BH particles swallowed\n",
709  Ntot_gas_swallowed, Ntot_BH_swallowed);
710 
711  int total_bh;
712  double total_mdoteddington;
713  double total_mass_holes, total_mdot;
714 
715  double Local_BH_mass = 0;
716  double Local_BH_Mdot = 0;
717  double Local_BH_Medd = 0;
718  int Local_BH_num = 0;
719  /* Compute total mass of black holes
720  * present by summing contents of black hole array*/
721  #pragma omp parallel for reduction(+ : Local_BH_num) reduction(+: Local_BH_mass) reduction(+: Local_BH_Mdot) reduction(+: Local_BH_Medd)
722  for(i = 0; i < SlotsManager->info[5].size; i ++)
723  {
724  if(BhP[i].SwallowID != (MyIDType) -1)
725  continue;
726  Local_BH_num++;
727  Local_BH_mass += BhP[i].Mass;
728  Local_BH_Mdot += BhP[i].Mdot;
729  Local_BH_Medd += BhP[i].Mdot/BhP[i].Mass;
730  }
731 
732  MPI_Reduce(&Local_BH_mass, &total_mass_holes, 1, MPI_DOUBLE, MPI_SUM, 0, MPI_COMM_WORLD);
733  MPI_Reduce(&Local_BH_Mdot, &total_mdot, 1, MPI_DOUBLE, MPI_SUM, 0, MPI_COMM_WORLD);
734  MPI_Reduce(&Local_BH_Medd, &total_mdoteddington, 1, MPI_DOUBLE, MPI_SUM, 0, MPI_COMM_WORLD);
735  MPI_Reduce(&Local_BH_num, &total_bh, 1, MPI_INT, MPI_SUM, 0, MPI_COMM_WORLD);
736 
737  if(FdBlackHoles)
738  {
739  /* convert to solar masses per yr */
740  double mdot_in_msun_per_year =
742 
743  total_mdoteddington *= 1.0 / ((4 * M_PI * GRAVITY * LIGHTCGS * PROTONMASS /
744  (0.1 * LIGHTCGS * LIGHTCGS * THOMPSON)) * units.UnitTime_in_s);
745 
746  fprintf(FdBlackHoles, "%g %d %g %g %g %g\n",
747  atime, total_bh, total_mass_holes, total_mdot, mdot_in_msun_per_year, total_mdoteddington);
748  fflush(FdBlackHoles);
749  }
750  walltime_measure("/BH/Info");
751 }
static void blackhole_accretion_postprocess(int n, TreeWalk *tw)
Definition: blackhole.c:905
static void blackhole_accretion_copy(int place, TreeWalkQueryBHAccretion *I, TreeWalk *tw)
Definition: blackhole.c:1492
static void blackhole_dynfric_reduce(int place, TreeWalkResultBHDynfric *remote, enum TreeWalkReduceMode mode, TreeWalk *tw)
Definition: blackhole.c:846
static void blackhole_feedback_postprocess(int n, TreeWalk *tw)
Definition: blackhole.c:1029
static void blackhole_dynfric_ngbiter(TreeWalkQueryBHDynfric *I, TreeWalkResultBHDynfric *O, TreeWalkNgbIterBHDynfric *iter, LocalTreeWalk *lv)
Definition: blackhole.c:866
static void blackhole_dynfric_copy(int place, TreeWalkQueryBHDynfric *I, TreeWalk *tw)
Definition: blackhole.c:859
static int blackhole_dynfric_haswork(int n, TreeWalk *tw)
Definition: blackhole.c:840
static void blackhole_feedback_ngbiter(TreeWalkQueryBHFeedback *I, TreeWalkResultBHFeedback *O, TreeWalkNgbIterBHFeedback *iter, LocalTreeWalk *lv)
Definition: blackhole.c:1278
static void blackhole_feedback_copy(int place, TreeWalkQueryBHFeedback *I, TreeWalk *tw)
Definition: blackhole.c:1516
static void blackhole_accretion_preprocess(int n, TreeWalk *tw)
Definition: blackhole.c:1017
static void blackhole_accretion_ngbiter(TreeWalkQueryBHAccretion *I, TreeWalkResultBHAccretion *O, TreeWalkNgbIterBHAccretion *iter, LocalTreeWalk *lv)
Definition: blackhole.c:1063
static void blackhole_dynfric_postprocess(int n, TreeWalk *tw)
Definition: blackhole.c:757
static int blackhole_feedback_haswork(int n, TreeWalk *tw)
Definition: blackhole.c:1509
static void blackhole_accretion_reduce(int place, TreeWalkResultBHAccretion *remote, enum TreeWalkReduceMode mode, TreeWalk *tw)
Definition: blackhole.c:1457
static void blackhole_feedback_reduce(int place, TreeWalkResultBHFeedback *remote, enum TreeWalkReduceMode mode, TreeWalk *tw)
Definition: blackhole.c:1543
static void collect_BH_info(int *ActiveBlackHoles, int NumActiveBlackHoles, struct BHPriv *priv, FILE *FdBlackholeDetails)
Definition: blackhole.c:423
double hubble_function(const Cosmology *CP, double a)
Definition: cosmology.c:58
void message(int where, const char *fmt,...)
Definition: endrun.c:175
#define mymalloc(name, size)
Definition: mymalloc.h:15
#define ta_malloc(name, type, nele)
Definition: mymalloc.h:25
#define ta_free(p)
Definition: mymalloc.h:28
#define myfree(x)
Definition: mymalloc.h:19
#define SOLAR_MASS
Definition: physconst.h:6
#define GRAVITY
Definition: physconst.h:5
#define PROTONMASS
Definition: physconst.h:21
#define THOMPSON
Definition: physconst.h:23
#define SEC_PER_YEAR
Definition: physconst.h:32
#define LIGHTCGS
Definition: physconst.h:18
static Cosmology * CP
Definition: power.c:27
struct slots_manager_type SlotsManager[1]
Definition: slotsmanager.c:7
#define BhP
Definition: slotsmanager.h:121
int * ActiveParticle
Definition: timestep.h:13
int64_t NumActiveParticle
Definition: timestep.h:12
double hubble
Definition: blackhole.c:194
int64_t * N_sph_swallowed
Definition: blackhole.c:198
MyFloat(* BH_SurroundingVel)[3]
Definition: blackhole.c:167
double a3inv
Definition: blackhole.c:193
MyFloat * BH_accreted_BHMass
Definition: blackhole.c:176
MyFloat * MinPot
Definition: blackhole.c:159
MyFloat * MgasEnc
Definition: blackhole.c:187
MyFloat * BH_accreted_Mass
Definition: blackhole.c:175
struct UnitSystem units
Definition: blackhole.c:195
MyFloat * BH_SurroundingDensity
Definition: blackhole.c:165
int64_t * N_BH_swallowed
Definition: blackhole.c:199
MyFloat * BH_SurroundingRmsVel
Definition: blackhole.c:168
MyFloat * BH_Entropy
Definition: blackhole.c:160
int * KEflag
Definition: blackhole.c:189
MyFloat(* BH_SurroundingGasVel)[3]
Definition: blackhole.c:161
double atime
Definition: blackhole.c:192
MyIDType * SPH_SwallowID
Definition: blackhole.c:157
MyFloat(* V1sumDM)[3]
Definition: blackhole.c:185
int * BH_SurroundingParticles
Definition: blackhole.c:166
MyFloat * BH_FeedbackWeightSum
Definition: blackhole.c:181
MyFloat * NumDM
Definition: blackhole.c:184
Cosmology * CP
Definition: blackhole.c:196
MyFloat * V2sumDM
Definition: blackhole.c:186
MyFloat * BH_accreted_Mtrack
Definition: blackhole.c:177
MyFloat(* BH_accreted_momentum)[3]
Definition: blackhole.c:172
int BH_DynFrictionMethod
Definition: blackhole.c:46
void * priv
Definition: treewalk.h:85
TreeWalkHasWorkFunction haswork
Definition: treewalk.h:100
int64_t WorkSetSize
Definition: treewalk.h:163
size_t ngbiter_type_elsize
Definition: treewalk.h:97
TreeWalkProcessFunction preprocess
Definition: treewalk.h:105
const ForceTree * tree
Definition: treewalk.h:88
TreeWalkProcessFunction postprocess
Definition: treewalk.h:104
TreeWalkNgbIterFunction ngbiter
Definition: treewalk.h:103
TreeWalkReduceResultFunction reduce
Definition: treewalk.h:102
const char * ev_label
Definition: treewalk.h:91
int * WorkSet
Definition: treewalk.h:161
TreeWalkVisitFunction visit
Definition: treewalk.h:99
int repeatdisallowed
Definition: treewalk.h:117
size_t result_type_elsize
Definition: treewalk.h:96
TreeWalkFillQueryFunction fill
Definition: treewalk.h:101
size_t query_type_elsize
Definition: treewalk.h:95
double UnitMass_in_g
Definition: unitsystem.h:8
double UnitTime_in_s
Definition: unitsystem.h:11
int64_t size
Definition: slotsmanager.h:12
struct slot_info info[6]
Definition: slotsmanager.h:112
#define MPIU_Barrier(comm)
Definition: system.h:103
#define MPI_INT64
Definition: system.h:12
int treewalk_visit_ngbiter(TreeWalkQueryBase *I, TreeWalkResultBase *O, LocalTreeWalk *lv)
Definition: treewalk.c:924
void treewalk_run(TreeWalk *tw, int *active_set, size_t size)
Definition: treewalk.c:619
void treewalk_build_queue(TreeWalk *tw, int *active_set, const size_t size, int may_have_garbage)
Definition: treewalk.c:394
int(* TreeWalkVisitFunction)(TreeWalkQueryBase *input, TreeWalkResultBase *output, LocalTreeWalk *lv)
Definition: treewalk.h:68
void(* TreeWalkReduceResultFunction)(const int j, TreeWalkResultBase *result, const enum TreeWalkReduceMode mode, TreeWalk *tw)
Definition: treewalk.h:76
void(* TreeWalkProcessFunction)(const int i, TreeWalk *tw)
Definition: treewalk.h:73
void(* TreeWalkNgbIterFunction)(TreeWalkQueryBase *input, TreeWalkResultBase *output, TreeWalkNgbIterBase *iter, LocalTreeWalk *lv)
Definition: treewalk.h:70
void(* TreeWalkFillQueryFunction)(const int j, TreeWalkQueryBase *query, TreeWalk *tw)
Definition: treewalk.h:75
uint64_t MyIDType
Definition: types.h:10
LOW_PRECISION MyFloat
Definition: types.h:19
#define walltime_measure(name)
Definition: walltime.h:8

References BHPriv::a3inv, ActiveParticles::ActiveParticle, BHPriv::atime, BHPriv::BH_accreted_BHMass, BHPriv::BH_accreted_Mass, BHPriv::BH_accreted_momentum, BHPriv::BH_accreted_Mtrack, BlackholeParams::BH_DynFrictionMethod, BHPriv::BH_Entropy, BHPriv::BH_FeedbackWeightSum, BHPriv::BH_SurroundingDensity, BHPriv::BH_SurroundingGasVel, BHPriv::BH_SurroundingParticles, BHPriv::BH_SurroundingRmsVel, BHPriv::BH_SurroundingVel, BhP, blackhole_accretion_copy(), blackhole_accretion_ngbiter(), blackhole_accretion_postprocess(), blackhole_accretion_preprocess(), blackhole_accretion_reduce(), blackhole_dynfric_copy(), blackhole_dynfric_haswork(), blackhole_dynfric_ngbiter(), blackhole_dynfric_postprocess(), blackhole_dynfric_reduce(), blackhole_feedback_copy(), blackhole_feedback_haswork(), blackhole_feedback_ngbiter(), blackhole_feedback_postprocess(), blackhole_feedback_reduce(), blackhole_params, collect_BH_info(), BHPriv::CP, CP, TreeWalk::ev_label, TreeWalk::fill, GRAVITY, TreeWalk::haswork, BHPriv::hubble, hubble_function(), slots_manager_type::info, BHPriv::KEflag, LIGHTCGS, message(), BHPriv::MgasEnc, BHPriv::MinPot, MPI_INT64, MPIU_Barrier, myfree, mymalloc, BHPriv::N_BH_swallowed, BHPriv::N_sph_swallowed, TreeWalk::ngbiter, TreeWalk::ngbiter_type_elsize, ActiveParticles::NumActiveParticle, BHPriv::NumDM, TreeWalk::postprocess, TreeWalk::preprocess, TreeWalk::priv, PROTONMASS, TreeWalk::query_type_elsize, TreeWalk::reduce, TreeWalk::repeatdisallowed, TreeWalk::result_type_elsize, SEC_PER_YEAR, slot_info::size, SlotsManager, SOLAR_MASS, BHPriv::SPH_SwallowID, ta_free, ta_malloc, THOMPSON, TreeWalk::tree, treewalk_build_queue(), treewalk_run(), treewalk_visit_ngbiter(), UnitSystem::UnitMass_in_g, BHPriv::units, UnitSystem::UnitTime_in_s, BHPriv::V1sumDM, BHPriv::V2sumDM, TreeWalk::visit, walltime_measure, TreeWalk::WorkSet, and TreeWalk::WorkSetSize.

Referenced by run().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ blackhole_make_one()

void blackhole_make_one ( int  index,
const double  atime 
)

Definition at line 1576 of file blackhole.c.

1576  {
1577  if(P[index].Type != 0)
1578  endrun(7772, "Only Gas turns into blackholes, what's wrong?");
1579 
1580  int child = index;
1581 
1582  /* Make the new particle a black hole: use all the P[i].Mass
1583  * so we don't have lots of low mass tracers.
1584  * If the BH seed mass is small this may lead to a mismatch
1585  * between the gas and BH mass. */
1586  child = slots_convert(child, 5, -1, PartManager, SlotsManager);
1587 
1588  BHP(child).base.ID = P[child].ID;
1589  /* The accretion mass should always be the seed black hole mass,
1590  * irrespective of the gravitational mass of the particle.*/
1592  BHP(child).Mass = bh_powerlaw_seed_mass(P[child].ID);
1593  else
1594  BHP(child).Mass = blackhole_params.SeedBlackHoleMass;
1595 
1596  BHP(child).Mseed = BHP(child).Mass;
1597  BHP(child).Mdot = 0;
1598  BHP(child).FormationTime = atime;
1599  BHP(child).SwallowID = (MyIDType) -1;
1600  BHP(child).Density = 0;
1601 
1602  /* It is important to initialize MinPotPos to the current position of
1603  * a BH to avoid drifting to unknown locations (0,0,0) immediately
1604  * after the BH is created. */
1605  int j;
1606  for(j = 0; j < 3; j++) {
1607  BHP(child).MinPotPos[j] = P[child].Pos[j];
1608  BHP(child).DFAccel[j] = 0;
1609  BHP(child).DragAccel[j] = 0;
1610  }
1611  BHP(child).JumpToMinPot = 0;
1612  BHP(child).CountProgs = 1;
1613 
1615  BHP(child).Mtrack = P[child].Mass;
1616  P[child].Mass = blackhole_params.SeedBHDynMass;
1617  }
1618  else{
1619  BHP(child).Mtrack = -1; /* This column is not used then. */
1620  }
1621  /* Initialize KineticFdbkEnergy, keep zero if BlackHoleKineticOn is not turned on */
1622  BHP(child).KineticFdbkEnergy = 0;
1623 }
static double bh_powerlaw_seed_mass(MyIDType ID)
Definition: blackhole.c:1560
void endrun(int where, const char *fmt,...)
Definition: endrun.c:147
struct part_manager_type PartManager[1]
Definition: partmanager.c:11
#define P
Definition: partmanager.h:88
int slots_convert(int parent, int ptype, int placement, struct part_manager_type *pman, struct slots_manager_type *sman)
Definition: slotsmanager.c:60
#define BHP(i)
Definition: slotsmanager.h:125
double SeedBHDynMass
Definition: blackhole.c:51
double SeedBlackHoleMass
Definition: blackhole.c:53
double MaxSeedBlackHoleMass
Definition: blackhole.c:54

References BHPriv::atime, bh_powerlaw_seed_mass(), BHP, blackhole_params, endrun(), BlackholeParams::MaxSeedBlackHoleMass, P, PartManager, BlackholeParams::SeedBHDynMass, BlackholeParams::SeedBlackHoleMass, slots_convert(), and SlotsManager.

Referenced by fof_seed_make_one().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ set_blackhole_params()

void set_blackhole_params ( ParameterSet ps)

Definition at line 257 of file blackhole.c.

258 {
259  int ThisTask;
260  MPI_Comm_rank(MPI_COMM_WORLD, &ThisTask);
261  if(ThisTask == 0) {
262  blackhole_params.BlackHoleAccretionFactor = param_get_double(ps, "BlackHoleAccretionFactor");
263  blackhole_params.BlackHoleEddingtonFactor = param_get_double(ps, "BlackHoleEddingtonFactor");
264 
265  blackhole_params.BlackHoleFeedbackFactor = param_get_double(ps, "BlackHoleFeedbackFactor");
266 
267  blackhole_params.BlackHoleFeedbackMethod = (enum BlackHoleFeedbackMethod) param_get_enum(ps, "BlackHoleFeedbackMethod");
268  blackhole_params.BlackHoleRepositionEnabled = param_get_int(ps, "BlackHoleRepositionEnabled");
269 
270  blackhole_params.BlackHoleKineticOn = param_get_int(ps,"BlackHoleKineticOn");
271  blackhole_params.BHKE_EddingtonThrFactor = param_get_double(ps, "BHKE_EddingtonThrFactor");
272  blackhole_params.BHKE_EddingtonMFactor = param_get_double(ps, "BHKE_EddingtonMFactor");
273  blackhole_params.BHKE_EddingtonMPivot = param_get_double(ps, "BHKE_EddingtonMPivot");
274  blackhole_params.BHKE_EddingtonMIndex = param_get_double(ps, "BHKE_EddingtonMIndex");
275  blackhole_params.BHKE_EffRhoFactor = param_get_double(ps, "BHKE_EffRhoFactor");
276  blackhole_params.BHKE_EffCap = param_get_double(ps, "BHKE_EffCap");
277  blackhole_params.BHKE_InjEnergyThr = param_get_double(ps, "BHKE_InjEnergyThr");
279  /***********************************************************************************/
280  blackhole_params.BH_DynFrictionMethod = param_get_int(ps, "BH_DynFrictionMethod");
281  blackhole_params.BH_DFBoostFactor = param_get_int(ps, "BH_DFBoostFactor");
282  blackhole_params.BH_DFbmax = param_get_double(ps, "BH_DFbmax");
283  blackhole_params.BH_DRAG = param_get_int(ps, "BH_DRAG");
284  blackhole_params.MergeGravBound = param_get_int(ps, "MergeGravBound");
285  blackhole_params.SeedBHDynMass = param_get_double(ps,"SeedBHDynMass");
286 
287  blackhole_params.SeedBlackHoleMass = param_get_double(ps, "SeedBlackHoleMass");
288  blackhole_params.MaxSeedBlackHoleMass = param_get_double(ps,"MaxSeedBlackHoleMass");
289  blackhole_params.SeedBlackHoleMassIndex = param_get_double(ps,"SeedBlackHoleMassIndex");
290  /***********************************************************************************/
291  }
292  MPI_Bcast(&blackhole_params, sizeof(struct BlackholeParams), MPI_BYTE, 0, MPI_COMM_WORLD);
293 }
BlackHoleFeedbackMethod
Definition: blackhole.h:8
double param_get_double(ParameterSet *ps, const char *name)
Definition: paramset.c:336
int param_get_int(ParameterSet *ps, const char *name)
Definition: paramset.c:368
int param_get_enum(ParameterSet *ps, const char *name)
Definition: paramset.c:378
double BHKE_EffRhoFactor
Definition: blackhole.c:39
enum BlackHoleFeedbackMethod BlackHoleFeedbackMethod
Definition: blackhole.c:30
double BHKE_EddingtonMFactor
Definition: blackhole.c:36
int BlackHoleKineticOn
Definition: blackhole.c:34
int BH_DFBoostFactor
Definition: blackhole.c:47
double BlackHoleAccretionFactor
Definition: blackhole.c:28
double BHKE_SfrCritOverDensity
Definition: blackhole.c:42
double BHKE_EddingtonMPivot
Definition: blackhole.c:37
double SeedBlackHoleMassIndex
Definition: blackhole.c:55
double BH_DFbmax
Definition: blackhole.c:48
double BHKE_EffCap
Definition: blackhole.c:40
int MergeGravBound
Definition: blackhole.c:44
double BlackHoleFeedbackFactor
Definition: blackhole.c:29
double BHKE_EddingtonThrFactor
Definition: blackhole.c:35
double BlackHoleEddingtonFactor
Definition: blackhole.c:31
double BHKE_InjEnergyThr
Definition: blackhole.c:41
double BHKE_EddingtonMIndex
Definition: blackhole.c:38
int ThisTask
Definition: test_exchange.c:23

References BlackholeParams::BH_DFbmax, BlackholeParams::BH_DFBoostFactor, BlackholeParams::BH_DRAG, BlackholeParams::BH_DynFrictionMethod, BlackholeParams::BHKE_EddingtonMFactor, BlackholeParams::BHKE_EddingtonMIndex, BlackholeParams::BHKE_EddingtonMPivot, BlackholeParams::BHKE_EddingtonThrFactor, BlackholeParams::BHKE_EffCap, BlackholeParams::BHKE_EffRhoFactor, BlackholeParams::BHKE_InjEnergyThr, BlackholeParams::BHKE_SfrCritOverDensity, blackhole_params, BlackholeParams::BlackHoleAccretionFactor, BlackholeParams::BlackHoleEddingtonFactor, BlackholeParams::BlackHoleFeedbackFactor, BlackholeParams::BlackHoleFeedbackMethod, BlackholeParams::BlackHoleKineticOn, BlackholeParams::BlackHoleRepositionEnabled, BlackholeParams::MaxSeedBlackHoleMass, BlackholeParams::MergeGravBound, param_get_double(), param_get_enum(), param_get_int(), BlackholeParams::SeedBHDynMass, BlackholeParams::SeedBlackHoleMass, BlackholeParams::SeedBlackHoleMassIndex, and ThisTask.

Referenced by read_parameter_file().

Here is the call graph for this function:
Here is the caller graph for this function: