CMS 3D CMS Logo

Public Member Functions | Private Attributes | Friends

hitfit::Fourvec_Constrainer Class Reference

Do a kinematic fit for a set of four-momenta, given a set of mass constraints. More...

#include <Fourvec_Constrainer.h>

List of all members.

Public Member Functions

void add_constraint (std::string s)
 Specify an additional constraint s for the problem. The format for s is described in the class description.
double constrain (Fourvec_Event &ev, double &m, double &sigm, Column_Vector &pullx, Column_Vector &pully)
 Do a constrained fit for event ev. Returns the requested mass and its uncertainty in m and sigm, and the pull quantities in pullx and pully. Returns the $\chi^{2}$, the value will be negative if the fit failed to converge.
 Fourvec_Constrainer (const Fourvec_Constrainer_Args &args)
 Constructor.
void mass_constraint (std::string s)
 Specify a combination of objects that will be returned by the constrain() method as mass. The format of s is the same as for normal constraints. The left-hand side specifies the mass to calculate, the right-hand side should be zero. This combination of objects will be called only once.

Private Attributes

const Fourvec_Constrainer_Args _args
std::vector< Constraint_constraints
std::vector< Constraint_mass_constraint

Friends

std::ostream & operator<< (std::ostream &s, const Fourvec_Constrainer &c)
 Output stream operator, print the content of this Fourvec_Constrainer to an output stream.

Detailed Description

Do a kinematic fit for a set of four-momenta, given a set of mass constraints.

The input is in a Fourvec_Event instance. This consists of a collection of objects, each of which has a four-momentum, an uncertainty on the four-momentum, and an integer label.

There are also a set of mass constraints, based on sums of objects with specified labels. A constraint can either require that the invariant mass of a set of objects is constant, or that the mass of two sets of objects be equal to each other. For example, the constraint

\[ (1~~2) = 80 \]

says that the sum of all objects with labels 1 and 2 should have a mass of 80, where the unit of mass is always in GeV. And the constraint

\[ (1~~2) = (3~~4) \]

says that the sum of all objects with labels 1 and 2 should have an invariant mass equal as the invariant mass of the sum of all objects with label 3 and 4.

All the objects are fixed to constant masses for the fit. These masses are attributes of the objects in the Fourvec_Event. This is done by scaling either the four-momentum's three-momentum or energy, depending on the setting of the use_e parameter.

If there is no neutrino present in the event, two additional constraints are automatically added for total $p_{T}$ balance, unless the parameter ignore_met has been set to TRUE.

When the fit completes, this object can compute an invariant mass out of some combination of the objects, including and uncertainty. The particular combination is specified through the method mass_constraint((); it gets a constraint string like normal with the LHS of the constraint i s the mass which will be calculater, and the RHS of the oncstrant should be zero.

Definition at line 230 of file Fourvec_Constrainer.h.


Constructor & Destructor Documentation

hitfit::Fourvec_Constrainer::Fourvec_Constrainer ( const Fourvec_Constrainer_Args args)

Constructor.

Parameters:
argsParameter settings for this instance.

Definition at line 203 of file Fourvec_Constrainer.cc.

  : _args (args)
{
}

Member Function Documentation

void hitfit::Fourvec_Constrainer::add_constraint ( std::string  s)

Specify an additional constraint s for the problem. The format for s is described in the class description.

Parameters:
sThe constraint to add.

Definition at line 215 of file Fourvec_Constrainer.cc.

References _constraints.

Referenced by hitfit::Constrained_Top::Constrained_Top(), and hitfit::Constrained_Z::Constrained_Z().

{
  _constraints.push_back (Constraint (s));
}
double hitfit::Fourvec_Constrainer::constrain ( Fourvec_Event ev,
double &  m,
double &  sigm,
Column_Vector pullx,
Column_Vector pully 
)

Do a constrained fit for event ev. Returns the requested mass and its uncertainty in m and sigm, and the pull quantities in pullx and pully. Returns the $\chi^{2}$, the value will be negative if the fit failed to converge.

Parameters:
evThe event to be fitted (input), and later after fitting (output).
mThe requested invariant mass to calculate.
sigmThe uncertainty in the requested invariant mass.
pullxPull quantities for well-measured variables.
pullyPull quantities for poorly-measured variables.
Input:
  • ev (before fit).
Output:
  • ev (after fit).
  • m.
  • sigm.
  • pullx
  • pully.
Returns:
The fit $\chi^{2}$, this value will be negative if the fit failed to converge.

Definition at line 1754 of file Fourvec_Constrainer.cc.

References _args, _constraints, _mass_constraint, hitfit::Fourvec_Constrainer_Args::chisq_constrainer_args(), hitfit::Fourvec_Event::has_neutrino(), hitfit::Fourvec_Constrainer_Args::ignore_met(), hitfit::Fourvec_Event::nobjs(), dttmaxenums::R, hitfit::Fourvec_Constrainer_Args::use_e(), x, and detailsBasic3DVector::y.

Referenced by hitfit::Constrained_Z::constrain(), and hitfit::Constrained_Top::constrain().

{
  adjust_fourvecs (ev, _args.use_e ());

  bool use_kt = ev.has_neutrino() || !_args.ignore_met();
  int nobjs = ev.nobjs ();
  int n_measured_vars = nobjs * 3;
  int n_unmeasured_vars = 0;
  int n_constraints = _constraints.size ();

  if (use_kt) {
    n_measured_vars += 2;

    if (ev.has_neutrino ())
      n_unmeasured_vars = 1;
    else
      n_constraints += 2;
  }

  Matrix G_i (n_measured_vars, n_measured_vars);
  Diagonal_Matrix Y (n_unmeasured_vars);
  Column_Vector x (n_measured_vars);
  Column_Vector y (n_unmeasured_vars);
  pack_event (ev, _args.use_e(), use_kt, x, y, G_i, Y);

  Column_Vector xm = x;
  Column_Vector ym = y;

  // ??? Should allow for using a different underlying fitter here.
  Chisq_Constrainer fitter (_args.chisq_constrainer_args());

  Fourvec_Constraint_Calculator cc (ev, _constraints, _args);

  Matrix Q;
  Matrix R;
  Matrix S;
  double chisq = fitter.fit (cc,
                             xm, x, ym, y, G_i, Y,
                             pullx, pully,
                             Q, R, S);

  unpack_event (ev, x, y, _args.use_e (), use_kt);

  calculate_mass (ev, _mass_constraint, _args, chisq, Q, R, S, m, sigm);

  return chisq;
}
void hitfit::Fourvec_Constrainer::mass_constraint ( std::string  s)

Specify a combination of objects that will be returned by the constrain() method as mass. The format of s is the same as for normal constraints. The left-hand side specifies the mass to calculate, the right-hand side should be zero. This combination of objects will be called only once.

Parameters:
sThe constraint defining the mass.

Definition at line 228 of file Fourvec_Constrainer.cc.

References _mass_constraint.

Referenced by hitfit::Constrained_Top::Constrained_Top().

{
  assert (_mass_constraint.size() == 0);
  _mass_constraint.push_back (Constraint (s));
}

Friends And Related Function Documentation

std::ostream& operator<< ( std::ostream &  s,
const Fourvec_Constrainer c 
) [friend]

Output stream operator, print the content of this Fourvec_Constrainer to an output stream.

Parameters:
sThe output stream to which to write.
cThe instance of Fourvec_Constrainer to be printed.

Definition at line 254 of file Fourvec_Constrainer.cc.

{
  s << "Constraints: (e_com = " << c._args.e_com() << ") ";
  if (c._args.use_e())
    s << "(E)";
  s << "\n";

  for (std::vector<Constraint>::size_type i=0; i < c._constraints.size(); i++)
    s << "  " << c._constraints[i] << "\n";

  if (c._mass_constraint.size() > 0) {
    s << "Mass constraint:\n";
    s << c._mass_constraint[0] << "\n";
  }
  return s;
}

Member Data Documentation

Parameter settings for this instance.

Definition at line 322 of file Fourvec_Constrainer.h.

Referenced by constrain(), and hitfit::operator<<().

The constraints for this problem.

Definition at line 328 of file Fourvec_Constrainer.h.

Referenced by add_constraint(), constrain(), and hitfit::operator<<().

The constraint giving the mass to be calculated. This should have no more than one entry.

Definition at line 336 of file Fourvec_Constrainer.h.

Referenced by constrain(), mass_constraint(), and hitfit::operator<<().