CMS 3D CMS Logo

Enumerations | Functions
secvtxMassTemplates.cc File Reference
#include "DataFormats/FWLite/interface/Handle.h"
#include "DataFormats/PatCandidates/interface/Jet.h"
#include "DataFormats/PatCandidates/interface/Muon.h"
#include "DataFormats/HepMCCandidate/interface/FlavorHistory.h"
#include "DataFormats/HepMCCandidate/interface/GenParticle.h"
#include "Math/GenVector/PxPyPzM4D.h"
#include "PhysicsTools/FWLite/interface/EventContainer.h"
#include "PhysicsTools/FWLite/interface/CommandLineParser.h"
#include "TROOT.h"
#include "TH2F.h"

Go to the source code of this file.

Enumerations

enum  { kNormalMode, kVqqMode, kLFMode, kWcMode }
 

Functions

void bookHistograms (fwlite::EventContainer &eventCont)
 
bool calcSampleName (fwlite::EventContainer &eventCont, string &sampleName)
 
int main (int argc, char *argv[])
 
void outputNameTagFunc (string &tag)
 

Enumeration Type Documentation

anonymous enum
Enumerator
kNormalMode 
kVqqMode 
kLFMode 
kWcMode 

Definition at line 21 of file secvtxMassTemplates.cc.

Function Documentation

void bookHistograms ( fwlite::EventContainer eventCont)

Definition at line 293 of file secvtxMassTemplates.cc.

References fwlite::EventContainer::add(), optutl::VariableMapCont::integerValue(), metsig::jet, kLFMode, kVqqMode, createfilelist::parser, fwlite::EventContainer::parser(), optutl::VariableMapCont::stringValue(), and GlobalPosition_Frontier_DevDB_cff::tag.

Referenced by main().

294 {
296  // First, come up with all possible base //
297  // names (E.g., Wbb, Wb2, etc.). //
299  CommandLineParser &parser = eventCont.parser();
300  CommandLineParser::SVec baseNameVec;
301  CommandLineParser::SVec beginningVec, endingVec;
302  switch ( parser.integerValue ("mode") )
303  {
304  case kVqqMode:
305  // We want Wbb, Wb2, .., Zbb, .. In this case, we completely
306  // ignore the sampleName that was passed in.
307  // Starts with
308  beginningVec.push_back ("X");
309  beginningVec.push_back ("W");
310  beginningVec.push_back ("Z");
311  // Ends with
312  endingVec.push_back( "bb" );
313  endingVec.push_back( "b2" );
314  endingVec.push_back( "cc" );
315  endingVec.push_back( "c2" );
316  for (CommandLineParser::SVecConstIter outerIter = beginningVec.begin();
317  beginningVec.end() != outerIter;
318  ++outerIter)
319  {
320  for (CommandLineParser::SVecConstIter innerIter = endingVec.begin();
321  endingVec.end() != innerIter;
322  ++innerIter)
323  {
324  baseNameVec.push_back( *outerIter + *innerIter);
325  } // for innerIter
326  } // for outerIter
327  break;
328  case kLFMode:
329  // just like the default case, except that we do have some
330  // heavy flavor pieces here, too.
331  baseNameVec.push_back(parser.stringValue ("sampleName") + "b3");
332  baseNameVec.push_back(parser.stringValue ("sampleName") + "c3");
333  // no break because to add just the name as well
334  default:
335  // We just want to use the sample name as it was given to us.
336  baseNameVec.push_back(parser.stringValue ("sampleName"));
337  } // for switch
338 
340  // Now the different tagging endings. //
342  CommandLineParser::SVec singleTagEndingVec, doubleTagEndingVec;
343  singleTagEndingVec.push_back ("_b");
344  singleTagEndingVec.push_back ("_c");
345  singleTagEndingVec.push_back ("_q");
346  doubleTagEndingVec.push_back ("_bb");
347  doubleTagEndingVec.push_back ("_cc");
348  doubleTagEndingVec.push_back ("_qq");
349  doubleTagEndingVec.push_back ("_bc");
350  doubleTagEndingVec.push_back ("_bq");
351  doubleTagEndingVec.push_back ("_cq");
352 
354  // Finally, let's put it all together. //
356  for (CommandLineParser::SVecConstIter baseIter = baseNameVec.begin();
357  baseNameVec.end() != baseIter;
358  ++baseIter)
359  {
361  // For each flavor, one jet/tag counting histogram. //
363  TString histName = *baseIter + "_jettag";
364  eventCont.add( new TH2F( histName, histName,
365  6, -0.5, 5.5,
366  3, -0.5, 2.5) );
367  for (int jet = 1; jet <= 5; ++jet)
368  {
369  for (int tag = 1; tag <= 2; ++tag)
370  {
372  // For each jet/tag, a single secvtx mass //
374  if (tag > jet) continue;
375  histName = *baseIter + Form ("_secvtxMass_%dj_%dt", jet, tag);
376  eventCont.add( new TH1F( histName, histName, 40, 0, 10) );
377  CommandLineParser::SVec *vecPtr = &singleTagEndingVec;
378  if (2 == tag)
379  {
380  vecPtr = &doubleTagEndingVec;
381  } // double tag
382  for (CommandLineParser::SVecConstIter tagIter = vecPtr->begin();
383  vecPtr->end() != tagIter;
384  ++tagIter)
385  {
387  // And different secvtx mass for each tag ending. //
389  histName = *baseIter + Form ("_secvtxMass_%dj_%dt", jet, tag)
390  + *tagIter;
391  eventCont.add( new TH1F( histName, histName, 40, 0, 10) );
392  } // for tagIter
393  } // for tag
394  } // for jet
395  } // for baseIter
396 }
std::string & stringValue(std::string key)
int & integerValue(std::string key)
void add(TH1 *histPtr, const std::string &directory="")
optutl::CommandLineParser & parser()
bool calcSampleName ( fwlite::EventContainer eventCont,
string &  sampleName 
)

Definition at line 399 of file secvtxMassTemplates.cc.

References funct::abs(), caloTruthProducer_cfi::genParticleCollection, fwlite::Handle< T >::getByLabel(), optutl::VariableMapCont::integerValue(), fwlite::Handle< T >::isValid(), kLFMode, kNormalMode, kWcMode, ALCARECOPromptCalibProdSiPixelAli0T_cff::mode, createfilelist::parser, fwlite::EventContainer::parser(), and optutl::VariableMapCont::stringValue().

Referenced by main().

400 {
401  // calculate sample name
402  CommandLineParser &parser = eventCont.parser();
403  sampleName = parser.stringValue ("sampleName");
404  int mode = parser.integerValue ("mode");
405 
407  // Normal Mode //
409  if (kNormalMode == mode)
410  {
411  // all we want is the sample name, so in this case we're done.
412  return true;
413  }
414  // Get the heavy flavor category
415  fwlite::Handle< unsigned int > heavyFlavorCategory;
416  heavyFlavorCategory.getByLabel (eventCont, "flavorHistoryFilter");
417  assert ( heavyFlavorCategory.isValid() );
418  int HFcat = (*heavyFlavorCategory);
419 
421  // Light Flavor Mode //
423  if (kLFMode == mode)
424  {
425  // Wqq
426  if (5 == HFcat)
427  {
428  sampleName += "b3";
429  } else if (6 == HFcat)
430  {
431  sampleName += "c3";
432  } else if (11 != HFcat)
433  {
434  // skip this event
435  return false;
436  } // else if ! 11
437  return true;
438  }
439 
441  // Wc Mode //
443  if (kWcMode == mode)
444  {
445  // Wc
446  if (4 != HFcat)
447  {
448  // skip this event
449  return false;
450  } // if not Wc
451  return true;
452  } // else if Wc
453 
455  // Vqq Mode //
457  // MadGraph (at least as CMS has implemented it) has this _lovely_
458  // feature that if the W or Z is far enough off-shell, it erases
459  // the W or Z from the event record. This means that in some
460  // number of cases, we won't be able to tell whether this is a W or
461  // Z event by looking for a W or Z in the GenParticle collection.
462  // (We'll eventually have to be more clever).
463  sampleName = "X";
465  genParticleCollection.getByLabel (eventCont, "genParticles");
466  assert ( genParticleCollection.isValid() );
467  // We don't know if it is a W, a Z, or neither
468  // Iterate over genParticles
469  const vector< reco::GenParticle>::const_iterator
470  kGenPartEnd = genParticleCollection->end();
471  for (vector< reco::GenParticle>::const_iterator gpIter =
472  genParticleCollection->begin();
473  gpIter != kGenPartEnd; ++gpIter )
474  {
475  if (gpIter->status() == 3 && std::abs(gpIter->pdgId()) == 23)
476  {
477  sampleName = "Z";
478  break;
479  }
480  if (gpIter->status() == 3 && std::abs(gpIter->pdgId()) == 24)
481  {
482  sampleName = "W";
483  break;
484  }
485  } // for gpIter
486  switch (HFcat)
487  {
488  // from:
489  // https://twiki.cern.ch/twiki/bin/viewauth/CMS/SWGuideFlavorHistory
490  // 1. W+bb with >= 2 jets from the ME (dr > 0.5)
491  // 2. W+b or W+bb with 1 jet from the ME
492  // 3. W+cc from the ME (dr > 0.5)
493  // 4. W+c or W+cc with 1 jet from the ME
494  // 5. W+bb with 1 jet from the parton shower (dr == 0.0)
495  // 6. W+cc with 1 jet from the parton shower (dr == 0.0)
496  // 7. W+bb with >= 2 partons but 1 jet from the ME (dr == 0.0)
497  // 8. W+cc with >= 2 partons but 1 jet from the ME (dr == 0.0)
498  // 9. W+bb with >= 2 partons but 2 jets from the PS (dr > 0.5)
499  // 10. W+cc with >= 2 partons but 2 jets from the PS (dr > 0.5)
500  // 11. Veto of all the previous (W+ light jets)
501  case 1:
502  sampleName += "bb";
503  break;
504  case 2:
505  // Sometimes this is referred to as 'b' (e.g., 'Wb'), but I
506  // am using the suffix '2' to keep this case clear for when
507  // we have charm (see below).
508  sampleName += "b2";
509  break;
510  case 3:
511  sampleName += "cc";
512  break;
513  case 4:
514  // We want to keep this case clear from real W + single charm
515  // produced (as opposed to two charm quarks produced and one
516  // goes down the beampipe), so we use 'c2' instead of 'c'.
517  sampleName += "c2";
518  break;
519  default:
520  // we don't want the rest of the cases. Return an empty
521  // string so we know.
522  return false;
523  } // switch HFcat
524  return true;
525 }
bool isValid() const
Definition: Handle.h:61
std::string & stringValue(std::string key)
void getByLabel(const P &iP, const char *iModuleLabel, const char *iProductInstanceLabel=0, const char *iProcessLabel=0)
Definition: Handle.h:91
int & integerValue(std::string key)
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
optutl::CommandLineParser & parser()
int main ( int  argc,
char *  argv[] 
)

Definition at line 53 of file secvtxMassTemplates.cc.

References funct::abs(), optutl::VariableMapCont::addOption(), fwlite::EventContainer::atEnd(), bookHistograms(), calcSampleName(), fwlite::Handle< T >::getByLabel(), fwlite::EventContainer::hist(), fwlite::Handle< T >::isValid(), min(), reco::TemplatedSecondaryVertexTagInfo< IPTI, VTX >::nVertices(), optutl::CommandLineParser::parseArguments(), createfilelist::parser, reco::TemplatedSecondaryVertexTagInfo< IPTI, VTX >::secondaryVertex(), optutl::VariableMapCont::stringValue(), caloDeepCSVTagInfos_cfi::svTagInfos, fwlite::EventContainer::toBegin(), reco::Vertex::tracks_begin(), and reco::Vertex::tracks_end().

54 {
56  // ////////////////////////// //
57  // // Command Line Options // //
58  // ////////////////////////// //
60 
61  // Tell people what this analysis code does and setup default options.
62  CommandLineParser parser ("Creates SecVtx Mass templates");
63 
65  // Add any command line options you would like here //
67  parser.addOption ("mode", CommandLineParser::kInteger,
68  "Normal(0), VQQ(1), LF(2), Wc(3)",
69  0);
70  parser.addOption ("sampleName", CommandLineParser::kString,
71  "Sample name (e.g., top, Wqq, etc.)");
72 
74  parser.stringValue ("outputFile") = "jetPt"; // .root added automatically
75 
76  // Parse the command line arguments
77  parser.parseArguments (argc, argv);
78 
80  // //////////////////////////// //
81  // // Create Event Container // //
82  // //////////////////////////// //
84 
85  // This object 'eventCont' is used both to get all information from the
86  // eventCont as well as to store histograms, etc.
87  // Parse the command line arguments
88  parser.parseArguments (argc, argv);
89 
91  // //////////////////////////// //
92  // // Create Event Container // //
93  // //////////////////////////// //
95 
96  fwlite::EventContainer eventCont (parser);
97 
99  // ////////////////////////////////// //
100  // // Begin Run // //
101  // // (e.g., book histograms, etc) // //
102  // ////////////////////////////////// //
104 
105  // Setup a style
106  gROOT->SetStyle ("Plain");
107 
108  // Book those histograms!
109  bookHistograms (eventCont);
110 
112  // //////////////// //
113  // // Event Loop // //
114  // //////////////// //
116 
117  for (eventCont.toBegin(); !eventCont.atEnd(); ++eventCont)
118  {
120  // Take What We Need From Event //
122  fwlite::Handle< vector< pat::Jet > > jetCollection;
123  jetCollection.getByLabel (eventCont, "selectedLayer1Jets");
124  assert ( jetCollection.isValid() );
125 
126  fwlite::Handle< vector< pat::Muon > > goodMuonCollection;
127  goodMuonCollection.getByLabel (eventCont, "goodLeptons");
128  assert ( goodMuonCollection.isValid() );
129 
130  // If we don't have any good leptons, don't bother
131  if (! goodMuonCollection->size() )
132  {
133  continue;
134  }
135 
136  // get the sample name for this event
137  string sampleName;
138  if ( ! calcSampleName (eventCont, sampleName) )
139  {
140  // We don't want this one.
141  continue;
142  }
143 
145  // Tagged Jets and Flavor Separator //
147  int numBottom = 0, numCharm = 0, numLight = 0;
148  int numTags = 0;
149  double sumVertexMass = 0.;
150  // Loop over the jets and find out which are tagged
151  const vector< pat::Jet >::const_iterator kJetEnd = jetCollection->end();
152  for (vector< pat::Jet >::const_iterator jetIter = jetCollection->begin();
153  kJetEnd != jetIter;
154  ++jetIter)
155  {
156  // Is this jet tagged and does it have a good secondary vertex
157  if( jetIter->bDiscriminator("simpleSecondaryVertexBJetTags") < 2.05 )
158  {
159  // This jet isn't tagged
160  continue;
161  }
163  = jetIter->tagInfoSecondaryVertex("secondaryVertex");
164  if ( svTagInfos->nVertices() <= 0 )
165  {
166  // Given that we are using simple secondary vertex
167  // tagging, I don't think this should ever happen.
168  // Maybe we should put a counter here just to check.
169  continue;
170  } // if we have no secondary verticies
171 
172  // count it
173  ++numTags;
174 
175  // What is the flavor of this jet
176  int jetFlavor = std::abs( jetIter->partonFlavour() );
177  if (5 == jetFlavor)
178  {
179  ++numBottom;
180  } // if bottom
181  else if (4 == jetFlavor)
182  {
183  ++numCharm;
184  } // if charm
185  else
186  {
187  ++numLight;
188  } // if light flavor
189 
191  // Calculate SecVtx Mass //
193  ROOT::Math::LorentzVector< ROOT::Math::PxPyPzM4D<double> > sumVec;
194  reco::CompositeCandidate vertexCand;
196  kEndTracks = svTagInfos->secondaryVertex(0).tracks_end();
197  for (reco::Vertex::trackRef_iterator trackIter =
198  svTagInfos->secondaryVertex(0).tracks_begin();
199  trackIter != kEndTracks;
200  ++trackIter )
201  {
202  const double kPionMass = 0.13957018;
203  ROOT::Math::LorentzVector< ROOT::Math::PxPyPzM4D<double> > vec;
204  vec.SetPx( (*trackIter)->px() );
205  vec.SetPy( (*trackIter)->py() );
206  vec.SetPz( (*trackIter)->pz() );
207  vec.SetM (kPionMass);
208  sumVec += vec;
209  } // for trackIter
210  sumVertexMass += sumVec.M();
211  if (2 == numTags)
212  {
213  // We've got enough. Stop.
214  break;
215  } // if we have enough tags
216  } // for jet
217 
219  // General Accounting //
221  int numJets = std::min( (int) jetCollection->size(), 5 );
222  eventCont.hist( sampleName + "_jettag")->Fill (numJets, numTags);
223 
224  // If we don't have any tags, don't bother going on
225  if ( ! numTags)
226  {
227  continue;
228  }
229 
231  // Calculate average SecVtx mass and //
232  // fill appropriate histograms. //
234  sumVertexMass /= numTags;
235  string whichtag = "";
236  if (1 == numTags)
237  {
238  // single tag
239  if (numBottom) whichtag = "_b";
240  else if (numCharm) whichtag = "_c";
241  else if (numLight) whichtag = "_q";
242  else whichtag = "_X";
243  } else {
244  // double tags
245  if (2 == numBottom) whichtag = "_bb";
246  else if (2 == numCharm) whichtag = "_cc";
247  else if (2 == numLight) whichtag = "_qq";
248  else if (numBottom && numCharm) whichtag = "_bc";
249  else if (numBottom && numLight) whichtag = "_bq";
250  else if (numCharm && numLight) whichtag = "_bq";
251  else whichtag = "_XX";
252  } // if two tags
253  string massName = sampleName
254  + Form("_secvtxMass_%dj_%dt", numJets, numTags);
255  eventCont.hist(massName)->Fill (sumVertexMass);
256  eventCont.hist(massName + whichtag)->Fill (sumVertexMass);
257  } // for eventCont
258 
260  // ////////////////// //
261  // // Clean Up Job // //
262  // ////////////////// //
264 
265  // Histograms will be automatically written to the root file
266  // specificed by command line options.
267 
268  // All done! Bye bye.
269  return 0;
270 }
trackRef_iterator tracks_end() const
last iterator over tracks
Definition: Vertex.cc:81
bool isValid() const
Definition: Handle.h:61
const VTX & secondaryVertex(unsigned int index) const
bool calcSampleName(fwlite::EventContainer &eventCont, string &sampleName)
void bookHistograms(fwlite::EventContainer &eventCont)
void getByLabel(const P &iP, const char *iModuleLabel, const char *iProductInstanceLabel=0, const char *iProcessLabel=0)
Definition: Handle.h:91
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
T min(T a, T b)
Definition: MathUtil.h:58
std::vector< TrackBaseRef >::const_iterator trackRef_iterator
The iteratator for the vector<TrackRef>
Definition: Vertex.h:37
trackRef_iterator tracks_begin() const
first iterator over tracks
Definition: Vertex.cc:76
void outputNameTagFunc ( string &  tag)

Definition at line 280 of file secvtxMassTemplates.cc.

281 {
282  // If you do not want to give you output filename any "tag" based
283  // on the command line options, simply do nothing here. This
284  // function is designed to be called by fwlite::EventContainer constructor.
285 
286  // if ( boolValue ("someCondition") )
287  // {
288  // tag += "_someCond";
289  // }
290 }