45 #include "TLorentzVector.h"
53 #include <Math/VectorUtil.h>
84 using namespace ROOT::Math::VectorUtil;
100 hMassNoCorr = fs->
make<TH1F>(
"hMassNoCorr",
"",100,100,300);
101 hMassCorFl0 = fs->
make<TH1F>(
"hMassCorFl0",
"",100,100,300);
102 hMassCorFlM = fs->
make<TH1F>(
"hMassCorFlM",
"",100,100,300);
103 hMassCorMix = fs->
make<TH1F>(
"hMassCorMix",
"",100,100,300);
109 cout <<
"[calcTopMass] analysing event " << iEvent.
id() << endl;
111 iEvent.
getByLabel (sourceByRefer_ , theJetPartonMatch );
113 cerr <<
"[calcTopMass] caught std::exception " << ce.what() << endl;
123 TLorentzVector jetsNoCorr[6];
124 TLorentzVector jetsCorFl0[6];
125 TLorentzVector jetsCorFlM[6];
126 TLorentzVector jetsCorMix[6];
128 bool isQuarkFound[6] = {
false};
130 for ( JetMatchedPartonsCollection::const_iterator
j = theJetPartonMatch->begin();
131 j != theJetPartonMatch->end();
140 int particPDG = thePhyDef.
get()->pdgId();
142 if( particPDG == 5 ) {
143 double bTcorr = bTopCorrector->
correction (theJet);
144 double tTcorr = tTopCorrector->
correction (theJet);
145 jetsNoCorr[0].SetPtEtaPhiM( theJet.pt() , theJet.eta() , theJet.phi() , theJet.mass() );
146 jetsCorFl0[0].SetPtEtaPhiM( theJet.pt()*bTcorr , theJet.eta() , theJet.phi() , 0 );
147 jetsCorFlM[0].SetPtEtaPhiM( theJet.pt()*bTcorr , theJet.eta() , theJet.phi() , bMass );
148 jetsCorMix[0].SetPtEtaPhiM( theJet.pt()*tTcorr , theJet.eta() , theJet.phi() , 0 );
149 isQuarkFound[0]=
true;
150 }
else if ( particPDG == -5 ) {
151 double bTcorr = bTopCorrector->
correction (theJet);
152 double tTcorr = tTopCorrector->
correction (theJet);
153 jetsNoCorr[3].SetPtEtaPhiM( theJet.pt() , theJet.eta() , theJet.phi() , theJet.mass() );
154 jetsCorFl0[3].SetPtEtaPhiM( theJet.pt()*bTcorr , theJet.eta() , theJet.phi() , 0 );
155 jetsCorFlM[3].SetPtEtaPhiM( theJet.pt()*bTcorr , theJet.eta() , theJet.phi() , bMass );
156 jetsCorMix[3].SetPtEtaPhiM( theJet.pt()*tTcorr , theJet.eta() , theJet.phi() , 0 );
157 isQuarkFound[3]=
true;
158 }
else if ( particPDG == 2 ) {
159 double qTcorr = qTopCorrector->
correction (theJet);
160 double tTcorr = tTopCorrector->
correction (theJet);
161 jetsNoCorr[1].SetPtEtaPhiM( theJet.pt() , theJet.eta() , theJet.phi() , theJet.mass() );
162 jetsCorFl0[1].SetPtEtaPhiM( theJet.pt()*qTcorr , theJet.eta() , theJet.phi() , 0 );
163 jetsCorFlM[1].SetPtEtaPhiM( theJet.pt()*qTcorr , theJet.eta() , theJet.phi() , qMass );
164 jetsCorMix[1].SetPtEtaPhiM( theJet.pt()*tTcorr , theJet.eta() , theJet.phi() , 0 );
165 isQuarkFound[1]=
true;
166 }
else if ( particPDG == 4 ) {
167 double cTcorr = cTopCorrector->
correction (theJet);
168 double tTcorr = tTopCorrector->
correction (theJet);
169 jetsNoCorr[1].SetPtEtaPhiM( theJet.pt() , theJet.eta() , theJet.phi() , theJet.mass() );
170 jetsCorFl0[1].SetPtEtaPhiM( theJet.pt()*cTcorr , theJet.eta() , theJet.phi() , 0 );
171 jetsCorFlM[1].SetPtEtaPhiM( theJet.pt()*cTcorr , theJet.eta() , theJet.phi() , cMass );
172 jetsCorMix[1].SetPtEtaPhiM( theJet.pt()*tTcorr , theJet.eta() , theJet.phi() , 0 );
173 isQuarkFound[1]=
true;
174 }
else if ( particPDG == -1 ) {
175 double qTcorr = qTopCorrector->
correction (theJet);
176 double tTcorr = tTopCorrector->
correction (theJet);
177 jetsNoCorr[2].SetPtEtaPhiM( theJet.pt() , theJet.eta() , theJet.phi() , theJet.mass() );
178 jetsCorFl0[2].SetPtEtaPhiM( theJet.pt()*qTcorr , theJet.eta() , theJet.phi() , 0 );
179 jetsCorFlM[2].SetPtEtaPhiM( theJet.pt()*qTcorr , theJet.eta() , theJet.phi() , qMass );
180 jetsCorMix[2].SetPtEtaPhiM( theJet.pt()*tTcorr , theJet.eta() , theJet.phi() , 0 );
181 isQuarkFound[2]=
true;
182 }
else if ( particPDG == -3 ) {
183 double qTcorr = qTopCorrector->
correction (theJet);
184 double tTcorr = tTopCorrector->
correction (theJet);
185 jetsNoCorr[2].SetPtEtaPhiM( theJet.pt() , theJet.eta() , theJet.phi() , theJet.mass() );
186 jetsCorFl0[2].SetPtEtaPhiM( theJet.pt()*qTcorr , theJet.eta() , theJet.phi() , 0 );
187 jetsCorFlM[2].SetPtEtaPhiM( theJet.pt()*qTcorr , theJet.eta() , theJet.phi() , qMass );
188 jetsCorMix[2].SetPtEtaPhiM( theJet.pt()*tTcorr , theJet.eta() , theJet.phi() , 0 );
189 isQuarkFound[2]=
true;
190 }
else if ( particPDG == -2 ) {
191 double qTcorr = qTopCorrector->
correction (theJet);
192 double tTcorr = tTopCorrector->
correction (theJet);
193 jetsNoCorr[4].SetPtEtaPhiM( theJet.pt() , theJet.eta() , theJet.phi() , theJet.mass() );
194 jetsCorFl0[4].SetPtEtaPhiM( theJet.pt()*qTcorr , theJet.eta() , theJet.phi() , 0 );
195 jetsCorFlM[4].SetPtEtaPhiM( theJet.pt()*qTcorr , theJet.eta() , theJet.phi() , qMass );
196 jetsCorMix[4].SetPtEtaPhiM( theJet.pt()*tTcorr , theJet.eta() , theJet.phi() , 0 );
197 isQuarkFound[4]=
true;
198 }
else if ( particPDG == -4 ) {
199 double qTcorr = qTopCorrector->
correction (theJet);
200 double tTcorr = tTopCorrector->
correction (theJet);
201 jetsNoCorr[4].SetPtEtaPhiM( theJet.pt() , theJet.eta() , theJet.phi() , theJet.mass() );
202 jetsCorFl0[4].SetPtEtaPhiM( theJet.pt()*qTcorr , theJet.eta() , theJet.phi() , 0 );
203 jetsCorFlM[4].SetPtEtaPhiM( theJet.pt()*qTcorr , theJet.eta() , theJet.phi() , qMass );
204 jetsCorMix[4].SetPtEtaPhiM( theJet.pt()*tTcorr , theJet.eta() , theJet.phi() , 0 );
205 isQuarkFound[4]=
true;
206 }
else if ( particPDG == 1 ) {
207 double qTcorr = qTopCorrector->
correction (theJet);
208 double tTcorr = tTopCorrector->
correction (theJet);
209 jetsNoCorr[5].SetPtEtaPhiM( theJet.pt() , theJet.eta() , theJet.phi() , theJet.mass() );
210 jetsCorFl0[5].SetPtEtaPhiM( theJet.pt()*qTcorr , theJet.eta() , theJet.phi() , 0 );
211 jetsCorFlM[5].SetPtEtaPhiM( theJet.pt()*qTcorr , theJet.eta() , theJet.phi() , qMass );
212 jetsCorMix[5].SetPtEtaPhiM( theJet.pt()*tTcorr , theJet.eta() , theJet.phi() , 0 );
213 isQuarkFound[5]=
true;
214 }
else if ( particPDG == 3 ) {
215 double cTcorr = cTopCorrector->
correction (theJet);
216 double tTcorr = tTopCorrector->
correction (theJet);
217 jetsNoCorr[5].SetPtEtaPhiM( theJet.pt() , theJet.eta() , theJet.phi() , theJet.mass() );
218 jetsCorFl0[5].SetPtEtaPhiM( theJet.pt()*cTcorr , theJet.eta() , theJet.phi() , 0 );
219 jetsCorFlM[5].SetPtEtaPhiM( theJet.pt()*cTcorr , theJet.eta() , theJet.phi() , cMass );
220 jetsCorMix[5].SetPtEtaPhiM( theJet.pt()*tTcorr , theJet.eta() , theJet.phi() , 0 );
221 isQuarkFound[5]=
true;
226 if( isQuarkFound[0] && isQuarkFound[1] && isQuarkFound[2] ) {
227 TLorentzVector *theTopPNoCorr =
new TLorentzVector( jetsNoCorr[0] + jetsNoCorr[1] + jetsNoCorr[2] );
228 TLorentzVector *theTopPCorFl0 =
new TLorentzVector( jetsCorFl0[0] + jetsCorFl0[1] + jetsCorFl0[2] );
229 TLorentzVector *theTopPCorFlM =
new TLorentzVector( jetsCorFlM[0] + jetsCorFlM[1] + jetsCorFlM[2] );
230 TLorentzVector *theTopPCorMix =
new TLorentzVector( jetsCorMix[0] + jetsCorMix[1] + jetsCorMix[2] );
231 hMassNoCorr->Fill( theTopPNoCorr->M() );
232 hMassCorFl0->Fill( theTopPCorFl0->M() );
233 hMassCorFlM->Fill( theTopPCorFlM->M() );
234 hMassCorMix->Fill( theTopPCorMix->M() );
237 if( isQuarkFound[3] && isQuarkFound[4] && isQuarkFound[5] ) {
238 TLorentzVector *theTopMNoCorr =
new TLorentzVector( jetsNoCorr[3] + jetsNoCorr[4] + jetsNoCorr[5] );
239 TLorentzVector *theTopMCorFl0 =
new TLorentzVector( jetsCorFl0[3] + jetsCorFl0[4] + jetsCorFl0[5] );
240 TLorentzVector *theTopMCorFlM =
new TLorentzVector( jetsCorFlM[3] + jetsCorFlM[4] + jetsCorFlM[5] );
241 TLorentzVector *theTopMCorMix =
new TLorentzVector( jetsCorMix[3] + jetsCorMix[4] + jetsCorMix[5] );
242 hMassNoCorr->Fill( theTopMNoCorr->M() );
243 hMassCorFl0->Fill( theTopMCorFl0->M() );
244 hMassCorFlM->Fill( theTopMCorFlM->M() );
245 hMassCorMix->Fill( theTopMCorMix->M() );
T getParameter(std::string const &) const
edm::Handle< reco::JetMatchedPartonsCollection > theJetPartonMatch
std::string m_cT_CorrectorName
virtual double correction(const LorentzVector &fJet) const =0
get correction using Jet information only
#define DEFINE_FWK_MODULE(type)
edm::InputTag sourceByRefer_
T * make(const Args &...args) const
make new ROOT object
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
bool isNonnull() const
Checks for non-null.
virtual void analyze(const edm::Event &iEvent, const edm::EventSetup &iSetup) override
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
calcTopMass(const edm::ParameterSet &)
std::string m_qT_CorrectorName
const GenParticleRef & physicsDefinitionParton() const
static const JetCorrector * getJetCorrector(const std::string &fName, const edm::EventSetup &fSetup)
retrieve corrector from the event setup. troughs exception if something is missing ...
std::string m_tT_CorrectorName
edm::InputTag sourcePartons_
T const * get() const
Returns C++ pointer to the item.
std::string m_bT_CorrectorName