CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_3/src/TopQuarkAnalysis/TopHitFit/src/Fourvec_Constrainer.cc

Go to the documentation of this file.
00001 //
00002 // $Id: Fourvec_Constrainer.cc,v 1.1 2011/05/26 09:47:00 mseidel Exp $
00003 //
00004 // File: src/Fourvec_Constrainer.cc
00005 // Purpose: Do a kinematic fit for a set of 4-vectors, given a set
00006 //          of mass constraints.
00007 // Created: Jul, 2000, sss, based on run 1 mass analysis code.
00008 //
00009 // CMSSW File      : src/Fourvec_Constrainer.cc
00010 // Original Author : Scott Stuart Snyder <snyder@bnl.gov> for D0
00011 // Imported to CMSSW by Haryo Sumowidagdo <Suharyo.Sumowidagdo@cern.ch>
00012 //
00013 
00014 
00041 #include "TopQuarkAnalysis/TopHitFit/interface/Fourvec_Constrainer.h"
00042 #include "TopQuarkAnalysis/TopHitFit/interface/Fourvec_Event.h"
00043 #include "TopQuarkAnalysis/TopHitFit/interface/Pair_Table.h"
00044 #include "TopQuarkAnalysis/TopHitFit/interface/Chisq_Constrainer.h"
00045 #include "TopQuarkAnalysis/TopHitFit/interface/matutil.h"
00046 #include "TopQuarkAnalysis/TopHitFit/interface/Defaults.h"
00047 #include <cmath>
00048 #include <iostream>
00049 
00050 
00051 using std::sqrt;
00052 using std::exp;
00053 using std::cos;
00054 using std::sin;
00055 using std::ostream;
00056 using std::vector;
00057 
00058 
00059 namespace hitfit {
00060 
00061 
00062 //*************************************************************************
00063 // Argument handling.
00064 //
00065 
00066 
00067 Fourvec_Constrainer_Args::Fourvec_Constrainer_Args (const Defaults& defs)
00068 //
00069 // Purpose: Constructor.
00070 //
00071 // Inputs:
00072 //   defs -        The Defaults instance from which to initialize.
00073 //
00074   : _use_e (defs.get_bool ("use_e")),
00075     _e_com (defs.get_float ("e_com")),
00076     _ignore_met (defs.get_bool ("ignore_met")),
00077     _chisq_constrainer_args (defs)
00078 {
00079 }
00080 
00081 
00082 bool Fourvec_Constrainer_Args::use_e () const
00083 //
00084 // Purpose: Return the use_e parameter.
00085 //          See the header for documentation.
00086 //
00087 {
00088   return _use_e;
00089 }
00090 
00091 
00092 double Fourvec_Constrainer_Args::e_com () const
00093 //
00094 // Purpose: Return the e_com parameter.
00095 //          See the header for documentation.
00096 //
00097 {
00098   return _e_com;
00099 }
00100 
00101 
00102 bool Fourvec_Constrainer_Args::ignore_met () const
00103 //
00104 // Purpose: Return the ignore_met parameter.
00105 //          See the header for documentation.
00106 //
00107 {
00108   return _ignore_met;
00109 }
00110 
00111 
00112 const Chisq_Constrainer_Args&
00113 Fourvec_Constrainer_Args::chisq_constrainer_args () const
00114 //
00115 // Purpose: Return the contained subobject parameters.
00116 //
00117 {
00118   return _chisq_constrainer_args;
00119 }
00120 
00121 
00122 //*************************************************************************
00123 // Variable layout.
00124 //
00125 // We need to map the quantities we fit onto the vectors of well- and
00126 // poorly-measured quantities.
00127 //
00128 
00129 //
00130 // The well-measured variables consist of three variables for each
00131 // object.  If we are using transverse momentum constraints,
00132 // these fill be followed by the two cartesian components of kt.
00133 //
00134 // Each object is represented by three variables: the momentum (or 1/p
00135 // if the muon flag was set), and the two spherical angles, phi and eta.
00136 // Here is how they're ordered.
00137 //
00141 typedef enum {
00142   p_offs = 0,
00143   phi_offs = 1,
00144   eta_offs = 2
00145 } Offsets;
00146 
00151 typedef enum {
00152   x_offs = 0,
00153   y_offs = 1
00154 } Kt_Offsets;
00155 
00156 //
00157 // If there is a neutrino, then it is at index 1 of the poorly-measured
00158 // set (otherwise, that set is empty).
00159 //
00164 typedef enum {
00165   nu_z = 1
00166 } Unmeasured_Variables;
00167 
00168 
00169 
00170 namespace {
00171 
00172 
00179 int obj_index (int i)
00180 //
00181 // Purpose: Return the starting variable index for object I.
00182 //
00183 // Inputs:
00184 //   i -           The object index.
00185 //
00186 // Returns:
00187 //   The index in the well-measured set of the first variable
00188 //   for object I.
00189 //
00190 {
00191   return i*3 + 1;
00192 }
00193 
00194 
00195 } // unnamed namespace
00196 
00197 
00198 //*************************************************************************
00199 // Object management.
00200 //
00201 
00202 
00203 Fourvec_Constrainer::Fourvec_Constrainer (const Fourvec_Constrainer_Args& args)
00204 //
00205 // Purpose: Constructor.
00206 //
00207 // Inputs:
00208 //   args -        The parameter settings for this instance.
00209 //
00210   : _args (args)
00211 {
00212 }
00213 
00214 
00215 void Fourvec_Constrainer::add_constraint (std::string s)
00216 //
00217 // Purpose: Specify an additional constraint S for the problem.
00218 //          The format for S is described in the header.
00219 //
00220 // Inputs:
00221 //   s -           The constraint to add.
00222 //
00223 {
00224   _constraints.push_back (Constraint (s));
00225 }
00226 
00227 
00228 void Fourvec_Constrainer::mass_constraint (std::string s)
00229 //
00230 // Purpose: Specify the combination of objects that will be returned by
00231 //          constrain() as the mass.  The format of S is the same as for
00232 //          normal constraints.  The LHS specifies the mass to calculate;
00233 //          the RHS should be zero.
00234 //          This should only be called once.
00235 //
00236 // Inputs:
00237 //   s -           The constraint defining the mass.
00238 //
00239 {
00240   assert (_mass_constraint.size() == 0);
00241   _mass_constraint.push_back (Constraint (s));
00242 }
00243 
00244 
00254 std::ostream& operator<< (std::ostream& s, const Fourvec_Constrainer& c)
00255 //
00256 // Purpose: Print the object to S.
00257 //
00258 // Inputs:
00259 //   s -           The stream to which to write.
00260 //   c -           The object to write.
00261 //
00262 // Returns:
00263 //   The stream S.
00264 //
00265 {
00266   s << "Constraints: (e_com = " << c._args.e_com() << ") ";
00267   if (c._args.use_e())
00268     s << "(E)";
00269   s << "\n";
00270 
00271   for (std::vector<Constraint>::size_type i=0; i < c._constraints.size(); i++)
00272     s << "  " << c._constraints[i] << "\n";
00273 
00274   if (c._mass_constraint.size() > 0) {
00275     s << "Mass constraint:\n";
00276     s << c._mass_constraint[0] << "\n";
00277   }
00278   return s;
00279 }
00280 
00281 
00282 //*************************************************************************
00283 // Event packing and unpacking.
00284 //
00285 
00286 
00287 namespace {
00288 
00289 
00300 void adjust_fourvecs (Fourvec_Event& ev,
00301                       bool use_e_flag)
00302 //
00303 // Purpose: For all objects in EV, adjust their 4-momenta
00304 //          to have their requested masses.
00305 //
00306 // Inputs:
00307 //   ev -          The event on which to operate.
00308 //   use_e_flag -  If true, keep E and scale 3-momentum.
00309 //                 If false, keep the 3-momentum and scale E.
00310 //
00311 {
00312   int nobjs = ev.nobjs ();
00313   for (int i=0; i < nobjs; i++) {
00314     const FE_Obj& obj = ev.obj (i);
00315     Fourvec p = obj.p;
00316     if (use_e_flag)
00317       adjust_p_for_mass (p, obj.mass);
00318     else
00319       adjust_e_for_mass (p, obj.mass);
00320     ev.set_obj_p (i, p);
00321   }
00322 }
00323 
00324 
00341 Fourvec get_p_eta_phi_vec (const Column_Vector& c,
00342                            int ndx,
00343                            const FE_Obj& obj,
00344                            bool use_e_flag)
00345 //
00346 // Purpose: Convert object NDX from its representation in the set
00347 //          of well-measured variables C to a 4-vector.
00348 //
00349 // Inputs:
00350 //   c -           The vector of well-measured variables.
00351 //   ndx -         The index of the object in which we're interested.
00352 //   obj -         The object from the Fourvec_Event.
00353 //   use_e_flag -  If true, we're using E as the fit variable, otherwise p.
00354 //
00355 // Returns:
00356 //   The object's 4-momentum.
00357 //
00358 {
00359   // Get the energy and momentum of the object.
00360   double e, p;
00361 
00362   if (use_e_flag) {
00363     // We're using E as a fit variable.  Get it directly.
00364     e = c(ndx + p_offs);
00365 
00366     // Take into account the muon case.
00367     if (obj.muon_p) e = 1/e;
00368 
00369     // Find the momentum given the energy.
00370     if (obj.mass == 0)
00371       p = e;
00372     else {
00373       double xx = e*e - obj.mass*obj.mass;
00374       if (xx >= 0)
00375     p = sqrt (xx);
00376       else
00377     p = 0;
00378     }
00379   }
00380   else {
00381     // We're using P as a fit variable.  Fetch it.
00382     p = c(ndx + p_offs);
00383 
00384     // Take into account the muon case.
00385     if (obj.muon_p) p = 1/p;
00386 
00387     // Find the energy given the momentum.
00388     e = (obj.mass == 0 ? p : sqrt (obj.mass*obj.mass + p*p));
00389   }
00390 
00391   // Get angular variables.
00392   double phi = c(ndx + phi_offs);
00393   double eta = c(ndx + eta_offs);
00394   if (fabs (eta) > 50) {
00395     // Protect against ridiculously large etas
00396     eta = eta > 0 ? 50 : -50;
00397   }
00398   double exp_eta = exp (eta);
00399   double iexp_eta = 1/exp_eta;
00400   double sin_theta = 2 / (exp_eta + iexp_eta);
00401   double cos_theta = (exp_eta - iexp_eta) / (exp_eta + iexp_eta);
00402 
00403   // Form the 4-momentum.
00404   return Fourvec (p * sin_theta * cos (phi),
00405                   p * sin_theta * sin (phi),
00406                   p * cos_theta,
00407                   e);
00408 }
00409 
00410 
00427 void set_p_eta_phi_vec (const FE_Obj& obj,
00428                         Column_Vector& c,
00429                         int ndx,
00430                         bool use_e_flag)
00431 //
00432 // Purpose: Initialize the variables in the well-measured set C describing
00433 //          object NDX from its Fourvec_Event representation OBJ.
00434 //
00435 // Inputs:
00436 //   obj -         The object from the Fourvec_Event.
00437 //   c -           The vector of well-measured variables.
00438 //   ndx -         The index of the object in which we're interested.
00439 //   use_e_flag -  If true, we're using E as the fit variable, otherwise p.
00440 //
00441 //
00442 {
00443   if (use_e_flag)
00444     c(ndx + p_offs) = obj.p.e();
00445   else
00446     c(ndx + p_offs) = obj.p.vect().mag();
00447 
00448   if (obj.muon_p) c(ndx + p_offs) = 1/c(ndx + p_offs);
00449   c(ndx + phi_offs) = obj.p.phi();
00450   c(ndx + eta_offs) = obj.p.pseudoRapidity();
00451 }
00452 
00453 
00463 Matrix error_matrix (double p_sig,
00464                      double phi_sig,
00465                      double eta_sig)
00466 //
00467 // Purpose: Set up the 3x3 error matrix for an object.
00468 //
00469 // Inputs:
00470 //   p_sig -       The momentum uncertainty.
00471 //   phi_sig -     The phi uncertainty.
00472 //   eta_sig -     The eta uncertainty.
00473 //
00474 // Returns:
00475 //   The object's error matrix.
00476 //
00477 {
00478   Matrix err (3, 3, 0);
00479   err(1+  p_offs, 1+  p_offs) = p_sig * p_sig;
00480   err(1+phi_offs, 1+phi_offs) = phi_sig * phi_sig;
00481   err(1+eta_offs, 1+eta_offs) = eta_sig * eta_sig;
00482   return err;
00483 }
00484 
00485 
00511 void pack_event (const Fourvec_Event& ev,
00512                  bool use_e_flag,
00513                  bool use_kt_flag,
00514                  Column_Vector& xm,
00515                  Column_Vector& ym,
00516                  Matrix& G_i,
00517                  Diagonal_Matrix& Y)
00518 //
00519 // Purpose: Take the information from a Fourvec_Event EV and pack
00520 //          it into vectors of well- and poorly-measured variables.
00521 //          Also set up the error matrices.
00522 //
00523 // Inputs:
00524 //   ev -          The event to pack.
00525 //   use_e_flag -  If true, we're using E as the fit variable, otherwise p.
00526 //   use_kt_flag - True if we're to pack kt variables.
00527 //
00528 // Outputs:
00529 //   xm -          Vector of well-measured variables.
00530 //   ym -          Vector of poorly-measured variables.
00531 //   G_i -         Error matrix for well-measured variables.
00532 //   Y -           Inverse error matrix for poorly-measured variables.
00533 //
00534 {
00535   // Number of objects in the event.
00536   int nobjs = ev.nobjs ();
00537 
00538   int n_measured_vars = nobjs * 3;
00539   if (use_kt_flag)
00540     n_measured_vars += 2;
00541 
00542   // Clear the error matrix.
00543   G_i = Matrix (n_measured_vars, n_measured_vars, 0);
00544 
00545   // Loop over objects.
00546   for (int i=0; i<nobjs; i++) {
00547     const FE_Obj& obj = ev.obj (i);
00548     int this_index = obj_index (i);
00549     set_p_eta_phi_vec (obj, xm, this_index, use_e_flag);
00550     G_i.sub (this_index, this_index, error_matrix (obj.p_error,
00551                                                    obj.phi_error,
00552                                                    obj.eta_error));
00553 
00554   }
00555 
00556   if (use_kt_flag) {
00557     // Set up kt.
00558     int kt_ndx = obj_index (nobjs);
00559     xm (kt_ndx+x_offs) = ev.kt().x();
00560     xm (kt_ndx+y_offs) = ev.kt().y();
00561 
00562     // And its error matrix.
00563     G_i(kt_ndx+x_offs, kt_ndx+x_offs) = ev.kt_x_error() * ev.kt_x_error();
00564     G_i(kt_ndx+y_offs, kt_ndx+y_offs) = ev.kt_y_error() * ev.kt_y_error();
00565     G_i(kt_ndx+x_offs, kt_ndx+y_offs) = ev.kt_xy_covar();
00566     G_i(kt_ndx+y_offs, kt_ndx+x_offs) = ev.kt_xy_covar();
00567   }
00568 
00569   // Handle a neutrino.
00570   if (ev.has_neutrino()) {
00571     ym(nu_z) = ev.nu().z();
00572     Y = Diagonal_Matrix (1, 0);
00573   }
00574 }
00575 
00576 
00596 void unpack_event (Fourvec_Event& ev,
00597                    const Column_Vector& x,
00598                    const Column_Vector& y,
00599                    bool use_e_flag,
00600                    bool use_kt_flag)
00601 //
00602 // Purpose: Update the contents of EV from the variable sets X and Y.
00603 //
00604 // Inputs:
00605 //   ev -          The event.
00606 //   x -           Vector of well-measured variables.
00607 //   use_e_flag -  If true, we're using E as the fit variable, otherwise p.
00608 //   use_kt_flag - True if we're too unpack kt variables.
00609 //
00610 // Outputs:
00611 //   ev -          The event after updating.
00612 //
00613 {
00614   // Do all the objects.
00615   Fourvec sum = 0;
00616   for (int j=0; j<ev.nobjs(); j++) {
00617     const FE_Obj& obj = ev.obj (j);
00618     ev.set_obj_p (j, get_p_eta_phi_vec (x, obj_index (j), obj, use_e_flag));
00619     sum += obj.p;
00620   }
00621 
00622 
00623   if (use_kt_flag) {
00624     int kt_ndx = obj_index (ev.nobjs());
00625     Fourvec kt = Fourvec (x(kt_ndx+x_offs), x(kt_ndx+y_offs), 0, 0);
00626     Fourvec nu = kt - sum;
00627     if (ev.has_neutrino()) {
00628       nu.setPz (y(nu_z));
00629       adjust_e_for_mass (nu, 0);
00630       ev.set_nu_p (nu);
00631     }
00632     else {
00633       adjust_e_for_mass (nu, 0);
00634       ev.set_x_p (nu);
00635     }
00636   }
00637 }
00638 
00639 
00640 } // unnamed namespace
00641 
00642 
00643 //*************************************************************************
00644 // Constraint evaluation.
00645 //
00646 
00647 
00648 namespace {
00649 
00650 
00679 double dot_and_gradient (const Fourvec& v1,
00680                          const Fourvec& v2,
00681                          bool use_e_flag,
00682                          double v1_x[3],
00683                          double v2_x[3],
00684                          bool& badflag)
00685 //
00686 // Purpose: Compute the dot product v1.v2 and its gradients wrt
00687 //          p, phi, and theta of each 4-vector.
00688 //
00689 // Inputs:
00690 //   v1 -          The first 4-vector in the dot product.
00691 //   v2 -          The second 4-vector in the dot product.
00692 //   use_e_flag -  If true, we're using E as the fit variable, otherwise p.
00693 //
00694 // Outputs:
00695 //   v1_x -        Gradients of the dot product wrt v1's p, phi, theta.
00696 //   v2_x -        Gradients of the dot product wrt v2's p, phi, theta.
00697 //   badflag -     Set to true for the singular case (vectors vanish).
00698 //
00699 // Returns:
00700 //   The dot product.
00701 //
00702 {
00703   // Calculate the dot product.
00704   double dot = v1 * v2;
00705 
00706   double p1 = v1.vect().mag();
00707   double p2 = v2.vect().mag();
00708   double e1 = v1.e();
00709   double e2 = v2.e();
00710   double pt1 = v1.vect().perp();
00711   double pt2 = v2.vect().perp();
00712 
00713   // Protect against the singular case.
00714   badflag = false;
00715   if (p1 == 0 || p2 == 0 || e1 == 0 || e2 == 0 || pt1 == 0 || pt2 == 0) {
00716     badflag = true;
00717     v1_x[p_offs] = v1_x[phi_offs] = v1_x[eta_offs] = 0;
00718     v2_x[p_offs] = v2_x[phi_offs] = v2_x[eta_offs] = 0;
00719     return false;
00720   }
00721 
00722   // Calculate the gradients.
00723   v1_x[p_offs] = (dot - v1.m2() * e2 / e1) / p1;
00724   v2_x[p_offs] = (dot - v2.m2() * e1 / e2) / p2;
00725 
00726   if (use_e_flag) {
00727     v1_x[p_offs] *= e1 / p1;
00728     v2_x[p_offs] *= e2 / p2;
00729   }
00730 
00731   v1_x[phi_offs] = v1(1)*v2(0) - v1(0)*v2(1);
00732   v2_x[phi_offs] = -v1_x[phi_offs];
00733 
00734   double fac = v1(0)*v2(0) + v1(1)*v2(1);
00735   v1_x[eta_offs] = pt1*v2(2) - v1(2)/pt1 * fac;
00736   v2_x[eta_offs] = pt2*v1(2) - v2(2)/pt2 * fac;
00737 
00738   return dot;
00739 }
00740 
00741 
00763 void addin_obj_gradient (int constraint_no,
00764                          int sign,
00765                          int base_index,
00766                          const double grad[],
00767                          Matrix& Bx)
00768 //
00769 // Purpose: Tally up the dot product gradients for an object
00770 //          into Bx.
00771 //
00772 // Inputs:
00773 //   constraint_no-The number of the constraint.
00774 //   base_index -  The index in the well-measured variable list
00775 //                 of the first variable for this object.
00776 //   sign -        The sign with which these gradients should be
00777 //                 added into Bx, either +1 or -1.  (I.e., which
00778 //                 side of the constraint equation.)
00779 //   grad -        The gradients for this object, vs. p, phi, theta.
00780 //   Bx -          The well-measured variable gradients.
00781 //
00782 // Outputs:
00783 //   Bx -          The well-measured variable gradients, updated.
00784 //
00785 {
00786   Bx(base_index + p_offs,   constraint_no) += sign * grad[p_offs];
00787   Bx(base_index + phi_offs, constraint_no) += sign * grad[phi_offs];
00788   Bx(base_index + eta_offs, constraint_no) += sign * grad[eta_offs];
00789 }
00790 
00791 
00819 void addin_nu_gradient (int constraint_no,
00820                         int sign,
00821                         int kt_ndx,
00822                         const double grad[],
00823                         Matrix& Bx, Matrix& By)
00824 //
00825 // Purpose: Tally up the dot product gradients for a neutrino
00826 //          into Bx and By.
00827 //
00828 // Inputs:
00829 //   constraint_no-The number of the constraint.
00830 //   sign -        The sign with which these gradients should be
00831 //                 added into Bx, either +1 or -1.  (I.e., which
00832 //                 side of the constraint equation.)
00833 //   kt_ndx -      The index of the kt variables in the variables array.
00834 //   grad -        The gradients for this object, vs. p, phi, theta.
00835 //   Bx -          The well-measured variable gradients.
00836 //   By -          The poorly-measured variable gradients.
00837 //
00838 // Outputs:
00839 //   Bx -          The well-measured variable gradients, updated.
00840 //   By -          The poorly-measured variable gradients, updated.
00841 //
00842 {
00843   Bx(kt_ndx+x_offs,constraint_no) += sign*grad[p_offs];  // Really p for now.
00844   Bx(kt_ndx+y_offs,constraint_no) += sign*grad[phi_offs];// Really phi for now.
00845   By(nu_z,         constraint_no) += sign*grad[eta_offs]; // Really theta ...
00846 }
00847 
00848 
00875 void addin_gradient (const Fourvec_Event& ev,
00876                      int constraint_no,
00877                      int sign,
00878                      int obj_no,
00879                      const double grad[],
00880                      Matrix& Bx,
00881                      Matrix& By)
00882 //
00883 // Purpose: Tally up the dot product gradients for an object (which may
00884 //          or may not be a neutrino) into Bx and By.
00885 //
00886 // Inputs:
00887 //   ev -          The event we're fitting.
00888 //   constraint_no-The number of the constraint.
00889 //   sign -        The sign with which these gradients should be
00890 //                 added into Bx, either +1 or -1.  (I.e., which
00891 //                 side of the constraint equation.)
00892 //   obj_no -      The number of the object.
00893 //   grad -        The gradients for this object, vs. p, phi, theta.
00894 //   Bx -          The well-measured variable gradients.
00895 //   By -          The poorly-measured variable gradients.
00896 //
00897 // Outputs:
00898 //   Bx -          The well-measured variable gradients, updated.
00899 //   By -          The poorly-measured variable gradients, updated.
00900 //
00901 {
00902   if (obj_no >= ev.nobjs()) {
00903     assert (obj_no == ev.nobjs());
00904     addin_nu_gradient (constraint_no, sign, obj_index (obj_no), grad, Bx, By);
00905   }
00906   else
00907     addin_obj_gradient (constraint_no, sign, obj_index (obj_no), grad, Bx);
00908 }
00909 
00910 
00941 void addin_gradients (const Fourvec_Event& ev,
00942                       int constraint_no,
00943                       int sign,
00944                       int i,
00945                       const double igrad[],
00946                       int j,
00947                       const double jgrad[],
00948                       Matrix& Bx,
00949                       Matrix& By)
00950 //
00951 // Purpose: Tally up the gradients from a single dot product into Bx and By.
00952 //
00953 // Inputs:
00954 //   ev -          The event we're fitting.
00955 //   constraint_no-The number of the constraint.
00956 //   sign -        The sign with which these gradients should be
00957 //                 added into Bx, either +1 or -1.  (I.e., which
00958 //                 side of the constraint equation.)
00959 //   i -           The number of the first object.
00960 //   igrad -       The gradients for the first object, vs. p, phi, theta.
00961 //   j -           The number of the second object.
00962 //   jgrad -       The gradients for the second object, vs. p, phi, theta.
00963 //   Bx -          The well-measured variable gradients.
00964 //   By -          The poorly-measured variable gradients.
00965 //
00966 // Outputs:
00967 //   Bx -          The well-measured variable gradients, updated.
00968 //   By -          The poorly-measured variable gradients, updated.
00969 //
00970 {
00971   addin_gradient (ev, constraint_no, sign, i, igrad, Bx, By);
00972   addin_gradient (ev, constraint_no, sign, j, jgrad, Bx, By);
00973 }
00974 
00975 
00988 void add_mass_terms (const Fourvec_Event& ev,
00989                      const vector<Constraint>& constraints,
00990                      Row_Vector& F)
00991 //
00992 // Purpose: Add the m^2 terms into the constraint values.
00993 //
00994 // Inputs:
00995 //   ev -          The event we're fitting.
00996 //   constraints - The list of constraints.
00997 //   F -           Vector of constraint values.
00998 //
00999 // Outputs:
01000 //   F -           Vector of constraint values, updated.
01001 //
01002 {
01003   for (std::vector<Constraint>::size_type i=0; i<constraints.size(); i++)
01004     F(i+1) += constraints[i].sum_mass_terms (ev);
01005 }
01006 
01007 
01041 bool calculate_mass_constraints (const Fourvec_Event& ev,
01042                                  const Pair_Table& pt,
01043                                  const vector<Constraint>& constraints,
01044                                  bool use_e_flag,
01045                                  Row_Vector& F,
01046                                  Matrix& Bx,
01047                                  Matrix& By)
01048 //
01049 // Purpose: Calculate the mass constraints and gradients.
01050 //          Note: At this stage, the gradients are calculated not
01051 //                quite with respect to the fit variables; instead, for
01052 //                all objects (including the neutrino) we calculate
01053 //                the gradients with respect to p, phi, theta.  They'll
01054 //                be converted via appropriate Jacobian transformations
01055 //                later.
01056 //
01057 // Inputs:
01058 //   ev -          The event we're fitting.
01059 //   pt -          Table of cached pair assignments.
01060 //   constraints - The list of constraints.
01061 //   use_e_flag -  If true, we're using E as the fit variable, otherwise p.
01062 //
01063 // Outputs:
01064 //   (nb. these should be passed in correctly dimensioned.)
01065 //   F -           Vector of constraint values.
01066 //   Bx -          The well-measured variable gradients.
01067 //   By -          The poorly-measured variable gradients.
01068 //
01069 {
01070   int npairs = pt.npairs ();
01071   for (int p=0; p < npairs; p++) {
01072     const Objpair& objpair = pt.get_pair (p);
01073     int i = objpair.i();
01074     int j = objpair.j();
01075     double igrad[3], jgrad[3];
01076     bool badflag = false;
01077     double dot = dot_and_gradient (ev.obj (i).p,
01078                                    ev.obj (j).p,
01079                                    use_e_flag,
01080                                    igrad,
01081                                    jgrad,
01082                                    badflag);
01083     if (badflag)
01084       return false;
01085 
01086     for (std::vector<Constraint>::size_type k=0; k < constraints.size(); k++)
01087       if (objpair.for_constraint (k)) {
01088     F(k+1) += objpair.for_constraint (k) * dot;
01089     addin_gradients (ev, k+1, objpair.for_constraint (k),
01090              i, igrad, j, jgrad, Bx, By);
01091       }
01092 
01093   }
01094 
01095   add_mass_terms (ev, constraints, F);
01096   return true;
01097 }
01098 
01099 
01122 void add_nuterm (unsigned ndx,
01123                  const Fourvec& v,
01124                  Matrix& Bx,
01125                  bool use_e_flag,
01126                  int kt_ndx)
01127 //
01128 // Purpose: Carry out the Jacobian transformation for
01129 //          (p_nu^x,p_nu_y) -> (kt^x, kt_y) for a given object.
01130 //
01131 // Inputs:
01132 //   ndx -         Index of the object for which to transform gradients.
01133 //   v -           The object's 4-momentum.
01134 //   Bx -          The well-measured variable gradients.
01135 //   use_e_flag -  If true, we're using E as the fit variable, otherwise p.
01136 //   kt_ndx -      The index of the kt variables in the variables array.
01137 //
01138 // Outputs:
01139 //   Bx -          The well-measured variable gradients, updated.
01140 //
01141 {
01142   double px = v.px();
01143   double py = v.py();
01144   double cot_theta = v.pz() / v.vect().perp();
01145 
01146   for (int j=1; j<=Bx.num_col(); j++) {
01147     double dxnu = Bx(kt_ndx+x_offs, j);
01148     double dynu = Bx(kt_ndx+y_offs, j);
01149 
01150     if (dxnu != 0 || dynu != 0) {
01151       double fac = 1 / v.vect().mag();
01152       if (use_e_flag)
01153     fac = v.e() * fac * fac;
01154       Bx(ndx +   p_offs, j) -= (px*dxnu + py*dynu) * fac;
01155       Bx(ndx + phi_offs, j) +=  py*dxnu - px*dynu;
01156       Bx(ndx + eta_offs, j) -= (px*dxnu + py*dynu) * cot_theta;
01157     }
01158   }
01159 }
01160 
01161 
01185 void convert_neutrino (const Fourvec_Event& ev,
01186                        bool use_e_flag,
01187                        Matrix& Bx,
01188                        Matrix& By)
01189 //
01190 // Purpose: Carry out the Jacobian transformations the neutrino.
01191 //          First, convert from spherical (p, phi, theta) coordinates
01192 //          to rectangular (x, y, z).  Then convert from neutrino pt
01193 //          components to kt components.
01194 //
01195 // Inputs:
01196 //   ev -          The event we're fitting.
01197 //   use_e_flag -  If true, we're using E as the fit variable, otherwise p.
01198 //   Bx -          The well-measured variable gradients.
01199 //   By -          The poorly-measured variable gradients.
01200 //
01201 // Outputs:
01202 //   Bx -          The well-measured variable gradients, updated.
01203 //   By -          The poorly-measured variable gradients, updated.
01204 //
01205 {
01206   int nconstraints = Bx.num_col ();
01207 
01208   const Fourvec& nu = ev.nu ();
01209 
01210   // convert neutrino from polar coordinates to rectangular.
01211   double pnu2  = nu.vect().mag2();  double pnu = sqrt (pnu2);
01212   double ptnu2 = nu.vect().perp2(); double ptnu = sqrt (ptnu2);
01213 
01214   // Doesn't matter whether we use E or P here, since nu is massless.
01215 
01216   double thfac = nu.z()/pnu2/ptnu;
01217   double fac[3][3];
01218   fac[0][0] = nu(0)/pnu;     fac[0][1] = nu(1)/pnu;     fac[0][2] = nu(2)/pnu;
01219   fac[1][0] = - nu(1)/ptnu2; fac[1][1] =   nu(0)/ptnu2; fac[1][2] = 0;
01220   fac[2][0] = nu(0)*thfac;   fac[2][1] = nu(1)*thfac;   fac[2][2] = -ptnu/pnu2;
01221 
01222   int kt_ndx = obj_index (ev.nobjs());
01223   for (int j=1; j<=nconstraints; j++) {
01224     double tmp1 = fac[0][0]*Bx(kt_ndx+x_offs,j) +
01225                   fac[1][0]*Bx(kt_ndx+y_offs,j) +
01226                   fac[2][0]*By(nu_z,j);
01227     Bx(kt_ndx+y_offs,j) = fac[0][1]*Bx(kt_ndx+x_offs,j) +
01228                           fac[1][1]*Bx(kt_ndx+y_offs,j) +
01229                           fac[2][1]*By(nu_z,j);
01230     By(nu_z,j) = fac[0][2]*Bx(kt_ndx+x_offs,j) + fac[2][2]*By(nu_z,j);
01231 
01232     Bx(kt_ndx+x_offs,j) = tmp1;
01233   }
01234 
01235   // Add nu terms.
01236   for (int j=0; j<ev.nobjs(); j++) {
01237     add_nuterm (obj_index (j), ev.obj(j).p, Bx, use_e_flag, kt_ndx);
01238   }
01239 }
01240 
01241 
01262 void calculate_kt_constraints (const Fourvec_Event& ev,
01263                                bool use_e_flag,
01264                                Row_Vector& F,
01265                                Matrix& Bx)
01266 //
01267 // Purpose: Calculate the overall kt constraints (kt = 0) for the case
01268 //          where there is no neutrino.
01269 //
01270 // Inputs:
01271 //   ev -          The event we're fitting.
01272 //   use_e_flag -  If true, we're using E as the fit variable, otherwise p.
01273 //   F -           Vector of constraint values.
01274 //   Bx -          The well-measured variable gradients.
01275 //
01276 // Outputs:
01277 //   F -           Vector of constraint values, updated.
01278 //   Bx -          The well-measured variable gradients, updated.
01279 //
01280 {
01281   Fourvec tmp = 0;
01282   int base = F.num_col() - 2;
01283   int nobjs = ev.nobjs();
01284   for (int j=0; j<nobjs; j++) {
01285     const Fourvec& obj = ev.obj (j).p;
01286     tmp += obj;
01287 
01288     int ndx = obj_index (j);
01289     double p = obj.vect().mag();
01290     double cot_theta = obj.z() / obj.vect().perp();
01291 
01292     Bx(ndx +   p_offs, base+1) = obj(0) / p;
01293     Bx(ndx + phi_offs, base+1) = -obj(1);
01294     Bx(ndx + eta_offs, base+1) = obj(0) * cot_theta;
01295 
01296     Bx(ndx +   p_offs, base+2) = obj(1) / p;
01297     Bx(ndx + phi_offs, base+2) = obj(0);
01298     Bx(ndx + eta_offs, base+2) = obj(1) * cot_theta;
01299 
01300     if (use_e_flag) {
01301       Bx(ndx +   p_offs, base+1) *= obj.e() / p;
01302       Bx(ndx +   p_offs, base+2) *= obj.e() / p;
01303     }
01304   }
01305 
01306   int kt_ndx = obj_index (nobjs);
01307   Bx(kt_ndx+x_offs, base+1) = -1;
01308   Bx(kt_ndx+y_offs, base+2) = -1;
01309 
01310   F(base+1) = tmp(0) - ev.kt().x ();
01311   F(base+2) = tmp(1) - ev.kt().y ();
01312 }
01313 
01314 
01328 void ddtheta_to_ddeta (int i, double cos_theta, Matrix& Bx)
01329 //
01330 // Purpose: Do the Jacobian transformation from theta -> eta
01331 //          for a single object.
01332 //
01333 // Inputs:
01334 //   i -           The index of the object.
01335 //   cos_theta -   cos(theta) for the object.
01336 //   Bx -          The well-measured variable gradients.
01337 //
01338 // Outputs:
01339 //   Bx -          The well-measured variable gradients, updated.
01340 //
01341 {
01342   double sin_theta = sqrt (1 - cos_theta * cos_theta);
01343   for (int j=1; j<=Bx.num_col(); j++)
01344     Bx(i,j) *= - sin_theta;   /* \sin\theta = 1 / \cosh\eta */
01345 }
01346 
01347 
01348 } // unnamed namespace
01349 
01350 
01376 class Fourvec_Constraint_Calculator
01377   : public Constraint_Calculator
01378 //
01379 // Purpose: Constraint evaluator.
01380 //
01381 {
01382 public:
01383   // Constructor, destructor.
01384   Fourvec_Constraint_Calculator (Fourvec_Event& ev,
01385                                  const vector<Constraint>& constraints,
01386                                  const Fourvec_Constrainer_Args& args);
01387   virtual ~Fourvec_Constraint_Calculator () {}
01388 
01389   // Evaluate constraints at the point described by X and Y (well-measured
01390   // and poorly-measured variables, respectively).  The results should
01391   // be stored in F.  BX and BY should be set to the gradients of F with
01392   // respect to X and Y, respectively.
01393   //
01394   // Return true if the point X, Y is accepted.
01395   // Return false if it is rejected (i.e., in an unphysical region).
01396   // The constraints need not be evaluated in that case.
01397   virtual bool eval (const Column_Vector& x,
01398                      const Column_Vector& y,
01399                      Row_Vector& F,
01400                      Matrix& Bx,
01401                      Matrix& By);
01402 
01403 
01404   // Calculate the constraint functions and gradients.
01405   bool calculate_constraints (Row_Vector& F,
01406                               Matrix& Bx,
01407                               Matrix& By) const;
01408 
01409 private:
01410   // The event we're fitting.
01411   Fourvec_Event& _ev;
01412 
01413   // Vector of constraints.
01414   const vector<Constraint>& _constraints;
01415 
01416   // Argument values.
01417   const Fourvec_Constrainer_Args& _args;
01418 
01419   // The pair table.
01420   Pair_Table _pt;
01421 };
01422 
01423 
01433 Fourvec_Constraint_Calculator::Fourvec_Constraint_Calculator
01434   (Fourvec_Event& ev,
01435    const vector<Constraint>& constraints,
01436    const Fourvec_Constrainer_Args& args)
01437 //
01438 // Purpose: Constructor.
01439 //
01440 // Inputs:
01441 //   ev -          The event we're fitting.
01442 //   constraints - The list of constraints.
01443 //   args -        The parameter settings.
01444 //
01445 //
01446   : Constraint_Calculator (constraints.size() +
01447                            ((ev.has_neutrino() || args.ignore_met()) ? 0 : 2)),
01448     _ev (ev),
01449     _constraints (constraints),
01450     _args (args),
01451     _pt (constraints, ev)
01452 
01453 {
01454 }
01455 
01456 
01471 bool
01472 Fourvec_Constraint_Calculator::calculate_constraints (Row_Vector& F,
01473                                                       Matrix& Bx,
01474                                                       Matrix& By) const
01475 //
01476 // Purpose: Calculate the constraint functions and gradients.
01477 //
01478 // Outputs:
01479 //   F -           Vector of constraint values.
01480 //   Bx -          The well-measured variable gradients.
01481 //   By -          The poorly-measured variable gradients.
01482 //
01483 {
01484   // Clear the matrices.
01485   Bx = Matrix (Bx.num_row(), Bx.num_col(), 0);
01486   By = Matrix (By.num_row(), By.num_col(), 0);
01487   F = Row_Vector (F.num_col(), 0);
01488 
01489   const double p_eps = 1e-10;
01490 
01491   if (_ev.has_neutrino() && _ev.nu().z() > _args.e_com()) {
01492     return false;
01493   }
01494 
01495   int nobjs = _ev.nobjs ();
01496 
01497   // Reject the point if any of the momenta get too small.
01498   for (int j=0; j<nobjs; j++) {
01499     if (_ev.obj(j).p.perp() <= p_eps || _ev.obj(j).p.e() <= p_eps) {
01500       return false;
01501     }
01502   }
01503 
01504   if (! calculate_mass_constraints (_ev, _pt, _constraints, _args.use_e(),
01505                                     F, Bx, By))
01506     return false;
01507 
01508   if (_ev.has_neutrino())
01509     convert_neutrino (_ev, _args.use_e(), Bx, By);
01510   else if (!_args.ignore_met())
01511   {
01512     /* kt constraints */
01513     calculate_kt_constraints (_ev, _args.use_e(), F, Bx);
01514   }
01515 
01516   /* convert d/dtheta to d/deta */
01517   for (int j=0; j<nobjs; j++) {
01518     ddtheta_to_ddeta (obj_index (j) + eta_offs,
01519                       _ev.obj(j).p.cosTheta(),
01520                       Bx);
01521   }
01522 
01523   /* handle muons */
01524   for (int j=0; j<nobjs; j++) {
01525     const FE_Obj& obj = _ev.obj (j);
01526     if (obj.muon_p) {
01527       // Again, E vs. P doesn't matter here since we assume muons to be massless.
01528       double pmu2 = obj.p.vect().mag2();
01529       int ndx = obj_index (j) + p_offs;
01530       for (int k=1; k<=Bx.num_col(); k++)
01531     Bx(ndx, k) = - Bx(ndx, k) * pmu2;
01532     }
01533   }
01534 
01535   return true;
01536 }
01537 
01538 
01568 bool Fourvec_Constraint_Calculator::eval (const Column_Vector& x,
01569                                           const Column_Vector& y,
01570                                           Row_Vector& F,
01571                                           Matrix& Bx,
01572                                           Matrix& By)
01573 //
01574 // Purpose: Evaluate constraints at the point described by X and Y
01575 //          (well-measured and poorly-measured variables, respectively).
01576 //          The results should be stored in F.  BX and BY should be set
01577 //          to the gradients of F with respect to X and Y, respectively.
01578 //
01579 // Inputs:
01580 //   x -           Vector of well-measured variables.
01581 //   y -           Vector of poorly-measured variables.
01582 //
01583 // Outputs:
01584 //   F -           Vector of constraint values.
01585 //   Bx -          The well-measured variable gradients.
01586 //   By -          The poorly-measured variable gradients.
01587 //
01588 {
01589   int nobjs = _ev.nobjs();
01590 
01591   const double p_eps = 1e-10;
01592   const double eta_max = 10;
01593 
01594   // Give up if we've gone into an obviously unphysical region.
01595   for (int j=0; j<nobjs; j++)
01596     if (x(obj_index (j) + p_offs) < p_eps ||
01597         fabs(x(obj_index (j) + eta_offs)) > eta_max) {
01598       return false;
01599     }
01600 
01601   unpack_event (_ev, x, y, _args.use_e(),
01602                 _ev.has_neutrino() || !_args.ignore_met());
01603 
01604   return calculate_constraints (F, Bx, By);
01605 }
01606 
01607 
01608 //*************************************************************************
01609 // Mass calculation.
01610 //
01611 
01612 
01613 namespace {
01614 
01615 
01633 double calculate_sigm (const Matrix& Q,
01634                        const Matrix& R,
01635                        const Matrix& S,
01636                        const Column_Vector& Bx,
01637                        const Column_Vector& By)
01638 //
01639 // Purpose: Do error propagation to find the uncertainty in the final mass.
01640 //
01641 // Inputs:
01642 //   Q -           The final error matrix for the well-measured variables.
01643 //   R -           The final error matrix for the poorly-measured variables.
01644 //   S -           The final cross error matrix for the two sets of variables.
01645 //   Bx -          The well-measured variable gradients.
01646 //   By -          The poorly-measured variable gradients.
01647 //
01648 {
01649   double sig2 = scalar (Bx.T() * Q * Bx);
01650 
01651   if (By.num_row() > 0) {
01652     sig2 += scalar (By.T() * R * By);
01653     sig2 += 2 * scalar (Bx.T() * S * By);
01654   }
01655 
01656   assert (sig2 >= 0);
01657   return sqrt (sig2);
01658 }
01659 
01660 
01683 void calculate_mass (Fourvec_Event& ev,
01684                      const vector<Constraint>& mass_constraint,
01685                      const Fourvec_Constrainer_Args& args,
01686                      double chisq,
01687                      const Matrix& Q,
01688                      const Matrix& R,
01689                      const Matrix& S,
01690                      double& m,
01691                      double& sigm)
01692 //
01693 // Purpose: Calculate the final requested mass and its uncertainty.
01694 //
01695 // Inputs:
01696 //   ev -          The event we're fitting.
01697 //   mass_constraint- The description of the mass we're to calculate.
01698 //   args -        Parameter settings.
01699 //   chisq -       The chisq from the fit.
01700 //   Q -           The final error matrix for the well-measured variables.
01701 //   R -           The final error matrix for the poorly-measured variables.
01702 //   S -           The final cross error matrix for the two sets of variables.
01703 //
01704 // Outputs:
01705 //   m -           The mass.
01706 //   sigm -        Its uncertainty.
01707 //
01708 {
01709   // Don't do anything if the mass wasn't specified.
01710   if (mass_constraint.size () == 0) {
01711     m = 0;
01712     sigm = 0;
01713     return;
01714   }
01715 
01716   // Do the constraint calculation.
01717   int n_measured_vars = ev.nobjs()*3 + 2;
01718   int n_unmeasured_vars = 0;
01719   if (ev.has_neutrino ()) n_unmeasured_vars = 1;
01720 
01721   Row_Vector F(1);
01722   Matrix Bx (n_measured_vars, 1);
01723   Matrix By (n_unmeasured_vars, 1);
01724 
01725   Fourvec_Constraint_Calculator cc (ev, mass_constraint, args);
01726   cc.calculate_constraints (F, Bx, By);
01727 
01728   // Calculate the mass.
01729   //assert (F(1) >= 0);
01730   if (F(1) >= 0.)
01731     m = sqrt (F(1) * 2);
01732   else {
01733     m = 0.;
01734     chisq = -100.;
01735   }
01736 
01737   // And the uncertainty.
01738   // We can only do this if the fit converged.
01739   if (chisq < 0)
01740     sigm = 0;
01741   else {
01742     //assert (F(1) > 0);
01743     Bx = Bx / (m);
01744     By = By / (m);
01745 
01746     sigm = calculate_sigm (Q, R, S, Bx, By);
01747   }
01748 }
01749 
01750 
01751 } // unnamed namespace
01752 
01753 
01754 double Fourvec_Constrainer::constrain (Fourvec_Event& ev,
01755                                        double& m,
01756                                        double& sigm,
01757                                        Column_Vector& pullx,
01758                                        Column_Vector& pully)
01759 //
01760 // Purpose: Do a constrained fit for EV.  Returns the requested mass and
01761 //          its error in M and SIGM, and the pull quantities in PULLX and
01762 //          PULLY.  Returns the chisq; this will be < 0 if the fit failed
01763 //          to converge.
01764 //
01765 // Inputs:
01766 //   ev -          The event we're fitting.
01767 //
01768 // Outputs:
01769 //   ev -          The fitted event.
01770 //   m -           Requested invariant mass.
01771 //   sigm -        Uncertainty on m.
01772 //   pullx -       Pull quantities for well-measured variables.
01773 //   pully -       Pull quantities for poorly-measured variables.
01774 //
01775 // Returns:
01776 //   The fit chisq, or < 0 if the fit didn't converge.
01777 //
01778 {
01779   adjust_fourvecs (ev, _args.use_e ());
01780 
01781   bool use_kt = ev.has_neutrino() || !_args.ignore_met();
01782   int nobjs = ev.nobjs ();
01783   int n_measured_vars = nobjs * 3;
01784   int n_unmeasured_vars = 0;
01785   int n_constraints = _constraints.size ();
01786 
01787   if (use_kt) {
01788     n_measured_vars += 2;
01789 
01790     if (ev.has_neutrino ())
01791       n_unmeasured_vars = 1;
01792     else
01793       n_constraints += 2;
01794   }
01795 
01796   Matrix G_i (n_measured_vars, n_measured_vars);
01797   Diagonal_Matrix Y (n_unmeasured_vars);
01798   Column_Vector x (n_measured_vars);
01799   Column_Vector y (n_unmeasured_vars);
01800   pack_event (ev, _args.use_e(), use_kt, x, y, G_i, Y);
01801 
01802   Column_Vector xm = x;
01803   Column_Vector ym = y;
01804 
01805   // ??? Should allow for using a different underlying fitter here.
01806   Chisq_Constrainer fitter (_args.chisq_constrainer_args());
01807 
01808   Fourvec_Constraint_Calculator cc (ev, _constraints, _args);
01809 
01810   Matrix Q;
01811   Matrix R;
01812   Matrix S;
01813   double chisq = fitter.fit (cc,
01814                              xm, x, ym, y, G_i, Y,
01815                              pullx, pully,
01816                              Q, R, S);
01817 
01818   unpack_event (ev, x, y, _args.use_e (), use_kt);
01819 
01820   calculate_mass (ev, _mass_constraint, _args, chisq, Q, R, S, m, sigm);
01821 
01822   return chisq;
01823 }
01824 
01825 
01826 } // namespace hitfit