CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
Constrained_Top.cc
Go to the documentation of this file.
1 //
2 //
3 // File: src/Constrained_Top.cc
4 // Purpose: Do kinematic fitting for a ttbar -> ljets event.
5 // Created: Jul, 2000, sss, based on run 1 mass analysis code.
6 //
7 // CMSSW File : src/Constrained_Top.cc
8 // Original Author : Scott Stuart Snyder <snyder@bnl.gov> for D0
9 // Imported to CMSSW by Haryo Sumowidagdo <Suharyo.Sumowidagdo@cern.ch>
10 //
11 
12 
42 #include <ostream>
43 #include <cassert>
44 #include <stdio.h>
45 
46 
47 using std::ostream;
48 
49 
50 namespace hitfit {
51 
52 
53 //*************************************************************************
54 // Argument handling.
55 //
56 
57 
59 //
60 // Purpose: Constructor.
61 //
62 // Inputs:
63 // defs - The Defaults instance from which to initialize.
64 //
65  : _bmass (defs.get_float ("bmass")),
66  _fourvec_constrainer_args (defs),
67  _equal_side(defs.get_bool("equal_side"))
68 
69 {
70 }
71 
72 
74 //
75 // Purpose: Return the bmass parameter.
76 // See the header for documentation.
77 //
78 {
79  return _bmass;
80 }
81 
82 
85 //
86 // Purpose: Return the contained subobject parameters.
87 //
88 {
90 }
91 
92 
94 //
95 // Purpose: Return the equal_side parameter
96 // See the header for documentation.
97 //
98 {
99  return _equal_side;
100 }
101 
102 
103 //*************************************************************************
104 
105 
107  double lepw_mass,
108  double hadw_mass,
109  double top_mass)
110 //
111 // Purpose: Constructor.
112 //
113 // Inputs:
114 // args - The parameter settings for this instance.
115 // lepw_mass - The mass to which the leptonic W should be constrained,
116 // or 0 to skip this constraint.
117 // hadw_mass - The mass to which the hadronic W should be constrained,
118 // or 0 to skip this constraint.
119 // top_mass - The mass to which the top quarks should be constrained,
120 // or 0 to skip this constraint.
121 //
122  : _args (args),
123  _constrainer (args.fourvec_constrainer_args())
124 {
125  char buf[256];
126  if (lepw_mass > 0) {
127  sprintf (buf, "(%d %d) = %f", nu_label, lepton_label, lepw_mass);
129  }
130 
131  if (hadw_mass > 0) {
132  sprintf (buf, "(%d %d) = %f", hadw1_label, hadw2_label, hadw_mass);
134  }
135 
136  if (args.equal_side()) {
137  sprintf (buf, "(%d %d %d) = (%d %d %d)",
141  }
142 
143  if (top_mass > 0) {
144  sprintf (buf, "(%d %d %d) = %f",
145  hadw1_label, hadw2_label, hadb_label, top_mass);
147  }
148  else {
149  sprintf (buf, "(%d %d %d) = 0", hadw1_label, hadw2_label, hadb_label);
151  }
152 }
153 
154 namespace {
155 
156 
171 FE_Obj make_fe_obj (const Lepjets_Event_Lep& obj, double mass, int type)
172 //
173 // Purpose: Helper to create an object to put into the Fourvec_Event.
174 //
175 // Inputs:
176 // obj - The input object.
177 // mass - The mass to which it should be constrained.
178 // type - The type to assign it.
179 //
180 // Returns:
181 // The constructed FE_Obj.
182 //
183 {
184  return FE_Obj (obj.p(), mass, type,
185  obj.p_sigma(), obj.eta_sigma(), obj.phi_sigma(),
186  obj.res().p_res().inverse());
187 }
188 
189 
208 void do_import (const Lepjets_Event& ev, double bmass, Fourvec_Event& fe)
209 //
210 // Purpose: Convert from a Lepjets_Event to a Fourvec_Event.
211 //
212 // Inputs:
213 // ev - The input event.
214 // bmass - The mass to which b-jets should be fixed.
215 //
216 // Outputs:
217 // fe - The initialized Fourvec_Event.
218 //
219 {
220  assert (ev.nleps() == 1);
221  fe.add (make_fe_obj (ev.lep(0), 0, lepton_label));
222 
223  bool saw_lepb = false;
224  bool saw_hadb = false;
225  for (std::vector<Lepjets_Event_Jet>::size_type j=0; j < ev.njets(); j++) {
226  if (ev.jet(j).type() == isr_label || ev.jet(j).type() == higgs_label)
227  continue;
228  double mass = 0;
229  if (ev.jet(j).type() == lepb_label && !saw_lepb) {
230  mass = bmass;
231  saw_lepb = true;
232  }
233  else if (ev.jet(j).type() == hadb_label && !saw_hadb) {
234  mass = bmass;
235  saw_hadb = true;
236  }
237  fe.add (make_fe_obj (ev.jet(j), mass, ev.jet(j).type()));
238  }
239 
240  fe.set_nu_p (ev.met());
241  Fourvec kt = ev.kt ();
242  fe.set_kt_error (ev.kt_res().sigma (kt.x()),
243  ev.kt_res().sigma (kt.y()),
244  0);
245 }
246 
247 
264 void do_export (const Fourvec_Event& fe, Lepjets_Event& ev)
265 //
266 // Purpose: Convert from a Fourvec_Event to a Lepjets_Event.
267 //
268 // Inputs:
269 // fe - The input event.
270 // ev - The original Lepjets_Event.
271 //
272 // Outputs:
273 // ev - The updated Lepjets_Event.
274 //
275 {
276  ev.lep(0).p() = fe.obj(0).p;
277  for (std::vector<Lepjets_Event_Jet>::size_type j=0, k=1; j < ev.njets(); j++) {
278  if (ev.jet(j).type() == isr_label || ev.jet(j).type() == higgs_label)
279  continue;
280  ev.jet(j).p() = fe.obj(k++).p;
281  }
282 
283  ev.met() = fe.nu();
284 }
285 
286 
287 } // unnamed namespace
288 
289 
291  double& mt,
292  double& sigmt,
293  Column_Vector& pullx,
294  Column_Vector& pully)
295 //
296 // Purpose: Do a constrained fit for EV. Returns the top mass and
297 // its error in MT and SIGMT, and the pull quantities in PULLX and
298 // PULLY. Returns the chisq; this will be < 0 if the fit failed
299 // to converge.
300 //
301 // Inputs:
302 // ev - The event we're fitting.
303 //
304 // Outputs:
305 // ev - The fitted event.
306 // mt - Requested invariant mass.
307 // sigmt - Uncertainty on mt.
308 // pullx - Pull quantities for well-measured variables.
309 // pully - Pull quantities for poorly-measured variables.
310 //
311 // Returns:
312 // The fit chisq, or < 0 if the fit didn't converge.
313 //
314 {
315  Fourvec_Event fe;
316  do_import (ev, _args.bmass (), fe);
317  double chisq = _constrainer.constrain (fe, mt, sigmt, pullx, pully);
318  do_export (fe, ev);
319 
320  return chisq;
321 }
322 
323 
333 std::ostream& operator<< (std::ostream& s, const Constrained_Top& ct)
334 //
335 // Purpose: Print the object to S.
336 //
337 // Inputs:
338 // s - The stream to which to write.
339 // ct - The object to write.
340 //
341 // Returns:
342 // The stream S.
343 //
344 {
345  return s << ct._constrainer;
346 }
347 
348 
349 } // namespace hitfit
type
Definition: HCALResponse.h:21
Hold on to parameters for the Constrained_Top class.
Define an abstract interface for getting parameter settings.
Constrained_Top(const Constrained_Top_Args &args, double lepw_mass, double hadw_mass, double top_mass)
Constructor, create an instance of the Constrained_Top object from the arguments object and the mass ...
Represent a lepton in an instance of Lepjets_Event class. This class hold the following information: ...
const Fourvec_Constrainer_Args & fourvec_constrainer_args() const
Represent an event for kinematic fitting as a collection of four-momenta. Each object is represented ...
Do a constrained kinematic fitting for a event.
bool inverse() const
Return the setting of the inverse flag.
Definition: Resolution.cc:171
CLHEP::HepVector Column_Vector
Definition: matutil.h:66
const int nu_label
Fourvec_Constrainer _constrainer
bool ev
uint16_t size_type
double eta_sigma() const
Return the uncertainty in pseudorapidity .
const Vector_Resolution & res() const
Return a constant reference to the resolution.
double phi_sigma() const
Return the uncertainty in azimuthal angle .
double p_sigma() const
Return the uncertainty in momentum or ( or if the lepton is a tracking object).
Represent a simple event consisting of lepton(s) and jet(s).
Represent a simple event consisting of lepton(s) and jet(s). An instance of this class holds a list o...
Definition: Lepjets_Event.h:66
CLHEP::HepLorentzVector Fourvec
Typedef for a HepLorentzVector.
Definition: fourvec.h:57
int j
Definition: DBlmapReader.cc:9
double constrain(Lepjets_Event &ev, double &mt, double &sigmt, Column_Vector &pullx, Column_Vector &pully)
Do a constrained fit of events. Returns the top mass and its error in mt and sigmt, and the pull quantities in pullx and pully. Returns the , this will be negative if the fit failed to converge.
const Resolution & p_res() const
Return a constant reference to the momentum resolution.
Constrained_Top_Args(const Defaults &defs)
int k[5][pyjets_maxn]
Do a constrained kinematic fit of a event.
string const
Definition: compareJSON.py:14
Fourvec & p()
Return a reference to the four-momentum.
Fourvec_Constrainer_Args _fourvec_constrainer_args
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...
const Constrained_Top_Args _args
Define an interface for getting parameter settings.
Definition: Defaults.h:61
std::ostream & operator<<(std::ostream &s, const Constraint_Intermed &ci)
Output stream operator, print the content of this Constraint_Intermed to an output stream...
void mass_constraint(std::string s)
Specify a combination of objects that will be returned by the constrain() method as mass...
Represent a single object in a Fourvec_Event, this is just a dumb data container. Each object in a Fo...
Definition: Fourvec_Event.h:94
void add_constraint(std::string s)
Specify an additional constraint s for the problem. The format for s is described in the class descri...
Hold on to parameters for the Fourvec_Constrainer class.
Represent an event for kinematic fitting as a collection of four-momenta.