CMS 3D CMS Logo

CalcTopMassExample.cc
Go to the documentation of this file.
1 //
2 // Example to calculated top mass from genjets corrected on the fly with parton calibration
3 // Top mass is calculated from:
4 // 1) Uncorrected jets
5 // 2) Corrected jets with flavour ttbar correction assuming jet mass = quark mass
6 // 3) The same as 3 assuming massless jets
7 // 4) Corrected jets with mixed ttbar correction assuming massless jets
8 //
9 // Author: Attilio Santocchia
10 // Date: 15.06.2009
11 // Tested on CMSSW_3_1_0
12 //
13 // user include files
22 
28 
30 
34 
39 
41 
42 #include "TFile.h"
43 #include "TH1.h"
44 #include "TF1.h"
45 #include "TLorentzVector.h"
46 
47 // system include files
48 #include <memory>
49 #include <string>
50 #include <iostream>
51 #include <vector>
52 #include <Math/VectorUtil.h>
53 #include <TMath.h>
54 
55 class calcTopMass : public edm::one::EDAnalyzer<edm::one::SharedResources> {
56 public:
57  explicit calcTopMass(const edm::ParameterSet &);
58  ~calcTopMass() override{};
59  void analyze(const edm::Event &iEvent, const edm::EventSetup &iSetup) override;
60 
61 private:
65 
70 
71  float bMass;
72  float cMass;
73  float qMass;
74 
75  TH1F *hMassNoCorr;
76  TH1F *hMassCorFl0;
77  TH1F *hMassCorFlM;
78  TH1F *hMassCorMix;
79 };
80 
81 using namespace std;
82 using namespace reco;
83 using namespace edm;
84 using namespace ROOT::Math::VectorUtil;
85 
87  sourceByRefer_ = iConfig.getParameter<InputTag>("srcByReference");
88  m_qT_CorrectorName = iConfig.getParameter<std::string>("qTopCorrector");
89  m_cT_CorrectorName = iConfig.getParameter<std::string>("cTopCorrector");
90  m_bT_CorrectorName = iConfig.getParameter<std::string>("bTopCorrector");
91  m_tT_CorrectorName = iConfig.getParameter<std::string>("tTopCorrector");
92 
93  bMass = 4.5;
94  cMass = 1.5;
95  qMass = 0.3;
96 
97  usesResource(TFileService::kSharedResource);
99  hMassNoCorr = fs->make<TH1F>("hMassNoCorr", "", 100, 100, 300);
100  hMassCorFl0 = fs->make<TH1F>("hMassCorFl0", "", 100, 100, 300);
101  hMassCorFlM = fs->make<TH1F>("hMassCorFlM", "", 100, 100, 300);
102  hMassCorMix = fs->make<TH1F>("hMassCorMix", "", 100, 100, 300);
103 }
104 
106  cout << "[calcTopMass] analysing event " << iEvent.id() << endl;
107  try {
108  iEvent.getByLabel(sourceByRefer_, theJetPartonMatch);
109  } catch (std::exception &ce) {
110  cerr << "[calcTopMass] caught std::exception " << ce.what() << endl;
111  return;
112  }
113 
114  // get all correctors from top events
115  const JetCorrector *qTopCorrector = JetCorrector::getJetCorrector(m_qT_CorrectorName, iSetup);
116  const JetCorrector *cTopCorrector = JetCorrector::getJetCorrector(m_cT_CorrectorName, iSetup);
117  const JetCorrector *bTopCorrector = JetCorrector::getJetCorrector(m_bT_CorrectorName, iSetup);
118  const JetCorrector *tTopCorrector = JetCorrector::getJetCorrector(m_tT_CorrectorName, iSetup);
119 
120  TLorentzVector jetsNoCorr[6];
121  TLorentzVector jetsCorFl0[6];
122  TLorentzVector jetsCorFlM[6];
123  TLorentzVector jetsCorMix[6];
124 
125  bool isQuarkFound[6] = {false};
126 
127  for (JetMatchedPartonsCollection::const_iterator j = theJetPartonMatch->begin(); j != theJetPartonMatch->end(); j++) {
128  const math::XYZTLorentzVector theJet = (*j).first.get()->p4();
129  const MatchedPartons aMatch = (*j).second;
130  const GenParticleRef &thePhyDef = aMatch.physicsDefinitionParton();
131 
132  if (thePhyDef.isNonnull()) {
133  int particPDG = thePhyDef.get()->pdgId();
134 
135  if (particPDG == 5) { //b from t
136  double bTcorr = bTopCorrector->correction(theJet);
137  double tTcorr = tTopCorrector->correction(theJet);
138  jetsNoCorr[0].SetPtEtaPhiM(theJet.pt(), theJet.eta(), theJet.phi(), theJet.mass());
139  jetsCorFl0[0].SetPtEtaPhiM(theJet.pt() * bTcorr, theJet.eta(), theJet.phi(), 0);
140  jetsCorFlM[0].SetPtEtaPhiM(theJet.pt() * bTcorr, theJet.eta(), theJet.phi(), bMass);
141  jetsCorMix[0].SetPtEtaPhiM(theJet.pt() * tTcorr, theJet.eta(), theJet.phi(), 0);
142  isQuarkFound[0] = true;
143  } else if (particPDG == -5) { //bbar from tbar
144  double bTcorr = bTopCorrector->correction(theJet);
145  double tTcorr = tTopCorrector->correction(theJet);
146  jetsNoCorr[3].SetPtEtaPhiM(theJet.pt(), theJet.eta(), theJet.phi(), theJet.mass());
147  jetsCorFl0[3].SetPtEtaPhiM(theJet.pt() * bTcorr, theJet.eta(), theJet.phi(), 0);
148  jetsCorFlM[3].SetPtEtaPhiM(theJet.pt() * bTcorr, theJet.eta(), theJet.phi(), bMass);
149  jetsCorMix[3].SetPtEtaPhiM(theJet.pt() * tTcorr, theJet.eta(), theJet.phi(), 0);
150  isQuarkFound[3] = true;
151  } else if (particPDG == 2) { //W+ from t
152  double qTcorr = qTopCorrector->correction(theJet);
153  double tTcorr = tTopCorrector->correction(theJet);
154  jetsNoCorr[1].SetPtEtaPhiM(theJet.pt(), theJet.eta(), theJet.phi(), theJet.mass());
155  jetsCorFl0[1].SetPtEtaPhiM(theJet.pt() * qTcorr, theJet.eta(), theJet.phi(), 0);
156  jetsCorFlM[1].SetPtEtaPhiM(theJet.pt() * qTcorr, theJet.eta(), theJet.phi(), qMass);
157  jetsCorMix[1].SetPtEtaPhiM(theJet.pt() * tTcorr, theJet.eta(), theJet.phi(), 0);
158  isQuarkFound[1] = true;
159  } else if (particPDG == 4) { //W+ from t
160  double cTcorr = cTopCorrector->correction(theJet);
161  double tTcorr = tTopCorrector->correction(theJet);
162  jetsNoCorr[1].SetPtEtaPhiM(theJet.pt(), theJet.eta(), theJet.phi(), theJet.mass());
163  jetsCorFl0[1].SetPtEtaPhiM(theJet.pt() * cTcorr, theJet.eta(), theJet.phi(), 0);
164  jetsCorFlM[1].SetPtEtaPhiM(theJet.pt() * cTcorr, theJet.eta(), theJet.phi(), cMass);
165  jetsCorMix[1].SetPtEtaPhiM(theJet.pt() * tTcorr, theJet.eta(), theJet.phi(), 0);
166  isQuarkFound[1] = true;
167  } else if (particPDG == -1) { //W+ from t
168  double qTcorr = qTopCorrector->correction(theJet);
169  double tTcorr = tTopCorrector->correction(theJet);
170  jetsNoCorr[2].SetPtEtaPhiM(theJet.pt(), theJet.eta(), theJet.phi(), theJet.mass());
171  jetsCorFl0[2].SetPtEtaPhiM(theJet.pt() * qTcorr, theJet.eta(), theJet.phi(), 0);
172  jetsCorFlM[2].SetPtEtaPhiM(theJet.pt() * qTcorr, theJet.eta(), theJet.phi(), qMass);
173  jetsCorMix[2].SetPtEtaPhiM(theJet.pt() * tTcorr, theJet.eta(), theJet.phi(), 0);
174  isQuarkFound[2] = true;
175  } else if (particPDG == -3) { //W+ from t
176  double qTcorr = qTopCorrector->correction(theJet);
177  double tTcorr = tTopCorrector->correction(theJet);
178  jetsNoCorr[2].SetPtEtaPhiM(theJet.pt(), theJet.eta(), theJet.phi(), theJet.mass());
179  jetsCorFl0[2].SetPtEtaPhiM(theJet.pt() * qTcorr, theJet.eta(), theJet.phi(), 0);
180  jetsCorFlM[2].SetPtEtaPhiM(theJet.pt() * qTcorr, theJet.eta(), theJet.phi(), qMass);
181  jetsCorMix[2].SetPtEtaPhiM(theJet.pt() * tTcorr, theJet.eta(), theJet.phi(), 0);
182  isQuarkFound[2] = true;
183  } else if (particPDG == -2) { //W- from tbar
184  double qTcorr = qTopCorrector->correction(theJet);
185  double tTcorr = tTopCorrector->correction(theJet);
186  jetsNoCorr[4].SetPtEtaPhiM(theJet.pt(), theJet.eta(), theJet.phi(), theJet.mass());
187  jetsCorFl0[4].SetPtEtaPhiM(theJet.pt() * qTcorr, theJet.eta(), theJet.phi(), 0);
188  jetsCorFlM[4].SetPtEtaPhiM(theJet.pt() * qTcorr, theJet.eta(), theJet.phi(), qMass);
189  jetsCorMix[4].SetPtEtaPhiM(theJet.pt() * tTcorr, theJet.eta(), theJet.phi(), 0);
190  isQuarkFound[4] = true;
191  } else if (particPDG == -4) { //W- from tbar
192  double qTcorr = qTopCorrector->correction(theJet);
193  double tTcorr = tTopCorrector->correction(theJet);
194  jetsNoCorr[4].SetPtEtaPhiM(theJet.pt(), theJet.eta(), theJet.phi(), theJet.mass());
195  jetsCorFl0[4].SetPtEtaPhiM(theJet.pt() * qTcorr, theJet.eta(), theJet.phi(), 0);
196  jetsCorFlM[4].SetPtEtaPhiM(theJet.pt() * qTcorr, theJet.eta(), theJet.phi(), qMass);
197  jetsCorMix[4].SetPtEtaPhiM(theJet.pt() * tTcorr, theJet.eta(), theJet.phi(), 0);
198  isQuarkFound[4] = true;
199  } else if (particPDG == 1) { //W- from tbar
200  double qTcorr = qTopCorrector->correction(theJet);
201  double tTcorr = tTopCorrector->correction(theJet);
202  jetsNoCorr[5].SetPtEtaPhiM(theJet.pt(), theJet.eta(), theJet.phi(), theJet.mass());
203  jetsCorFl0[5].SetPtEtaPhiM(theJet.pt() * qTcorr, theJet.eta(), theJet.phi(), 0);
204  jetsCorFlM[5].SetPtEtaPhiM(theJet.pt() * qTcorr, theJet.eta(), theJet.phi(), qMass);
205  jetsCorMix[5].SetPtEtaPhiM(theJet.pt() * tTcorr, theJet.eta(), theJet.phi(), 0);
206  isQuarkFound[5] = true;
207  } else if (particPDG == 3) { //W- from tbar
208  double cTcorr = cTopCorrector->correction(theJet);
209  double tTcorr = tTopCorrector->correction(theJet);
210  jetsNoCorr[5].SetPtEtaPhiM(theJet.pt(), theJet.eta(), theJet.phi(), theJet.mass());
211  jetsCorFl0[5].SetPtEtaPhiM(theJet.pt() * cTcorr, theJet.eta(), theJet.phi(), 0);
212  jetsCorFlM[5].SetPtEtaPhiM(theJet.pt() * cTcorr, theJet.eta(), theJet.phi(), cMass);
213  jetsCorMix[5].SetPtEtaPhiM(theJet.pt() * tTcorr, theJet.eta(), theJet.phi(), 0);
214  isQuarkFound[5] = true;
215  }
216  }
217  }
218 
219  if (isQuarkFound[0] && isQuarkFound[1] && isQuarkFound[2]) {
220  TLorentzVector *theTopPNoCorr = new TLorentzVector(jetsNoCorr[0] + jetsNoCorr[1] + jetsNoCorr[2]);
221  TLorentzVector *theTopPCorFl0 = new TLorentzVector(jetsCorFl0[0] + jetsCorFl0[1] + jetsCorFl0[2]);
222  TLorentzVector *theTopPCorFlM = new TLorentzVector(jetsCorFlM[0] + jetsCorFlM[1] + jetsCorFlM[2]);
223  TLorentzVector *theTopPCorMix = new TLorentzVector(jetsCorMix[0] + jetsCorMix[1] + jetsCorMix[2]);
224  hMassNoCorr->Fill(theTopPNoCorr->M());
225  hMassCorFl0->Fill(theTopPCorFl0->M());
226  hMassCorFlM->Fill(theTopPCorFlM->M());
227  hMassCorMix->Fill(theTopPCorMix->M());
228  }
229 
230  if (isQuarkFound[3] && isQuarkFound[4] && isQuarkFound[5]) {
231  TLorentzVector *theTopMNoCorr = new TLorentzVector(jetsNoCorr[3] + jetsNoCorr[4] + jetsNoCorr[5]);
232  TLorentzVector *theTopMCorFl0 = new TLorentzVector(jetsCorFl0[3] + jetsCorFl0[4] + jetsCorFl0[5]);
233  TLorentzVector *theTopMCorFlM = new TLorentzVector(jetsCorFlM[3] + jetsCorFlM[4] + jetsCorFlM[5]);
234  TLorentzVector *theTopMCorMix = new TLorentzVector(jetsCorMix[3] + jetsCorMix[4] + jetsCorMix[5]);
235  hMassNoCorr->Fill(theTopMNoCorr->M());
236  hMassCorFl0->Fill(theTopMCorFl0->M());
237  hMassCorFlM->Fill(theTopMCorFlM->M());
238  hMassCorMix->Fill(theTopMCorMix->M());
239  }
240 }
241 
static const std::string kSharedResource
Definition: TFileService.h:76
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
edm::Handle< reco::JetMatchedPartonsCollection > theJetPartonMatch
std::string m_cT_CorrectorName
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
edm::InputTag sourceByRefer_
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:238
double correction(const LorentzVector &fJet) const
get correction using Jet information only
Definition: JetCorrector.h:46
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
int iEvent
Definition: GenABIO.cc:224
void analyze(const edm::Event &iEvent, const edm::EventSetup &iSetup) override
~calcTopMass() override
const GenParticleRef & physicsDefinitionParton() const
calcTopMass(const edm::ParameterSet &)
std::string m_qT_CorrectorName
static const JetCorrector * getJetCorrector(const std::string &fName, const edm::EventSetup &fSetup)
retrieve corrector from the event setup. troughs exception if something is missing ...
Definition: JetCorrector.cc:48
std::string m_tT_CorrectorName
fixed size matrix
HLT enums.
edm::InputTag sourcePartons_
T const * get() const
Returns C++ pointer to the item.
Definition: Ref.h:232
std::string m_bT_CorrectorName