/******************************************************************************* * This file is part of SWIFT. * Coypright (c) 2018 Matthieu Schaller (schaller@strw.leidenuniv.nl) * * This program is free software: you can redistribute it and/or modify * it under the terms of the GNU Lesser General Public License as published * by the Free Software Foundation, either version 3 of the License, or * (at your option) any later version. * * This program is distributed in the hope that it will be useful, * but WITHOUT ANY WARRANTY; without even the implied warranty of * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the * GNU General Public License for more details. * * You should have received a copy of the GNU Lesser General Public License * along with this program. If not, see . * ******************************************************************************/ #ifndef SWIFT_EAGLE_BLACK_HOLES_PROPERTIES_H #define SWIFT_EAGLE_BLACK_HOLES_PROPERTIES_H #include "chemistry.h" #include "hydro_properties.h" struct black_holes_props { /* ----- Basic neighbour search properties ------ */ /*! Resolution parameter */ float eta_neighbours; /*! Target weightd number of neighbours (for info only)*/ float target_neighbours; /*! Smoothing length tolerance */ float h_tolerance; /*! Tolerance on neighbour number (for info only)*/ float delta_neighbours; /*! Maximal number of iterations to converge h */ int max_smoothing_iterations; /*! Maximal change of h over one time-step */ float log_max_h_change; /* ----- Properties of the accretion model ------ */ /*! Maximal fraction of the Eddington rate allowed. */ float f_Edd; /*! Radiative efficiency of the black holes. */ float epsilon_r; /*! Feedback coupling efficiency of the black holes. */ float epsilon_f; /* ---- Properties of the feedback model ------- */ /*! Temperature increase induced by AGN feedback (Kelvin) */ float AGN_delta_T_desired; }; /** * @brief Initialise the black hole properties from the parameter file. * * For the basic black holes neighbour finding properties we use the * defaults from the hydro scheme if the users did not provide specific * values. * * @param bp The #black_holes_props. * @param phys_const The physical constants in the internal unit system. * @param us The internal unit system. * @param params The parsed parameters. * @param hydro_props The already read-in properties of the hydro scheme. * @param cosmo The cosmological model. */ INLINE static void black_holes_props_init(struct black_holes_props *bp, const struct phys_const *phys_const, const struct unit_system *us, struct swift_params *params, const struct hydro_props *hydro_props, const struct cosmology *cosmo) { /* Read in the basic neighbour search properties or default to the hydro ones if the user did not provide any different values */ /* Kernel properties */ bp->eta_neighbours = parser_get_opt_param_float( params, "BlackHoles:resolution_eta", hydro_props->eta_neighbours); /* Tolerance for the smoothing length Newton-Raphson scheme */ bp->h_tolerance = parser_get_opt_param_float(params, "BlackHoles:h_tolerance", hydro_props->h_tolerance); /* Get derived properties */ bp->target_neighbours = pow_dimension(bp->eta_neighbours) * kernel_norm; const float delta_eta = bp->eta_neighbours * (1.f + bp->h_tolerance); bp->delta_neighbours = (pow_dimension(delta_eta) - pow_dimension(bp->eta_neighbours)) * kernel_norm; /* Number of iterations to converge h */ bp->max_smoothing_iterations = parser_get_opt_param_int(params, "BlackHoles:max_ghost_iterations", hydro_props->max_smoothing_iterations); /* Time integration properties */ const float max_volume_change = parser_get_opt_param_float(params, "BlackHoles:max_volume_change", -1); if (max_volume_change == -1) bp->log_max_h_change = hydro_props->log_max_h_change; else bp->log_max_h_change = logf(powf(max_volume_change, hydro_dimension_inv)); /* Accretion parameters ---------------------------------- */ bp->f_Edd = parser_get_param_float(params, "EAGLEAGN:max_eddington_fraction"); bp->epsilon_r = parser_get_param_float(params, "EAGLEAGN:radiative_efficiency"); bp->epsilon_f = parser_get_param_float(params, "EAGLEAGN:coupling_efficiency"); /* Feedback parameters ---------------------------------- */ bp->AGN_delta_T_desired = parser_get_param_float(params, "EAGLEAGN:AGN_delta_T_K"); } #endif /* SWIFT_EAGLE_BLACK_HOLES_PROPERTIES_H */