CMS 3D CMS Logo

List of all members | 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>

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. More...
 
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. More...
 
 Fourvec_Constrainer (const Fourvec_Constrainer_Args &args)
 Constructor. More...
 
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. More...
 

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. More...
 

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 220 of file Fourvec_Constrainer.h.

Constructor & Destructor Documentation

◆ Fourvec_Constrainer()

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

Constructor.

Parameters
argsParameter settings for this instance.

Definition at line 174 of file Fourvec_Constrainer.cc.

181  : _args(args) {}
const Fourvec_Constrainer_Args _args

Member Function Documentation

◆ add_constraint()

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 183 of file Fourvec_Constrainer.cc.

References _constraints, and alignCSCRings::s.

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

191  {
192  _constraints.push_back(Constraint(s));
193  }
std::vector< Constraint > _constraints

◆ constrain()

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 1612 of file Fourvec_Constrainer.cc.

References _args, _constraints, _mass_constraint, gpuPixelDoublets::cc, hitfit::Fourvec_Constrainer_Args::chisq_constrainer_args(), makeMEIFBenchmarkPlots::ev, hitfit::Fourvec_Constrainer_Args::ignore_met(), visualization-live-secondInstance_cfg::m, dttmaxenums::R, hitfit::Fourvec_Constrainer_Args::use_e(), x, and beamSpotPI::Y.

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

1633  {
1634  adjust_fourvecs(ev, _args.use_e());
1635 
1636  bool use_kt = ev.has_neutrino() || !_args.ignore_met();
1637  int nobjs = ev.nobjs();
1638  int n_measured_vars = nobjs * 3;
1639  int n_unmeasured_vars = 0;
1640 
1641  if (use_kt) {
1642  n_measured_vars += 2;
1643 
1644  if (ev.has_neutrino())
1645  n_unmeasured_vars = 1;
1646  }
1647 
1648  Matrix G_i(n_measured_vars, n_measured_vars);
1649  Diagonal_Matrix Y(n_unmeasured_vars);
1650  Column_Vector x(n_measured_vars);
1651  Column_Vector y(n_unmeasured_vars);
1652  pack_event(ev, _args.use_e(), use_kt, x, y, G_i, Y);
1653 
1654  Column_Vector xm = x;
1655  Column_Vector ym = y;
1656 
1657  // ??? Should allow for using a different underlying fitter here.
1658  Chisq_Constrainer fitter(_args.chisq_constrainer_args());
1659 
1660  Fourvec_Constraint_Calculator cc(ev, _constraints, _args);
1661 
1662  Matrix Q;
1663  Matrix R;
1664  Matrix S;
1665  double chisq = fitter.fit(cc, xm, x, ym, y, G_i, Y, pullx, pully, Q, R, S);
1666 
1667  unpack_event(ev, x, y, _args.use_e(), use_kt);
1668 
1669  calculate_mass(ev, _mass_constraint, _args, chisq, Q, R, S, m, sigm);
1670 
1671  return chisq;
1672  }
uint32_t cc[maxCellsPerHit]
Definition: gpuFishbone.h:49
CLHEP::HepVector Column_Vector
Definition: matutil.h:63
CLHEP::HepMatrix Matrix
Definition: matutil.h:62
const Fourvec_Constrainer_Args _args
std::vector< Constraint > _mass_constraint
CLHEP::HepDiagMatrix Diagonal_Matrix
Definition: matutil.h:64
std::vector< Constraint > _constraints
const Chisq_Constrainer_Args & chisq_constrainer_args() const

◆ mass_constraint()

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 195 of file Fourvec_Constrainer.cc.

References _mass_constraint, cms::cuda::assert(), and alignCSCRings::s.

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

206  {
207  assert(_mass_constraint.empty());
208  _mass_constraint.push_back(Constraint(s));
209  }
assert(be >=bs)
std::vector< Constraint > _mass_constraint

Friends And Related Function Documentation

◆ operator<<

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 220 of file Fourvec_Constrainer.cc.

231  {
232  s << "Constraints: (e_com = " << c._args.e_com() << ") ";
233  if (c._args.use_e())
234  s << "(E)";
235  s << "\n";
236 
237  for (std::vector<Constraint>::size_type i = 0; i < c._constraints.size(); i++)
238  s << " " << c._constraints[i] << "\n";
239 
240  if (!c._mass_constraint.empty()) {
241  s << "Mass constraint:\n";
242  s << c._mass_constraint[0] << "\n";
243  }
244  return s;
245  }
uint16_t size_type
if(threadIdxLocalY==0 &&threadIdxLocalX==0)

Member Data Documentation

◆ _args

const Fourvec_Constrainer_Args hitfit::Fourvec_Constrainer::_args
private

Parameter settings for this instance.

Definition at line 306 of file Fourvec_Constrainer.h.

Referenced by constrain().

◆ _constraints

std::vector<Constraint> hitfit::Fourvec_Constrainer::_constraints
private

The constraints for this problem.

Definition at line 312 of file Fourvec_Constrainer.h.

Referenced by add_constraint(), and constrain().

◆ _mass_constraint

std::vector<Constraint> hitfit::Fourvec_Constrainer::_mass_constraint
private

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

Definition at line 320 of file Fourvec_Constrainer.h.

Referenced by constrain(), and mass_constraint().