#define ETA_S_CORR 4.0 #define ETA_BH_CORR 4.0 #define DTMAXPOWER -3.0 #define DTMINPOWER -36.0 #include #include #include #include #include #include #include #include #include #include #include #include #include "black_holes.h" #include "config.h" #include "double3.h" #include "external.h" #include "grape6.h" #include "io.h" #ifdef ETICS #include "grapite.h" #endif //chrono::steady_clock::time_point walltime_start; namespace std::chrono { struct Timer { void start() { t_start = steady_clock::now(); } void stop() { t_stop = steady_clock::now(); time = duration_cast(t_stop - t_start).count()*1E-9; } double time; // seconds steady_clock::time_point t_start, t_stop; }; } std::chrono::Timer timer; class Calc_self_grav { public: Calc_self_grav(const int N, const int n_loc, const int clusterid, const int npipe, const double eps) : g6_calls(0), n_loc(n_loc), clusterid(clusterid), npipe(npipe), eps2(eps*eps) { h2.assign(N, eps2); pot_loc.resize(N); acc_loc.resize(N); jrk_loc.resize(N); } void operator()(const double t, const int n_act, std::vector &ind_act, std::vector &x_act, std::vector &v_act, std::vector& pot, double3 acc[], double3 jrk[]) { g6_set_ti(clusterid, t); for (int i=0; i h2; std::vector pot_loc; // the _loc variables are for this node only. std::vector acc_loc, jrk_loc; }; void calc_ext_grav(std::vector &external_gravity_components, int n, const std::vector &x, const std::vector &v, double *pot, double3 *acc, double3* jrk) // TODO should just be a class that has this pointer array as a member { std::fill(pot, pot+n, 0.); for (auto component : external_gravity_components) { if (component->is_active) component->apply(n, x, v, pot, acc, jrk); } } void energy_contr(const double time_cur, const double timesteps, const double n_act_sum, const double g6_calls, int N, const std::vector &m, const std::vector &x, const std::vector &v, const std::vector &pot, double pot_ext[]) { double E_pot = 0; for (int i=0; i 0) for (int i=0; i m; std::vector x, v; if (is_hdf5(config.input_file_name)) { #ifndef HAS_HDF5 fprintf(stderr, "ERROR: input file is in HDF5 format, but the code was compiled without HDF5 support\n"); return -1; #endif h5_read(config.input_file_name, &diskstep, &N, &time_cur, m.data(), x.data(), v.data()); } else ascii_read(config.input_file_name, diskstep, N, time_cur, m, x, v); std::vector ind(N); std::iota(begin(ind), end(ind), 0); double3 *a = new double3[N], *adot = new double3[N]; std::vector pot(N); double *pot_ext = new double[N], *t = new double[N], *dt = new double[N]; /* data for active particles */ int n_act; std::vector ind_act(N); std::vector pot_act_new(N); std::vector x_act_new(N), v_act_new(N); double *pot_act_ext = new double[N]; double3 *a_act_new = new double3[N], *adot_act_new = new double3[N]; double eps = config.eps; double eta = config.eta; double t_end = config.t_end; double dt_disk = config.dt_disk; double dt_contr = config.dt_contr; double dt_bh = config.dt_bh; if (myRank == rootRank) { printf("\n"); printf("Begin the calculation of phi-GRAPE program on %03d processors\n", n_proc); printf("\n"); printf("N = %07d \t eps = %.6E \n", N, eps); printf("t_beg = %.6E \t t_end = %.6E \n", time_cur, t_end); printf("dt_disk = %.6E \t dt_contr = %.6E \n", dt_disk, dt_contr); printf("dt_bh = %.6E \n", dt_bh); printf("eta = %.6E \n", eta); printf("\n"); if ((diskstep == 0) && (time_cur == 0)) { FILE *out = fopen("contr.dat", "w"); fclose(out); if (config.live_smbh_output && (config.live_smbh_count > 0)) { out = fopen("bh.dat", "w"); fclose(out); } if ((config.live_smbh_neighbor_output) && (config.live_smbh_count > 0)) { out = fopen("bh_neighbors.dat", "w"); fclose(out); } } } /* if (myRank == rootRank) */ double normalization_mass=1, normalization_length=1, normalization_velocity=1; if (config.ext_units_physical) { normalization_mass = 1/config.unit_mass; normalization_length = 1000/config.unit_length; normalization_velocity = 1.52484071426404437233e+01*sqrt(config.unit_length/config.unit_mass); } Plummer ext_bulge(config.ext_m_bulge*normalization_mass, config.ext_b_bulge*normalization_length); Miyamoto_Nagai ext_disk(config.ext_m_disk*normalization_mass, config.ext_a_disk*normalization_length, config.ext_b_disk*normalization_length); Plummer ext_halo_plummer(config.ext_m_halo_plummer*normalization_mass, config.ext_b_halo_plummer*normalization_length); Logarithmic_halo ext_log_halo(config.ext_log_halo_v*normalization_velocity, config.ext_log_halo_r*normalization_length); Dehnen ext_dehnen(config.ext_dehnen_m*normalization_mass, config.ext_dehnen_r*normalization_length, config.ext_dehnen_gamma); std::vector external_gravity_components; external_gravity_components.push_back(&ext_bulge); external_gravity_components.push_back(&ext_disk); external_gravity_components.push_back(&ext_halo_plummer); external_gravity_components.push_back(&ext_log_halo); external_gravity_components.push_back(&ext_dehnen); bool has_external_gravity = false; for (auto component : external_gravity_components) { if (component->is_active) { has_external_gravity = true; break; } } if (has_external_gravity) printf("External Potential: \n\n"); if (ext_bulge.is_active) printf("m_bulge = %.4E b_bulge = %.4E\n", config.ext_m_bulge*normalization_mass, config.ext_b_bulge*normalization_length); if (ext_disk.is_active) printf("m_disk = %.4E a_disk = %.4E b_disk = %.4E\n", config.ext_m_disk*normalization_mass, config.ext_a_disk*normalization_length, config.ext_b_disk*normalization_length); if (ext_halo_plummer.is_active) printf("m_halo = %.4E b_halo = %.4E\n", config.ext_m_halo_plummer*normalization_mass, config.ext_b_halo_plummer*normalization_length); if (ext_log_halo.is_active) printf("v_halo = %.6E r_halo = %.6E \n", config.ext_log_halo_v*normalization_velocity, config.ext_log_halo_r*normalization_length); if (ext_dehnen.is_active) printf("m_ext = %.6E r_ext = %.6E \t g_ext = %.3E \n", config.ext_dehnen_m*normalization_mass, config.ext_dehnen_r*normalization_length, config.ext_dehnen_gamma); printf("\n"); fflush(stdout); double eta_bh = eta/ETA_BH_CORR; const double dt_min = pow(2, DTMINPOWER); const double dt_max = std::min({dt_disk, dt_contr, pow(2, DTMAXPOWER)}); double t_disk = dt_disk*(1+floor(time_cur/dt_disk)); double t_contr = dt_contr*(1+floor(time_cur/dt_contr)); double t_bh = dt_bh*(1+floor(time_cur/dt_bh)); if (myRank == rootRank) { printf("t_disk = %.6E t_contr = %.6E t_bh = %.6E \n", t_disk, t_contr, t_bh); printf("\n"); fflush(stdout); } /* if (myRank == rootRank) */ std::fill(t, t+N, time_cur); std::fill(dt, dt+N, dt_min); /* some local settings for G6a boards */ int clusterid, numGPU; if (config.devices_per_node==0) { MPI_Comm shmcomm; MPI_Comm_split_type(MPI_COMM_WORLD, MPI_COMM_TYPE_SHARED, 0, MPI_INFO_NULL, &shmcomm); MPI_Comm_size(shmcomm, &numGPU); MPI_Comm_rank(shmcomm, &clusterid); } else { numGPU = config.devices_per_node; clusterid = myRank % numGPU; } printf("Rank of the processor %03d : Number of GPUs %01d : Cluster ID %01d \n", myRank, numGPU, clusterid); fflush(stdout); /* init the local GRAPEs */ g6_open(clusterid); int npipe = g6_npipes(); g6_set_tunit(51); g6_set_xunit(51); #ifdef ETICS grapite_set_dev_exec_threshold(config.grapite_dev_exec_threshold); bool grapite_active_search_flag = config.grapite_active_search; #else bool grapite_active_search_flag = false; #endif int n_loc = N/n_proc; Calc_self_grav calc_self_grav(N, n_loc, clusterid, npipe, eps); Active_search active_search(myRank, n_proc, n_loc, N, grapite_active_search_flag); Black_hole_physics black_hole_physics; if (config.live_smbh_count == 1) black_hole_physics = Black_hole_physics(m[0], 0, myRank, rootRank); else if (config.live_smbh_count == 2) black_hole_physics = Black_hole_physics(m[0], m[1], myRank, rootRank); if (config.binary_smbh_pn) { black_hole_physics.set_post_newtonian(config.pn_c, config.pn_usage.data()); if (config.pn_usage[6]) black_hole_physics.set_spins(config.smbh1_spin.data(), config.smbh2_spin.data()); } black_hole_physics.set_softening(eps, config.live_smbh_custom_eps); Binary_smbh_influence_sphere_output binary_smbh_influence_sphere_output(config.binary_smbh_influence_radius_factor, N, m, x, v, pot, dt); Write_bh_nb_data write_bh_nb_data(config.live_smbh_neighbor_number, config.live_smbh_count, N, m, x, v); #ifdef ETICS grapite_read_particle_tags(N, config.grapite_mask_file_name.c_str(), myRank, n_loc); grapite_set_dt_exp(config.dt_scf); grapite_set_t_exp(time_cur); #endif /* load the nj particles to the G6 */ for (int k=0; k= 0) { double3 xcm, vcm, xdc, vdc; grapite_calc_center(N, m, (double(*)[3])x, (double(*)[3])v, xcm, vcm, xdc, vdc); x[grapite_cep_index] = xdc; v[grapite_cep_index] = vdc; grapite_update_cep(time_cur, xdc, vdc, zeros, zeros); } #endif /* define the all particles as a active on all the processors for the first time grav calc. */ calc_self_grav(time_cur, N, ind, x, v, pot, a, adot); if (config.live_smbh_count == 2) { black_hole_physics.set_xv(x[0], x[1], v[0], v[1]); if (config.live_smbh_custom_eps >= 0) black_hole_physics.adjust_softening(pot[0], pot[1], a[0], a[1], adot[0], adot[1]); if (config.binary_smbh_pn) black_hole_physics.adjust_post_newtonian(dt[0], a[0], a[1], adot[0], adot[1]); } calc_ext_grav(external_gravity_components, N, x, v, pot_ext, a, adot); /* Energy control... */ if (myRank == rootRank) energy_contr(time_cur, timesteps, n_act_sum, calc_self_grav.g6_calls, N, m, x, v, pot, pot_ext); #ifdef ETICS if (config.etics_dump_coeffs && (diskstep==0)) { char out_fname[256]; sprintf(out_fname, "coeffs.%06d.%02d.dat", 0, myRank); grapite_dump(out_fname, 2); } if (grapite_cep_index >= 0) { double3 xcm, vcm, xdc, vdc; grapite_calc_center(N, m, (double(*)[3])x, (double(*)[3])v, xcm, vcm, xdc, vdc); x[grapite_cep_index] = xdc; v[grapite_cep_index] = vdc; grapite_update_cep(time_cur, xdc, vdc, a[grapite_cep_index], adot[grapite_cep_index]); } #endif /* Define initial timestep for all particles on all nodes */ for (int j=0; j dt_max) dt_tmp = dt_max; dt[j] = dt_tmp; if (config.dt_min_warning && (myRank == 0)) { if (dt[j] == dt_min) { printf("!!! Warning0: dt = dt_min = %.6E \t ind = %07d \n", dt[j], ind[j]); fflush(stdout); } } } /* j */ if (config.live_smbh_count > 0) { double min_dt = *std::min_element(dt, dt+N); for (int i=0; i0)) { int act_def_grapite_bh_count = 0; int i_bh[n_bh]; for (int i=0; i= 0) black_hole_physics.adjust_softening(pot_act_new[i_bh1], pot_act_new[i_bh2], a_act_new[i_bh1], a_act_new[i_bh2], adot_act_new[i_bh1], adot_act_new[i_bh2]); if (config.binary_smbh_pn) black_hole_physics.adjust_post_newtonian(dt[i_bh1], a_act_new[i_bh1], a_act_new[i_bh2], adot_act_new[i_bh1], adot_act_new[i_bh2]); } /* Calculate gravity on active particles due to external forces */ calc_ext_grav(external_gravity_components, n_act, x_act_new, v_act_new, pot_act_ext, a_act_new, adot_act_new); /* correct the active particles positions etc... on all the nodes */ double min_dt = dt_max; for (int i=0; i 0) && (ind_act[i] < config.live_smbh_count)) eta_curr = eta_bh; else eta_curr = eta; double dt_new = aarseth_step(eta_curr, dt_tmp, a_act_new[i], adot_act_new[i], a2, a3); //TODO the below should be moved to a function if (dt_new < dt_min) dt_tmp = dt_min; if ((dt_new < dt_tmp) && (dt_new > dt_min)) { int power = log(dt_new)/log(2.0) - 1; dt_tmp = pow(2.0, (double)power); // TODO why is this casting needed here? } if ((dt_new > 2*dt_tmp) && (fmod(min_t, 2*dt_tmp) == 0) && (2*dt_tmp <= dt_max)) { dt_tmp *= 2; } if (config.dt_min_warning && (myRank == 0)) { if (dt_tmp == dt_min) { printf("!!! Warning1: dt_act = dt_min = %.6E \t ind_act = %07d \n", dt_tmp, ind_act[i]); fflush(stdout); } } if (dt_tmp < min_dt) min_dt = dt_tmp; x[j_act] = x_act_new[i]; v[j_act] = v_act_new[i]; t[j_act] = min_t; dt[j_act] = dt_tmp; pot[j_act] = pot_act_new[i]; pot_ext[j_act] = pot_act_ext[i]; a[j_act] = a_act_new[i]; adot[j_act] = adot_act_new[i]; } /* i */ /* define the min. dt over all the act. part. and set it also for the BH... */ if (config.live_smbh_count > 0) { if (config.live_smbh_count>=1) dt[0] = min_dt; if (config.live_smbh_count==2) dt[1] = min_dt; } if (config.binary_smbh_influence_sphere_output && (myRank == rootRank)) { //TODO clean up this mass. We don't want to have all these _act arrays; they are only needed for this lame function. binary_smbh_influence_sphere_output(ind_act, n_act, timesteps, time_cur); } /* load the new values for active particles to the local GRAPE's */ for (int i=0; i= t_bh) { if (myRank == rootRank) { /* Write BH data... */ if (config.live_smbh_output) black_hole_physics.write_bh_data(time_cur, m, x, v, pot, a, adot, dt); /* Write BH NB data... */ if (config.live_smbh_neighbor_output) write_bh_nb_data(time_cur); } /* if (myRank == rootRank) */ t_bh += dt_bh; } /* if (time_cur >= t_bh) */ if (time_cur >= t_contr) { if (myRank == rootRank) { energy_contr(time_cur, timesteps, n_act_sum, calc_self_grav.g6_calls, N, m, x, v, pot, pot_ext); /* write cont data */ if (config.output_hdf5) h5_write("data.con", diskstep, N, time_cur, m, x, v, pot, a, adot, 0, true); else ascii_write("data.con", diskstep, N, time_cur, m, x, v, 16); } /* if (myRank == rootRank) */ #ifdef ETICS // We are /inside/ a control step, so all particles must be // synchronized; we can safely calculate their density centre. The // acceleration and jerk currently in the memory are for the // predicted position of the CEP, by calling grapite_calc_center we // "correct" the position and velocity, but not the gravity at that // point. if (grapite_cep_index >= 0) { double3 xcm, vcm, xdc, vdc; grapite_calc_center(N, m, (double(*)[3])x, (double(*)[3])v, xcm, vcm, xdc, vdc); x[grapite_cep_index] = xdc; v[grapite_cep_index] = vdc; grapite_update_cep(time_cur, xdc, vdc, a[grapite_cep_index], adot[grapite_cep_index]); } #endif t_contr += dt_contr; } /* if (time_cur >= t_contr) */ if (time_cur >= t_disk) { char out_fname[256]; diskstep++; if (myRank == rootRank) { sprintf(out_fname, "%06d", diskstep); if (config.output_hdf5) h5_write(std::string(out_fname) + ".h5", diskstep, N, time_cur, m, x, v, pot, a, adot, config.output_extra_mode, config.output_hdf5_double_precision); else ascii_write(std::string(out_fname) + ".dat", diskstep, N, time_cur, m, x, v, config.output_ascii_precision); } /* if (myRank == rootRank) */ #ifdef ETICS if (config.etics_dump_coeffs) { sprintf(out_fname, "coeffs.%06d.%02d.dat", diskstep, myRank); grapite_dump(out_fname, 2); } #endif t_disk += dt_disk; } /* if (time_cur >= t_disk) */ } /* while (time_cur < t_end) */ /* close the local GRAPEs */ timer.stop(); g6_close(clusterid); double g6_calls_sum; MPI_Reduce(&calc_self_grav.g6_calls, &g6_calls_sum, 1, MPI_DOUBLE, MPI_SUM, rootRank, MPI_COMM_WORLD); if (myRank == rootRank) { /* Write some output for the timestep annalize... */ printf("\n"); printf("timesteps = %.0f Total sum of integrated part. = %.0f g6_calls on all nodes = %.0f \n", timesteps, n_act_sum, g6_calls_sum); printf("\n"); printf("Real Speed = %.3f GFlops \n", 57.0*N*n_act_sum/(timer.time)/1.0E+09); fflush(stdout); } /* if (myRank == rootRank) */ delete[] a; delete[] adot; delete[] pot_ext; delete[] t; delete[] dt; delete[] a_act_new; delete[] adot_act_new; delete[] pot_act_ext; /* Finalize the MPI work */ MPI_Finalize(); return 0; }