CMS 3D CMS Logo

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 
41 #include <ostream>
42 #include <cassert>
43 #include <cstdio>
44 
45 using std::ostream;
46 
47 namespace hitfit {
48 
49  //*************************************************************************
50  // Argument handling.
51  //
52 
54  //
55  // Purpose: Constructor.
56  //
57  // Inputs:
58  // defs - The Defaults instance from which to initialize.
59  //
60  : _bmass(defs.get_float("bmass")),
61  _fourvec_constrainer_args(defs),
62  _equal_side(defs.get_bool("equal_side"))
63 
64  {}
65 
67  //
68  // Purpose: Return the bmass parameter.
69  // See the header for documentation.
70  //
71  {
72  return _bmass;
73  }
74 
76  //
77  // Purpose: Return the contained subobject parameters.
78  //
79  {
81  }
82 
84  //
85  // Purpose: Return the equal_side parameter
86  // See the header for documentation.
87  //
88  {
89  return _equal_side;
90  }
91 
92  //*************************************************************************
93 
94  Constrained_Top::Constrained_Top(const Constrained_Top_Args& args, double lepw_mass, double hadw_mass, double top_mass)
95  //
96  // Purpose: Constructor.
97  //
98  // Inputs:
99  // args - The parameter settings for this instance.
100  // lepw_mass - The mass to which the leptonic W should be constrained,
101  // or 0 to skip this constraint.
102  // hadw_mass - The mass to which the hadronic W should be constrained,
103  // or 0 to skip this constraint.
104  // top_mass - The mass to which the top quarks should be constrained,
105  // or 0 to skip this constraint.
106  //
107  : _args(args), _constrainer(args.fourvec_constrainer_args()) {
108  char buf[256];
109  if (lepw_mass > 0) {
110  sprintf(buf, "(%d %d) = %f", nu_label, lepton_label, lepw_mass);
112  }
113 
114  if (hadw_mass > 0) {
115  sprintf(buf, "(%d %d) = %f", hadw1_label, hadw2_label, hadw_mass);
117  }
118 
119  if (args.equal_side()) {
120  sprintf(buf, "(%d %d %d) = (%d %d %d)", nu_label, lepton_label, lepb_label, hadw1_label, hadw2_label, hadb_label);
122  }
123 
124  if (top_mass > 0) {
125  sprintf(buf, "(%d %d %d) = %f", hadw1_label, hadw2_label, hadb_label, top_mass);
127  } else {
128  sprintf(buf, "(%d %d %d) = 0", hadw1_label, hadw2_label, hadb_label);
130  }
131  }
132 
133  namespace {
134 
149  FE_Obj make_fe_obj(const Lepjets_Event_Lep& obj, double mass, int type)
150  //
151  // Purpose: Helper to create an object to put into the Fourvec_Event.
152  //
153  // Inputs:
154  // obj - The input object.
155  // mass - The mass to which it should be constrained.
156  // type - The type to assign it.
157  //
158  // Returns:
159  // The constructed FE_Obj.
160  //
161  {
162  return FE_Obj(obj.p(), mass, type, obj.p_sigma(), obj.eta_sigma(), obj.phi_sigma(), obj.res().p_res().inverse());
163  }
164 
183  void do_import(const Lepjets_Event& ev, double bmass, Fourvec_Event& fe)
184  //
185  // Purpose: Convert from a Lepjets_Event to a Fourvec_Event.
186  //
187  // Inputs:
188  // ev - The input event.
189  // bmass - The mass to which b-jets should be fixed.
190  //
191  // Outputs:
192  // fe - The initialized Fourvec_Event.
193  //
194  {
195  assert(ev.nleps() == 1);
196  fe.add(make_fe_obj(ev.lep(0), 0, lepton_label));
197 
198  bool saw_lepb = false;
199  bool saw_hadb = false;
200  for (std::vector<Lepjets_Event_Jet>::size_type j = 0; j < ev.njets(); j++) {
201  if (ev.jet(j).type() == isr_label || ev.jet(j).type() == higgs_label)
202  continue;
203  double mass = 0;
204  if (ev.jet(j).type() == lepb_label && !saw_lepb) {
205  mass = bmass;
206  saw_lepb = true;
207  } else if (ev.jet(j).type() == hadb_label && !saw_hadb) {
208  mass = bmass;
209  saw_hadb = true;
210  }
211  fe.add(make_fe_obj(ev.jet(j), mass, ev.jet(j).type()));
212  }
213 
214  fe.set_nu_p(ev.met());
215  Fourvec kt = ev.kt();
216  fe.set_kt_error(ev.kt_res().sigma(kt.x()), ev.kt_res().sigma(kt.y()), 0);
217  }
218 
235  void do_export(const Fourvec_Event& fe, Lepjets_Event& ev)
236  //
237  // Purpose: Convert from a Fourvec_Event to a Lepjets_Event.
238  //
239  // Inputs:
240  // fe - The input event.
241  // ev - The original Lepjets_Event.
242  //
243  // Outputs:
244  // ev - The updated Lepjets_Event.
245  //
246  {
247  ev.lep(0).p() = fe.obj(0).p;
248  for (std::vector<Lepjets_Event_Jet>::size_type j = 0, k = 1; j < ev.njets(); j++) {
249  if (ev.jet(j).type() == isr_label || ev.jet(j).type() == higgs_label)
250  continue;
251  ev.jet(j).p() = fe.obj(k++).p;
252  }
253 
254  ev.met() = fe.nu();
255  }
256 
257  } // unnamed namespace
258 
260  Lepjets_Event& ev, double& mt, double& sigmt, Column_Vector& pullx, Column_Vector& pully)
261  //
262  // Purpose: Do a constrained fit for EV. Returns the top mass and
263  // its error in MT and SIGMT, and the pull quantities in PULLX and
264  // PULLY. Returns the chisq; this will be < 0 if the fit failed
265  // to converge.
266  //
267  // Inputs:
268  // ev - The event we're fitting.
269  //
270  // Outputs:
271  // ev - The fitted event.
272  // mt - Requested invariant mass.
273  // sigmt - Uncertainty on mt.
274  // pullx - Pull quantities for well-measured variables.
275  // pully - Pull quantities for poorly-measured variables.
276  //
277  // Returns:
278  // The fit chisq, or < 0 if the fit didn't converge.
279  //
280  {
281  Fourvec_Event fe;
282  do_import(ev, _args.bmass(), fe);
283  double chisq = _constrainer.constrain(fe, mt, sigmt, pullx, pully);
284  do_export(fe, ev);
285 
286  return chisq;
287  }
288 
298  std::ostream& operator<<(std::ostream& s, const Constrained_Top& ct)
299  //
300  // Purpose: Print the object to S.
301  //
302  // Inputs:
303  // s - The stream to which to write.
304  // ct - The object to write.
305  //
306  // Returns:
307  // The stream S.
308  //
309  {
310  return s << ct._constrainer;
311  }
312 
313 } // namespace hitfit
Constrained_Top.h
Do a constrained kinematic fit of a event.
writedatasetfile.args
args
Definition: writedatasetfile.py:18
hitfit::Fourvec_Constrainer::mass_constraint
void mass_constraint(std::string s)
Specify a combination of objects that will be returned by the constrain() method as mass....
Definition: Fourvec_Constrainer.cc:195
hitfit::Column_Vector
CLHEP::HepVector Column_Vector
Definition: matutil.h:63
hitfit
Definition: Base_Constrainer.h:43
hitfit::Fourvec
CLHEP::HepLorentzVector Fourvec
Typedef for a HepLorentzVector.
Definition: fourvec.h:55
hitfit::hadb_label
Definition: Lepjets_Event_Jet.h:57
Fourvec_Event.h
Represent an event for kinematic fitting as a collection of four-momenta.
hitfit::lepb_label
Definition: Lepjets_Event_Jet.h:56
hitfit::Constrained_Top::Constrained_Top
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 ...
Definition: Constrained_Top.cc:94
hitfit::Constrained_Top_Args::equal_side
bool equal_side() const
Definition: Constrained_Top.cc:83
cms::cuda::assert
assert(be >=bs)
Defaults.h
Define an abstract interface for getting parameter settings.
watchdog.const
const
Definition: watchdog.py:83
hitfit::Constrained_Top_Args::Constrained_Top_Args
Constrained_Top_Args(const Defaults &defs)
Definition: Constrained_Top.cc:53
hitfit::Fourvec_Event
Represent an event for kinematic fitting as a collection of four-momenta. Each object is represented ...
Definition: Fourvec_Event.h:193
hitfit::hadw1_label
Definition: Lepjets_Event_Jet.h:58
hitfit::FE_Obj
Represent a single object in a Fourvec_Event, this is just a dumb data container. Each object in a Fo...
Definition: Fourvec_Event.h:90
Lepjets_Event.h
Represent a simple event consisting of lepton(s) and jet(s).
alignCSCRings.s
s
Definition: alignCSCRings.py:92
trigger::size_type
uint16_t size_type
Definition: TriggerTypeDefs.h:18
hitfit::Lepjets_Event
Represent a simple event consisting of lepton(s) and jet(s). An instance of this class holds a list o...
Definition: Lepjets_Event.h:62
hitfit::hadw2_label
Definition: Lepjets_Event_Jet.h:59
hitfit::Constrained_Top_Args::_fourvec_constrainer_args
Fourvec_Constrainer_Args _fourvec_constrainer_args
Definition: Constrained_Top.h:104
dqmdumpme.k
k
Definition: dqmdumpme.py:60
hitfit::operator<<
std::ostream & operator<<(std::ostream &s, const Constraint_Intermed &ci)
Output stream operator, print the content of this Constraint_Intermed to an output stream.
Definition: Constraint_Intermed.cc:268
hitfit::Constrained_Top_Args::fourvec_constrainer_args
const Fourvec_Constrainer_Args & fourvec_constrainer_args() const
Definition: Constrained_Top.cc:75
getGTfromDQMFile.obj
obj
Definition: getGTfromDQMFile.py:32
hitfit::Constrained_Top_Args
Hold on to parameters for the Constrained_Top class.
Definition: Constrained_Top.h:52
type
type
Definition: SiPixelVCal_PayloadInspector.cc:37
TtSemiLepEvtBuilder_cfi.sigmt
sigmt
Definition: TtSemiLepEvtBuilder_cfi.py:48
gainCalibHelper::gainCalibPI::type
type
Definition: SiPixelGainCalibHelper.h:39
hitfit::Fourvec_Constrainer::constrain
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,...
Definition: Fourvec_Constrainer.cc:1612
hitfit::Constrained_Top::constrain
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,...
Definition: Constrained_Top.cc:259
visDQMUpload.buf
buf
Definition: visDQMUpload.py:154
hitfit::higgs_label
Definition: Lepjets_Event_Jet.h:60
hitfit::Constrained_Top_Args::bmass
double bmass() const
Definition: Constrained_Top.cc:66
hitfit::nu_label
const int nu_label
Definition: Fourvec_Event.h:178
hitfit::Fourvec_Constrainer::add_constraint
void add_constraint(std::string s)
Specify an additional constraint s for the problem. The format for s is described in the class descri...
Definition: Fourvec_Constrainer.cc:183
hitfit::Constrained_Top::_constrainer
Fourvec_Constrainer _constrainer
Definition: Constrained_Top.h:198
hitfit::Constrained_Top_Args::_bmass
double _bmass
Definition: Constrained_Top.h:98
ev
bool ev
Definition: Hydjet2Hadronizer.cc:95
EgHLTOffHistBins_cfi.mass
mass
Definition: EgHLTOffHistBins_cfi.py:34
hitfit::Constrained_Top::_args
const Constrained_Top_Args _args
Definition: Constrained_Top.h:192
hitfit::Defaults
Define an interface for getting parameter settings.
Definition: Defaults.h:57
TtSemiLepEvtBuilder_cfi.mt
mt
Definition: TtSemiLepEvtBuilder_cfi.py:47
hitfit::Constrained_Top
Do a constrained kinematic fitting for a event.
Definition: Constrained_Top.h:120
hitfit::Fourvec_Constrainer_Args
Hold on to parameters for the Fourvec_Constrainer class.
Definition: Fourvec_Constrainer.h:93
dqmiolumiharvest.j
j
Definition: dqmiolumiharvest.py:66
hitfit::isr_label
Definition: Lepjets_Event_Jet.h:55
hitfit::lepton_label
Definition: Lepjets_Event_Lep.h:52
hitfit::Lepjets_Event_Lep
Represent a lepton in an instance of Lepjets_Event class. This class hold the following information:
Definition: Lepjets_Event_Lep.h:64
hitfit::Constrained_Top_Args::_equal_side
bool _equal_side
Definition: Constrained_Top.h:110