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 // $Id: Constrained_Top.cc,v 1.1 2011/05/26 09:46:59 mseidel Exp $
3 //
4 // File: src/Constrained_Top.cc
5 // Purpose: Do kinematic fitting for a ttbar -> ljets event.
6 // Created: Jul, 2000, sss, based on run 1 mass analysis code.
7 //
8 // CMSSW File : src/Constrained_Top.cc
9 // Original Author : Scott Stuart Snyder <snyder@bnl.gov> for D0
10 // Imported to CMSSW by Haryo Sumowidagdo <Suharyo.Sumowidagdo@cern.ch>
11 //
12 
13 
43 #include <ostream>
44 #include <cassert>
45 #include <stdio.h>
46 
47 
48 using std::ostream;
49 
50 
51 namespace hitfit {
52 
53 
54 //*************************************************************************
55 // Argument handling.
56 //
57 
58 
60 //
61 // Purpose: Constructor.
62 //
63 // Inputs:
64 // defs - The Defaults instance from which to initialize.
65 //
66  : _bmass (defs.get_float ("bmass")),
67  _fourvec_constrainer_args (defs),
68  _equal_side(defs.get_bool("equal_side"))
69 
70 {
71 }
72 
73 
75 //
76 // Purpose: Return the bmass parameter.
77 // See the header for documentation.
78 //
79 {
80  return _bmass;
81 }
82 
83 
86 //
87 // Purpose: Return the contained subobject parameters.
88 //
89 {
91 }
92 
93 
95 //
96 // Purpose: Return the equal_side parameter
97 // See the header for documentation.
98 //
99 {
100  return _equal_side;
101 }
102 
103 
104 //*************************************************************************
105 
106 
108  double lepw_mass,
109  double hadw_mass,
110  double top_mass)
111 //
112 // Purpose: Constructor.
113 //
114 // Inputs:
115 // args - The parameter settings for this instance.
116 // lepw_mass - The mass to which the leptonic W should be constrained,
117 // or 0 to skip this constraint.
118 // hadw_mass - The mass to which the hadronic W should be constrained,
119 // or 0 to skip this constraint.
120 // top_mass - The mass to which the top quarks should be constrained,
121 // or 0 to skip this constraint.
122 //
123  : _args (args),
124  _constrainer (args.fourvec_constrainer_args())
125 {
126  char buf[256];
127  if (lepw_mass > 0) {
128  sprintf (buf, "(%d %d) = %f", nu_label, lepton_label, lepw_mass);
130  }
131 
132  if (hadw_mass > 0) {
133  sprintf (buf, "(%d %d) = %f", hadw1_label, hadw2_label, hadw_mass);
135  }
136 
137  if (args.equal_side()) {
138  sprintf (buf, "(%d %d %d) = (%d %d %d)",
142  }
143 
144  if (top_mass > 0) {
145  sprintf (buf, "(%d %d %d) = %f",
146  hadw1_label, hadw2_label, hadb_label, top_mass);
148  }
149  else {
150  sprintf (buf, "(%d %d %d) = 0", hadw1_label, hadw2_label, hadb_label);
152  }
153 }
154 
155 namespace {
156 
157 
172 FE_Obj make_fe_obj (const Lepjets_Event_Lep& obj, double mass, int type)
173 //
174 // Purpose: Helper to create an object to put into the Fourvec_Event.
175 //
176 // Inputs:
177 // obj - The input object.
178 // mass - The mass to which it should be constrained.
179 // type - The type to assign it.
180 //
181 // Returns:
182 // The constructed FE_Obj.
183 //
184 {
185  return FE_Obj (obj.p(), mass, type,
186  obj.p_sigma(), obj.eta_sigma(), obj.phi_sigma(),
187  obj.res().p_res().inverse());
188 }
189 
190 
209 void do_import (const Lepjets_Event& ev, double bmass, Fourvec_Event& fe)
210 //
211 // Purpose: Convert from a Lepjets_Event to a Fourvec_Event.
212 //
213 // Inputs:
214 // ev - The input event.
215 // bmass - The mass to which b-jets should be fixed.
216 //
217 // Outputs:
218 // fe - The initialized Fourvec_Event.
219 //
220 {
221  assert (ev.nleps() == 1);
222  fe.add (make_fe_obj (ev.lep(0), 0, lepton_label));
223 
224  bool saw_lepb = false;
225  bool saw_hadb = false;
226  for (std::vector<Lepjets_Event_Jet>::size_type j=0; j < ev.njets(); j++) {
227  if (ev.jet(j).type() == isr_label || ev.jet(j).type() == higgs_label)
228  continue;
229  double mass = 0;
230  if (ev.jet(j).type() == lepb_label && !saw_lepb) {
231  mass = bmass;
232  saw_lepb = true;
233  }
234  else if (ev.jet(j).type() == hadb_label && !saw_hadb) {
235  mass = bmass;
236  saw_hadb = true;
237  }
238  fe.add (make_fe_obj (ev.jet(j), mass, ev.jet(j).type()));
239  }
240 
241  fe.set_nu_p (ev.met());
242  Fourvec kt = ev.kt ();
243  fe.set_kt_error (ev.kt_res().sigma (kt.x()),
244  ev.kt_res().sigma (kt.y()),
245  0);
246 }
247 
248 
265 void do_export (const Fourvec_Event& fe, Lepjets_Event& ev)
266 //
267 // Purpose: Convert from a Fourvec_Event to a Lepjets_Event.
268 //
269 // Inputs:
270 // fe - The input event.
271 // ev - The original Lepjets_Event.
272 //
273 // Outputs:
274 // ev - The updated Lepjets_Event.
275 //
276 {
277  ev.lep(0).p() = fe.obj(0).p;
278  for (std::vector<Lepjets_Event_Jet>::size_type j=0, k=1; j < ev.njets(); j++) {
279  if (ev.jet(j).type() == isr_label || ev.jet(j).type() == higgs_label)
280  continue;
281  ev.jet(j).p() = fe.obj(k++).p;
282  }
283 
284  ev.met() = fe.nu();
285 }
286 
287 
288 } // unnamed namespace
289 
290 
292  double& mt,
293  double& sigmt,
294  Column_Vector& pullx,
295  Column_Vector& pully)
296 //
297 // Purpose: Do a constrained fit for EV. Returns the top mass and
298 // its error in MT and SIGMT, and the pull quantities in PULLX and
299 // PULLY. Returns the chisq; this will be < 0 if the fit failed
300 // to converge.
301 //
302 // Inputs:
303 // ev - The event we're fitting.
304 //
305 // Outputs:
306 // ev - The fitted event.
307 // mt - Requested invariant mass.
308 // sigmt - Uncertainty on mt.
309 // pullx - Pull quantities for well-measured variables.
310 // pully - Pull quantities for poorly-measured variables.
311 //
312 // Returns:
313 // The fit chisq, or < 0 if the fit didn't converge.
314 //
315 {
317  do_import (ev, _args.bmass (), fe);
318  double chisq = _constrainer.constrain (fe, mt, sigmt, pullx, pully);
319  do_export (fe, ev);
320 
321  return chisq;
322 }
323 
324 
334 std::ostream& operator<< (std::ostream& s, const Constrained_Top& ct)
335 //
336 // Purpose: Print the object to S.
337 //
338 // Inputs:
339 // s - The stream to which to write.
340 // ct - The object to write.
341 //
342 // Returns:
343 // The stream S.
344 //
345 {
346  return s << ct._constrainer;
347 }
348 
349 
350 } // namespace hitfit
type
Definition: HCALResponse.h:22
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:172
CLHEP::HepVector Column_Vector
Definition: matutil.h:67
const int nu_label
Fourvec_Constrainer _constrainer
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:67
CLHEP::HepLorentzVector Fourvec
Typedef for a HepLorentzVector.
Definition: fourvec.h:58
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]
double & fe
Definition: VDTMath.h:196
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
dictionary 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:62
tuple mass
Definition: scaleCards.py:27
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:95
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.