// horizon_Jacobian.cc -- evaluate Jacobian matrix of LHS function H(h) // $Id$ // // <<>> // create_Jacobian - create a Jacobian matrix of the appropriate type // horizon_Jacobian - top-level driver to compute the Jacobian /// add_ghost_zone_Jacobian - add ghost zone dependencies to Jacobian // #include #include #include #include #include "util_Table.h" #include "cctk.h" #include "cctk_Arguments.h" #include "stdc.h" #include "config.hh" #include "../jtutil/util.hh" #include "../jtutil/array.hh" #include "../jtutil/cpm_map.hh" #include "../jtutil/linear_map.hh" using jtutil::error_exit; #include "../util/coords.hh" #include "../util/grid.hh" #include "../util/fd_grid.hh" #include "../util/patch.hh" #include "../util/patch_edge.hh" #include "../util/patch_interp.hh" #include "../util/ghost_zone.hh" #include "../util/patch_system.hh" #include "../elliptic/Jacobian.hh" #include "gfn.hh" #include "AHFinderDirect.hh" //****************************************************************************** // // ***** prototypes for functions local to this file ***** // namespace { void add_ghost_zone_Jacobian(const patch_system& ps, const patch& xp, const ghost_zone& xmgz, int x_II, int xm_irho, int xm_isigma, fp mol, Jacobian& Jac); } //****************************************************************************** // // This function decodes the Jacobian_type parameter and creates the // appropriate type of Jacobian matrix object. // Jacobian& create_Jacobian(const patch_system& ps, const char Jacobian_type[]) { if (STRING_EQUAL(Jacobian_type, "dense matrix")) then return *new dense_Jacobian(ps); else CCTK_VWarn(-1, __LINE__, __FILE__, CCTK_THORNSTRING, "unknown Jacobian_type=\"%s\"!", Jacobian_type); /*NOTREACHED*/ } //****************************************************************************** // // This function computes the Jacobian matrix of the LHS function H(h) // from the Jacobian coefficient (angular) gridfns. The computation is // done a Jacobian row at a time, using equation (25) of my 1996 apparent // horizon finding paper. // // Inputs (angular gridfns, on ghosted grid): // h # shape of trial surface // partial_H_wrt_partial_d_h, # Jacobian coefficients // partial_H_wrt_partial_dd_h, // // Outputs: // The Jacobian matrix is stored in the Jacobian object Jac. // void horizon_Jacobian(patch_system& ps, Jacobian& Jac) { CCTK_VInfo(CCTK_THORNSTRING, " horizon Jacobian"); ps.compute_synchronize_Jacobian(); Jac.zero(); for (int xpn = 0 ; xpn < ps.N_patches() ; ++xpn) { patch& xp = ps.ith_patch(xpn); for (int x_irho = xp.min_irho() ; x_irho <= xp.max_irho() ; ++x_irho) { for (int x_isigma = xp.min_isigma() ; x_isigma <= xp.max_isigma() ; ++x_isigma) { // // compute the Jacobian row for this grid point, i.e. // partial H(this point x, Jacobian row II) // --------------------------------------------- // partial h(other points y, Jacobian column JJ) // // FIXME FIXME: we still have to take into account the // position-dependence of the coefficients, // cf the difference between J[3H(h)] and J[2H(h)]; // here we're sort of computing the former, but // Newton's method really wants the latter // // Jacobian row index const int II = ps.gpn_of_patch_irho_isigma(xp, x_irho, x_isigma); // Jacobian coefficients for this point const fp Jacobian_coeff_rho = xp.gridfn(nominal_gfns::gfn__partial_H_wrt_partial_d_h_1, x_irho, x_isigma); const fp Jacobian_coeff_sigma = xp.gridfn(nominal_gfns::gfn__partial_H_wrt_partial_d_h_2, x_irho, x_isigma); const fp Jacobian_coeff_rho_rho = xp.gridfn(nominal_gfns::gfn__partial_H_wrt_partial_dd_h_11, x_irho, x_isigma); const fp Jacobian_coeff_rho_sigma = xp.gridfn(nominal_gfns::gfn__partial_H_wrt_partial_dd_h_12, x_irho, x_isigma); const fp Jacobian_coeff_sigma_sigma = xp.gridfn(nominal_gfns::gfn__partial_H_wrt_partial_dd_h_22, x_irho, x_isigma); // partial_rho, partial_rho_rho for (int m_irho = xp.molecule_min_m() ; m_irho <= xp.molecule_max_m() ; ++m_irho) { const int xm_irho = x_irho + m_irho; const fp Jac_rho = Jacobian_coeff_rho * xp.partial_rho_coeff(m_irho); const fp Jac_rho_rho = Jacobian_coeff_rho_rho * xp.partial_rho_rho_coeff(m_irho); const fp Jac_sum = Jac_rho + Jac_rho_rho; if (xp.is_in_nominal_grid(xm_irho, x_isigma)) then Jac(II, xp,xm_irho,x_isigma) += Jac_sum; else add_ghost_zone_Jacobian (ps, xp, xp.minmax_rho_ghost_zone(m_irho < 0), II, xm_irho, x_isigma, Jac_sum, Jac); } // partial_sigma, partial_sigma_sigma for (int m_isigma = xp.molecule_min_m() ; m_isigma <= xp.molecule_max_m() ; ++m_isigma) { const int xm_isigma = x_isigma + m_isigma; const fp Jac_sigma = Jacobian_coeff_sigma * xp.partial_sigma_coeff(m_isigma); const fp Jac_sigma_sigma = Jacobian_coeff_sigma_sigma * xp.partial_sigma_sigma_coeff(m_isigma); const fp Jac_sum = Jac_sigma + Jac_sigma_sigma; if (xp.is_in_nominal_grid(x_irho, xm_isigma)) then Jac(II, xp,x_irho,xm_isigma) += Jac_sum; else add_ghost_zone_Jacobian (ps, xp, xp.minmax_sigma_ghost_zone(m_isigma < 0), II, x_irho, xm_isigma, Jac_sum, Jac); } // partial_rho_sigma for (int m_irho = xp.molecule_min_m() ; m_irho <= xp.molecule_max_m() ; ++m_irho) { for (int m_isigma = xp.molecule_min_m() ; m_isigma <= xp.molecule_max_m() ; ++m_isigma) { const int xm_irho = x_irho + m_irho; const int xm_isigma = x_isigma + m_isigma; const fp Jac_rho_sigma = Jacobian_coeff_rho_sigma * xp.partial_rho_sigma_coeff(m_irho, m_isigma); if (xp.is_in_nominal_grid(xm_irho, xm_isigma)) then Jac(II, xp,xm_irho,xm_isigma) += Jac_rho_sigma; else add_ghost_zone_Jacobian (ps, xp, xp.corner_ghost_zone_containing_point(m_irho < 0, m_isigma < 0, xm_irho, xm_isigma), II, xm_irho, xm_isigma, Jac_rho_sigma, Jac); } } } } } } //****************************************************************************** // // This function adds the ghost-zone Jacobian dependency contributions // for a single ghost-zone point, to a Jacobian matrix. // // Arguments: // ps = The patch system. // xp = The patch containing the center point of the molecule. // xmgz = If the x+m point is in a ghost zone, this must be that ghost zone. // If the x+m point is not in a ghost zone, this argument is ignored. // x_II = The Jacobian row of the x point. // xm_(irho,isigma) = The coordinates (in xp) of the x+m point of the molecule. // mol = The molecule coefficient. // Jac = The Jacobian matrix. // namespace { void add_ghost_zone_Jacobian(const patch_system& ps, const patch& xp, const ghost_zone& xmgz, int x_II, int xm_irho, int xm_isigma, fp mol, Jacobian& Jac) { const patch_edge& xme = xmgz.my_edge(); const int xm_iperp = xme.iperp_of_irho_isigma(xm_irho, xm_isigma); const int xm_ipar = xme. ipar_of_irho_isigma(xm_irho, xm_isigma); // on what other points ym does this molecule point xm depend // via the patch_system::synchronize() operation? const patch& ymp = ps.synchronize_Jacobian_y_patch(xmgz); const patch_edge& yme = ps.synchronize_Jacobian_y_edge (xmgz); const int min_ym_ipar_m = ps.synchronize_Jacobian_min_y_ipar_m(xmgz); const int max_ym_ipar_m = ps.synchronize_Jacobian_max_y_ipar_m(xmgz); const int ym_iperp = ps.synchronize_Jacobian_y_iperp(xmgz, xm_iperp); const int ym_ipar_posn = ps.synchronize_Jacobian_y_ipar_posn (xmgz, xm_iperp, xm_ipar); // add the Jacobian contributions from the ym points for (int ym_ipar_m = min_ym_ipar_m ; ym_ipar_m <= max_ym_ipar_m ; ++ym_ipar_m) { const int ym_ipar = ym_ipar_posn + ym_ipar_m; const int ym_irho = yme. irho_of_iperp_ipar(ym_iperp,ym_ipar); const int ym_isigma = yme.isigma_of_iperp_ipar(ym_iperp,ym_ipar); const int JJ = ps.gpn_of_patch_irho_isigma(ymp, ym_irho, ym_isigma); const fp sync_Jac = ps.synchronize_Jacobian(xmgz, xm_iperp, xm_ipar, ym_ipar_m); Jac(x_II,JJ) += mol*sync_Jac; } } }