// find_horizons.cc -- top level driver for finding apparent horizons // $Header$ // // <<>> // <<>> // AHFinderDirect_find_horizons - top-level driver to find apparent horizons /// /// find_horizon - find a horizon /// print_summary_of_which_horizons_found /// /// do_evaluate_expansions /// do_test_expansion_Jacobian /// #include #include #include #include "util_String.h" #include "cctk.h" #include "cctk_Arguments.h" #include "config.h" #include "stdc.h" #include "../jtutil/util.hh" #include "../jtutil/array.hh" #include "../jtutil/cpm_map.hh" #include "../jtutil/linear_map.hh" #include "../patch/coords.hh" #include "../patch/grid.hh" #include "../patch/fd_grid.hh" #include "../patch/patch.hh" #include "../patch/patch_edge.hh" #include "../patch/patch_interp.hh" #include "../patch/ghost_zone.hh" #include "../patch/patch_system.hh" #include "../elliptic/Jacobian.hh" #include "../gr/gfns.hh" #include "../gr/gr.hh" #include "horizon_sequence.hh" #include "BH_diagnostics.hh" #include "driver.hh" // all the code in this file is inside this namespace namespace AHFinderDirect { using jtutil::error_exit; //****************************************************************************** // // ***** access to persistent data ***** // extern struct state state; //****************************************************************************** // // ***** prototypes for functions local to this file // namespace { void print_summary_of_which_horizons_found (int N_horizons, const struct AH_data* const* AH_data_array); void do_evaluate_expansions(int my_proc, int N_horizons, horizon_sequence& hs, struct AH_data* const AH_data_array[], const struct cactus_grid_info& cgi, const struct geometry_info& gi, const struct IO_info& IO_info, const struct error_info& error_info, const struct verbose_info& verbose_info, int timer_handle); void do_test_expansion_Jacobians(int my_proc, int N_horizons, struct AH_data* const AH_data_array[], const struct cactus_grid_info& cgi, const struct geometry_info& gi, struct Jacobian_info& Jac_info, bool test_all_Jacobian_compute_methods, const struct IO_info& IO_info, const struct error_info& error_info, const struct verbose_info& verbose_info, int timer_handle); } //****************************************************************************** // // This function is called by the Cactus scheduler to find the apparent // horizon or horizons in the current slice. // extern "C" void AHFinderDirect_find_horizons(CCTK_ARGUMENTS) { DECLARE_CCTK_ARGUMENTS // only try to find horizons every find_every time steps if (! state.find_now(cctk_iteration)) then return; // *** NO-OP RETURN *** if (state.timer_handle >= 0) then CCTK_TimerResetI(state.timer_handle); const int my_proc = state.my_proc; const int N_horizons = state.N_horizons; horizon_sequence& hs = *state.my_hs; const bool active_flag = hs.has_genuine_horizons(); const bool broadcast_horizon_shape = state.always_broadcast_horizon_shape || state.mask_info.set_mask_for_any_horizon || state.store_info_in_SS_vars_for_any_horizon; struct cactus_grid_info& cgi = state.cgi; const struct geometry_info& gi = state.gi; struct Jacobian_info& Jac_info = state.Jac_info; struct IO_info& IO_info = state.IO_info; const struct error_info& error_info = state.error_info; const struct verbose_info& verbose_info = state.verbose_info; // what are the semantics of the Cactus gxx variables? (these may // change from one call to another, so we have to re-check each time) if (CCTK_Equals(state.metric_type, "physical")) then cgi.use_Cactus_conformal_metric = false; else if (CCTK_Equals(state.metric_type, "static conformal")) then cgi.use_Cactus_conformal_metric = (*conformal_state > 0); else CCTK_VWarn(FATAL_ERROR, __LINE__, __FILE__, CCTK_THORNSTRING, "AHFinderDirect_find_horizons(): unknown ADMBase::metric_type=\"%s\"!", state.metric_type); /*NOTREACHED*/ // get the Cactus time step and decide if we want to output h and/or Theta now IO_info.time_iteration = cctk_iteration; IO_info.time = cctk_time; IO_info.output_h = (IO_info.output_h_every > 0) && ((IO_info.time_iteration % IO_info.output_h_every) == 0); IO_info.output_Theta = (IO_info.output_Theta_every > 0) && ((IO_info.time_iteration % IO_info.output_Theta_every) == 0); // set initial guess for any (genuine) horizons that need it, // i.e. for any (genuine) horizons where we didn't find the horizon previously for (int hn = hs.init_hn() ; hs.is_genuine() ; hn = hs.next_hn()) { assert( state.AH_data_array[hn] != NULL ); struct AH_data& AH_data = *state.AH_data_array[hn]; if (AH_data.found_flag) then AH_data.initial_find_flag = false; else { patch_system& ps = *AH_data.ps_ptr; setup_initial_guess(ps, AH_data.initial_guess_info, IO_info, hn, N_horizons, verbose_info); if (active_flag && IO_info.output_initial_guess) then output_gridfn(ps, gfns::gfn__h, IO_info, IO_info.h_base_file_name, hn, verbose_info .print_algorithm_highlights); AH_data.initial_find_flag = true; } } // // now the main horizon finding (or other computation) // switch (state.method) { case method__evaluate_expansions: do_evaluate_expansions(my_proc, N_horizons, *state.my_hs, state.AH_data_array, cgi, gi, IO_info, error_info, verbose_info, state.timer_handle); break; case method__test_expansion_Jacobians: do_test_expansion_Jacobians(my_proc, N_horizons, state.AH_data_array, cgi, gi, Jac_info, state.test_all_Jacobian_compute_methods, IO_info, error_info, verbose_info, state.timer_handle); break; case method__find_horizons: { if (state.timer_handle >= 0) then CCTK_TimerStartI(state.timer_handle); Newton(cctkGH, state.N_procs, state.N_active_procs, my_proc, *state.my_hs, state.AH_data_array, cgi, gi, Jac_info, state.solver_info, IO_info, state.BH_diagnostics_info, broadcast_horizon_shape, error_info, verbose_info, state.isb); if (state.timer_handle >= 0) then CCTK_TimerStopI(state.timer_handle); if (verbose_info.print_physics_highlights && !verbose_info.print_physics_details) then print_summary_of_which_horizons_found(N_horizons, state.AH_data_array); break; } default: CCTK_VWarn(FATAL_ERROR, __LINE__, __FILE__, CCTK_THORNSTRING, "\n" " find_horizons(): unknown method=(int)%d!\n" " (this should never happen!)" , int(state.method)); /*NOTREACHED*/ } if (state.timer_handle >= 0) then { CCTK_VInfo(CCTK_THORNSTRING, "timer stats for computation:"); CCTK_TimerPrintDataI(state.timer_handle, -1); } } //****************************************************************************** // // This function prints (using CCTK_VInfo()) a 1-line summary of // which AHs were / were not found. // namespace { void print_summary_of_which_horizons_found (int N_horizons, const struct AH_data* const* AH_data_array) { const int hn_buffer_size = 10; // buffer for single "%d" for hn char hn_buffer[hn_buffer_size]; const int msg_buffer_size = 200; // buffer for the entire line char msg_buffer[msg_buffer_size]; msg_buffer[0] = '\0'; Util_Strlcpy(msg_buffer, "found horizon(s) ", msg_buffer_size);; bool already_found_flag = false; // running inclusive-or of found_flag for (int hn = 1 ; hn <= N_horizons ; ++hn) { assert(AH_data_array[hn] != NULL); bool found_flag = AH_data_array[hn]->found_flag; if (found_flag) then { if (already_found_flag) then Util_Strlcat(msg_buffer, ",", msg_buffer_size); already_found_flag = true; snprintf(hn_buffer, hn_buffer_size, "%d", hn); Util_Strlcat(msg_buffer, hn_buffer, msg_buffer_size); } } if (already_found_flag) then { // we found one or more horizons Util_Strlcat(msg_buffer, " of ", msg_buffer_size); snprintf(hn_buffer, hn_buffer_size, "%d", N_horizons); Util_Strlcat(msg_buffer, hn_buffer, msg_buffer_size); } else Util_Strlcpy(msg_buffer, "no horizons found", msg_buffer_size); CCTK_VInfo(CCTK_THORNSTRING, msg_buffer); } } //****************************************************************************** //****************************************************************************** //****************************************************************************** // // This function implements AHFinderDirect::method == "horizon function": // On processor #0 it evaluates the Theta(h) function for each apparent // horizon (and does any I/O desired); on other processors it does N_horizons // dummy evaluations on horizon #0. // // Note that if we decide to output h, we output it *after* any Theta(h) // evaluation or horizon finding has been done, to ensure that all the // ghost zones are filled in in case we need to print them. // // Arguments: // timer_handle = a valid Cactus timer handle if we want to time the // apparent horizon process, or -ve to skip this // (we only time the computation, not the file I/O) // namespace { void do_evaluate_expansions(int my_proc, int N_horizons, horizon_sequence& hs, struct AH_data* const AH_data_array[], const struct cactus_grid_info& cgi, const struct geometry_info& gi, const struct IO_info& IO_info, const struct error_info& error_info, const struct verbose_info& verbose_info, int timer_handle) { const bool active_flag = (my_proc == 0); if (active_flag) then { assert( hs.N_horizons() == N_horizons ); assert( hs.my_N_horizons() == N_horizons ); for (int hn = hs.init_hn() ; hs.is_genuine() ; hn = hs.next_hn()) { assert( AH_data_array[hn] != NULL ); struct AH_data& AH_data = *AH_data_array[hn]; patch_system& ps = *AH_data.ps_ptr; if (timer_handle >= 0) then CCTK_TimerStartI(timer_handle); jtutil::norm Theta_norms; const bool Theta_ok = expansion(&ps, -AH_data.surface_expansion, cgi, gi, error_info, true,// initial eval false, // no Jacobian coeffs true, // yes, print msgs &Theta_norms); if (timer_handle >= 0) then CCTK_TimerStopI(timer_handle); if (IO_info.output_h) then output_gridfn(ps, gfns::gfn__h, IO_info, IO_info.h_base_file_name, hn, verbose_info.print_algorithm_details); if (Theta_ok) then { CCTK_VInfo(CCTK_THORNSTRING, " Theta(h) rms-norm %.2e, infinity-norm %.2e", Theta_norms.rms_norm(), Theta_norms.infinity_norm()); if (IO_info.output_Theta) then output_gridfn(ps, gfns::gfn__Theta, IO_info, IO_info .Theta_base_file_name, hn, verbose_info .print_algorithm_details); } } } else { for (int i = 0 ; i < N_horizons ; ++i) { expansion(NULL, 0.0, // dummy computation cgi, gi, error_info, true); // initial evaluation } } } } //****************************************************************************** // // This function implements // AHFinderDirect::method == "test expansion Jacobians": // On processor #0 it computes and prints the Jacobian matrix J[Theta(h)] // function for horizon #1; on other processors it does dummy Jacobian // computations. // // The Jacobian computation may optionally be done in several different // ways, in which case all the resulting Jacobian matrices are printed, // as are their differences. Alternatively, only // the numerical perturbation computation may be done/printed. // // Arguments: // timer_handle = a valid Cactus timer handle if we want to time the // apparent horizon process, or -ve to skip this // (we only time the computation, not the file I/O) // test_all_Jacobian_compute_methods // = true ==> Test all known methods of computing the Jacobian // matrix, and print all the resulting Jacobian matrices // and their differences. // false ==> Just test/print the numerical perturbation calculation. // (This may be useful if one or more of the other methods // is broken.) // namespace { void do_test_expansion_Jacobians(int my_proc, int N_horizons, struct AH_data* const AH_data_array[], const struct cactus_grid_info& cgi, const struct geometry_info& gi, struct Jacobian_info& Jac_info, bool test_all_Jacobian_compute_methods, const struct IO_info& IO_info, const struct error_info& error_info, const struct verbose_info& verbose_info, int timer_handle) { const bool active_flag = (my_proc == 0); assert(N_horizons >= 1); const bool print_msg_flag = true; const int hn = 1; struct AH_data* const AH_data_ptr = active_flag ? AH_data_array[hn] : NULL; patch_system* const ps_ptr = active_flag ? AH_data_ptr->ps_ptr : NULL; const fp add_to_expansion = active_flag ? -AH_data_ptr->surface_expansion : 0.0; // // numerical-perturbation Jacobian // Jacobian* Jac_NP_ptr = active_flag ? AH_data_ptr->Jac_ptr : NULL; expansion(ps_ptr, add_to_expansion, cgi, gi, error_info, true); // initial evaluation Jac_info.Jacobian_compute_method = Jacobian__numerical_perturbation; expansion_Jacobian(ps_ptr, Jac_NP_ptr, add_to_expansion, cgi, gi, Jac_info, error_info, true, // initial evaluation print_msg_flag); Jacobian* Jac_SD_FDdr_ptr = NULL; if (test_all_Jacobian_compute_methods) then { // symbolic differentiation with finite diff d/dr Jac_SD_FDdr_ptr = active_flag ? new_Jacobian(Jac_info.Jacobian_store_solve_method, *ps_ptr, verbose_info.print_algorithm_details) : NULL; expansion(ps_ptr, add_to_expansion, cgi, gi, error_info, true, // initial evaluation true); // compute SD Jacobian coeffs Jac_info.Jacobian_compute_method = Jacobian__symbolic_diff_with_FD_dr; expansion_Jacobian(ps_ptr, Jac_SD_FDdr_ptr, add_to_expansion, cgi, gi, Jac_info, error_info, true, // initial evaluation print_msg_flag); } if (active_flag) then output_Jacobians(*ps_ptr, Jac_NP_ptr, Jac_SD_FDdr_ptr, IO_info, IO_info.Jacobian_base_file_name, hn, print_msg_flag); } } //****************************************************************************** } // namespace AHFinderDirect