CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
Constrained_Z.cc
Go to the documentation of this file.
1 //
2 //
3 // File: Constrained_Z.cc
4 // Purpose: Do kinematic fitting for a (Z->ll)+jets event.
5 // Created: Apr, 2004, sss
6 //
7 // CMSSW File : src/Constrained_Z.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 <stdio.h>
44 
45 
46 namespace hitfit {
47 
48 
49 //*************************************************************************
50 // Argument handling.
51 //
52 
53 
55 //
56 // Purpose: Constructor.
57 //
58 // Inputs:
59 // defs - The Defaults instance from which to initialize.
60 //
61  : _zmass (defs.get_float ("zmass")),
62  _fourvec_constrainer_args (defs)
63 {
64 }
65 
66 
68 //
69 // Purpose: Return the zmass parameter.
70 // See the header for documentation.
71 //
72 {
73  return _zmass;
74 }
75 
76 
79 //
80 // Purpose: Return the contained subobject parameters.
81 //
82 {
84 }
85 
86 
87 //*************************************************************************
88 
89 
91 //
92 // Purpose: Constructor.
93 //
94 // Inputs:
95 // args - The parameter settings for this instance.
96 //
97  : _args (args),
98  _constrainer (args.fourvec_constrainer_args())
99 {
100  char buf[256];
101  sprintf (buf, "(%d) = %f", lepton_label, _args.zmass());
103 }
104 
105 
106 namespace {
107 
108 
123 FE_Obj make_fe_obj (const Lepjets_Event_Lep& obj, double mass, int type)
124 //
125 // Purpose: Helper to create an object to put into the Fourvec_Event.
126 //
127 // Inputs:
128 // obj - The input object.
129 // mass - The mass to which it should be constrained.
130 // type - The type to assign it.
131 //
132 // Returns:
133 // The constructed FE_Obj.
134 //
135 {
136  return FE_Obj (obj.p(), mass, type,
137  obj.p_sigma(), obj.eta_sigma(), obj.phi_sigma(),
138  obj.res().p_res().inverse());
139 }
140 
141 
157 void do_import (const Lepjets_Event& ev, Fourvec_Event& fe)
158 //
159 // Purpose: Convert from a Lepjets_Event to a Fourvec_Event.
160 //
161 // Inputs:
162 // ev - The input event.
163 //
164 // Outputs:
165 // fe - The initialized Fourvec_Event.
166 //
167 {
168  assert (ev.nleps() >= 2);
169  fe.add (make_fe_obj (ev.lep(0), 0, lepton_label));
170  fe.add (make_fe_obj (ev.lep(1), 0, lepton_label));
171 
172  for (std::vector<Lepjets_Event_Jet>::size_type j=0; j < ev.njets(); j++)
173  fe.add (make_fe_obj (ev.jet(j), 0, isr_label));
174 
175  Fourvec kt = ev.kt ();
176  fe.set_kt_error (ev.kt_res().sigma (kt.x()),
177  ev.kt_res().sigma (kt.y()),
178  0);
179  fe.set_x_p (ev.met());
180 }
181 
182 
199 void do_export (const Fourvec_Event& fe, Lepjets_Event& ev)
200 //
201 // Purpose: Convert from a Fourvec_Event to a Lepjets_Event.
202 //
203 // Inputs:
204 // fe - The input event.
205 // ev - The original Lepjets_Event.
206 //
207 // Outputs:
208 // ev - The updated Lepjets_Event.
209 //
210 {
211  ev.lep(0).p() = fe.obj(0).p;
212  ev.lep(1).p() = fe.obj(1).p;
213 
214  for (std::vector<Lepjets_Event_Jet>::size_type j=0, k=1; j < ev.njets(); j++)
215  ev.jet(j).p() = fe.obj(k++).p;
216 
217  Fourvec nu;
218  ev.met() = nu;
219 }
220 
221 
222 } // unnamed namespace
223 
224 
226 //
227 // Purpose: Do a constrained fit for EV.
228 // Returns the pull quantities in PULL.
229 // Returns the chisq; this will be < 0 if the fit failed
230 // to converge.
231 //
232 // Inputs:
233 // ev - The event we're fitting.
234 //
235 // Outputs:
236 // ev - The fitted event.
237 // pull - Pull quantities for well-measured variables.
238 //
239 // Returns:
240 // The fit chisq, or < 0 if the fit didn't converge.
241 //
242 {
243  Fourvec_Event fe;
244  do_import (ev, fe);
245  Column_Vector pully;
246  double m, sigm;
247  double chisq = _constrainer.constrain (fe, m, sigm, pull, pully);
248  do_export (fe, ev);
249 
250  return chisq;
251 }
252 
253 
261 std::ostream& operator<< (std::ostream& s, const Constrained_Z& cz)
262 //
263 // Purpose: Print the object to S.
264 //
265 // Inputs:
266 // s - The stream to which to write.
267 // cz - The object to write.
268 //
269 // Returns:
270 // The stream S.
271 //
272 {
273  return s << cz._constrainer;
274 }
275 
276 
277 } // namespace hitfit
type
Definition: HCALResponse.h:21
const Constrained_Z_Args & _args
Define an abstract interface for getting parameter settings.
Represent a lepton in an instance of Lepjets_Event class. This class hold the following information: ...
Represent an event for kinematic fitting as a collection of four-momenta. Each object is represented ...
bool inverse() const
Return the setting of the inverse flag.
Definition: Resolution.cc:171
CLHEP::HepVector Column_Vector
Definition: matutil.h:66
uint16_t size_type
double eta_sigma() const
Return the uncertainty in pseudorapidity .
Do a constrained kinematic fitting for a event.
Do a constrained kinematic fit of a event.
const Vector_Resolution & res() const
Return a constant reference to the resolution.
double phi_sigma() const
Return the uncertainty in azimuthal angle .
Fourvec_Constrainer _constrainer
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, Column_Vector &pull)
Do a constrained fit of event. Returns the pull quantities in pull. 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.
int k[5][pyjets_maxn]
Fourvec_Constrainer_Args _fourvec_constrainer_args
Definition: Constrained_Z.h:98
Constrained_Z(const Constrained_Z_Args &args)
Constructor, create an instance of the Constrained_Z object from the argument object.
Hold on to parameters for the Constrained_Z class.
Definition: Constrained_Z.h:55
string const
Definition: compareJSON.py:14
Fourvec & p()
Return a reference to the four-momentum.
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 Fourvec_Constrainer_Args & fourvec_constrainer_args() const
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...
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
Constrained_Z_Args(const Defaults &defs)
Constructor, initialize from a Defaults object.
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.