Go to the documentation of this file.00001
00002
00003
00004
00005
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022 #include <memory>
00023
00024
00025 #include "FWCore/Framework/interface/Frameworkfwd.h"
00026 #include "FWCore/Framework/interface/EDProducer.h"
00027
00028 #include "FWCore/Framework/interface/Event.h"
00029 #include "FWCore/Framework/interface/MakerMacros.h"
00030
00031 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00032
00033 #include "RecoJets/JetPlusTracks/plugins/JetPlusTrackProducerAA.h"
00034 #include "DataFormats/JetReco/interface/CaloJetCollection.h"
00035 #include "DataFormats/JetReco/interface/CaloJet.h"
00036 #include "DataFormats/JetReco/interface/JPTJetCollection.h"
00037 #include "DataFormats/JetReco/interface/JPTJet.h"
00038 #include "DataFormats/TrackReco/interface/TrackFwd.h"
00039 #include "DataFormats/TrackReco/interface/Track.h"
00040 #include "DataFormats/JetReco/interface/Jet.h"
00041 #include "DataFormats/VertexReco/interface/VertexFwd.h"
00042 #include "DataFormats/VertexReco/interface/Vertex.h"
00043 #include "DataFormats/Math/interface/deltaPhi.h"
00044 #include "DataFormats/Math/interface/deltaR.h"
00045
00046 #include <string>
00047
00048 using namespace std;
00049 using namespace jpt;
00050
00051
00052
00053
00054
00055
00056
00057
00058
00059
00060
00061
00062
00063 JetPlusTrackProducerAA::JetPlusTrackProducerAA(const edm::ParameterSet& iConfig)
00064 {
00065
00066 src = iConfig.getParameter<edm::InputTag>("src");
00067 alias = iConfig.getUntrackedParameter<string>("alias");
00068 mTracks = iConfig.getParameter<edm::InputTag> ("tracks");
00069 srcPVs_ = iConfig.getParameter<edm::InputTag>("srcPVs");
00070 vectorial_ = iConfig.getParameter<bool>("VectorialCorrection");
00071 useZSP = iConfig.getParameter<bool>("UseZSP");
00072 std::string tq = iConfig.getParameter<std::string>("TrackQuality");
00073 trackQuality_ = reco::TrackBase::qualityByName(tq);
00074 mConeSize = iConfig.getParameter<double> ("coneSize");
00075 mJPTalgo = new JetPlusTrackCorrector(iConfig);
00076 mZSPalgo = new ZSPJPTJetCorrector(iConfig);
00077
00078 produces<reco::JPTJetCollection>().setBranchAlias(alias);
00079
00080 }
00081
00082
00083 JetPlusTrackProducerAA::~JetPlusTrackProducerAA()
00084 {
00085
00086
00087
00088
00089 }
00090
00091
00092
00093
00094
00095
00096
00097 void
00098 JetPlusTrackProducerAA::produce(edm::Event& iEvent, const edm::EventSetup& iSetup)
00099 {
00100 using namespace edm;
00101
00102
00103
00104
00105
00106
00107 edm::Handle <edm::View <reco::CaloJet> > jets_h;
00108 iEvent.getByLabel (src, jets_h);
00109
00110 edm::Handle <reco::TrackCollection> tracks_h;
00111 iEvent.getByLabel (mTracks, tracks_h);
00112
00113 std::vector <reco::TrackRef> fTracks;
00114 fTracks.reserve (tracks_h->size());
00115 for (unsigned i = 0; i < tracks_h->size(); ++i) {
00116 fTracks.push_back (reco::TrackRef (tracks_h, i));
00117 }
00118
00119
00120 std::auto_ptr<reco::JPTJetCollection> pOut(new reco::JPTJetCollection());
00121
00122 reco::JPTJetCollection tmpColl;
00123
00124 for (unsigned i = 0; i < jets_h->size(); ++i) {
00125
00126 const reco::CaloJet* oldjet = &(*(jets_h->refAt(i)));
00127
00128 reco::CaloJet corrected = *oldjet;
00129
00130
00131
00132 double factorZSP = 1.;
00133 if(useZSP) factorZSP = mZSPalgo->correction(corrected, iEvent, iSetup);
00134
00135 corrected.scaleEnergy (factorZSP);
00136
00137
00138
00139 double scaleJPT = 1.;
00140
00141 math::XYZTLorentzVector p4;
00142
00143 if ( !vectorial_ ) {
00144
00145 scaleJPT = mJPTalgo->correction ( corrected, *oldjet, iEvent, iSetup );
00146 p4 = math::XYZTLorentzVector( corrected.px()*scaleJPT,
00147 corrected.py()*scaleJPT,
00148 corrected.pz()*scaleJPT,
00149 corrected.energy()*scaleJPT );
00150 } else {
00151 scaleJPT = mJPTalgo->correction( corrected, *oldjet, iEvent, iSetup, p4 );
00152 }
00153
00154
00155 jpt::MatchedTracks pions;
00156 jpt::MatchedTracks muons;
00157 jpt::MatchedTracks elecs;
00158
00159 bool ok = mJPTalgo->matchTracks( *oldjet,
00160 iEvent,
00161 iSetup,
00162 pions,
00163 muons,
00164 elecs );
00165
00166
00167 reco::JPTJet::Specific specific;
00168
00169 if(ok) {
00170 specific.pionsInVertexInCalo = pions.inVertexInCalo_;
00171 specific.pionsInVertexOutCalo = pions.inVertexOutOfCalo_;
00172 specific.pionsOutVertexInCalo = pions.outOfVertexInCalo_;
00173 specific.muonsInVertexInCalo = muons.inVertexInCalo_;
00174 specific.muonsInVertexOutCalo = muons.inVertexOutOfCalo_;
00175 specific.muonsOutVertexInCalo = muons.outOfVertexInCalo_;
00176 specific.elecsInVertexInCalo = elecs.inVertexInCalo_;
00177 specific.elecsInVertexOutCalo = elecs.inVertexOutOfCalo_;
00178 specific.elecsOutVertexInCalo = elecs.outOfVertexInCalo_;
00179 }
00180
00181
00182 edm::RefToBase<reco::Jet> myjet = (edm::RefToBase<reco::Jet>)jets_h->refAt(i);
00183 specific.theCaloJetRef = myjet;
00184 specific.mZSPCor = factorZSP;
00185 specific.mResponseOfChargedWithEff = (float)mJPTalgo->getResponseOfChargedWithEff();
00186 specific.mResponseOfChargedWithoutEff = (float)mJPTalgo->getResponseOfChargedWithoutEff();
00187 specific.mSumPtOfChargedWithEff = (float)mJPTalgo->getSumPtWithEff();
00188 specific.mSumPtOfChargedWithoutEff = (float)mJPTalgo->getSumPtWithoutEff();
00189 specific.mSumEnergyOfChargedWithEff = (float)mJPTalgo->getSumEnergyWithEff();
00190 specific.mSumEnergyOfChargedWithoutEff = (float)mJPTalgo->getSumEnergyWithoutEff();
00191 specific.mChargedHadronEnergy = (float)mJPTalgo->getSumEnergyWithoutEff();
00192
00193 double deR2Tr = 0.;
00194 double deEta2Tr = 0.;
00195 double dePhi2Tr = 0.;
00196 double Zch = 0.;
00197 double Pout2 = 0.;
00198 double Pout = 0.;
00199 double denominator_tracks = 0.;
00200 int ntracks = 0;
00201
00202 for( reco::TrackRefVector::const_iterator it = pions.inVertexInCalo_.begin(); it != pions.inVertexInCalo_.end(); it++) {
00203 double deR = deltaR((*it)->eta(), (*it)->phi(), p4.eta(), p4.phi());
00204 double deEta = (*it)->eta() - p4.eta();
00205 double dePhi = deltaPhi((*it)->phi(), p4.phi());
00206 if((**it).ptError()/(**it).pt() < 0.1) {
00207 deR2Tr = deR2Tr + deR*deR*(*it)->pt();
00208 deEta2Tr = deEta2Tr + deEta*deEta*(*it)->pt();
00209 dePhi2Tr = dePhi2Tr + dePhi*dePhi*(*it)->pt();
00210 denominator_tracks = denominator_tracks + (*it)->pt();
00211 Zch = Zch +
00212 ((*it)->px()*p4.Px()+(*it)->py()*p4.Py()+(*it)->pz()*p4.Pz())/(p4.P()*p4.P());
00213 Pout2 = Pout2 + (**it).p()*(**it).p() - (Zch*p4.P())*(Zch*p4.P());
00214 ntracks++;
00215 }
00216 }
00217 for( reco::TrackRefVector::const_iterator it = muons.inVertexInCalo_.begin(); it != muons.inVertexInCalo_.end(); it++) {
00218 double deR = deltaR((*it)->eta(), (*it)->phi(), p4.eta(), p4.phi());
00219 double deEta = (*it)->eta() - p4.eta();
00220 double dePhi = deltaPhi((*it)->phi(), p4.phi());
00221 if((**it).ptError()/(**it).pt() < 0.1) {
00222 deR2Tr = deR2Tr + deR*deR*(*it)->pt();
00223 deEta2Tr = deEta2Tr + deEta*deEta*(*it)->pt();
00224 dePhi2Tr = dePhi2Tr + dePhi*dePhi*(*it)->pt();
00225 denominator_tracks = denominator_tracks + (*it)->pt();
00226 Zch = Zch +
00227 ((*it)->px()*p4.Px()+(*it)->py()*p4.Py()+(*it)->pz()*p4.Pz())/(p4.P()*p4.P());
00228 Pout2 = Pout2 + (**it).p()*(**it).p() - (Zch*p4.P())*(Zch*p4.P());
00229 ntracks++;
00230 }
00231 }
00232 for( reco::TrackRefVector::const_iterator it = elecs.inVertexInCalo_.begin(); it != elecs.inVertexInCalo_.end(); it++) {
00233 double deR = deltaR((*it)->eta(), (*it)->phi(), p4.eta(), p4.phi());
00234 double deEta = (*it)->eta() - p4.eta();
00235 double dePhi = deltaPhi((*it)->phi(), p4.phi());
00236 if((**it).ptError()/(**it).pt() < 0.1) {
00237 deR2Tr = deR2Tr + deR*deR*(*it)->pt();
00238 deEta2Tr = deEta2Tr + deEta*deEta*(*it)->pt();
00239 dePhi2Tr = dePhi2Tr + dePhi*dePhi*(*it)->pt();
00240 denominator_tracks = denominator_tracks + (*it)->pt();
00241 Zch = Zch +
00242 ((*it)->px()*p4.Px()+(*it)->py()*p4.Py()+(*it)->pz()*p4.Pz())/(p4.P()*p4.P());
00243 Pout2 = Pout2 + (**it).p()*(**it).p() - (Zch*p4.P())*(Zch*p4.P());
00244 ntracks++;
00245 }
00246 }
00247
00248 if(ntracks > 0) {
00249 Pout = sqrt(fabs(Pout2))/ntracks;
00250 Zch = Zch/ntracks;
00251 }
00252
00253
00254 if (denominator_tracks!=0){
00255 deR2Tr = deR2Tr/denominator_tracks;
00256 deEta2Tr= deEta2Tr/denominator_tracks;
00257 dePhi2Tr= dePhi2Tr/denominator_tracks;
00258 }
00259
00260 specific.R2momtr = deR2Tr;
00261 specific.Eta2momtr = deEta2Tr;
00262 specific.Phi2momtr = dePhi2Tr;
00263 specific.Pout = Pout;
00264 specific.Zch = Zch;
00265
00266
00267
00268 reco::Particle::Point vertex_=reco::Jet::Point(0,0,0);
00269
00270
00271 edm::Handle<reco::VertexCollection> pvCollection;
00272 iEvent.getByLabel(srcPVs_,pvCollection);
00273 if ( pvCollection.isValid() && pvCollection->size()>0 ) vertex_=pvCollection->begin()->position();
00274
00275 reco::JPTJet fJet(p4, vertex_, specific, corrected.getJetConstituents());
00276
00277
00278
00279
00280 tmpColl.push_back(fJet);
00281
00282 }
00283
00284
00285 reco::TrackRefVector trBgOutOfVertex = calculateBGtracksJet(tmpColl,fTracks);
00286
00287 for(reco::JPTJetCollection::iterator ij=tmpColl.begin(); ij!=tmpColl.end(); ij++)
00288 {
00289
00290 double factorPU = mJPTalgo->correctAA(*ij,trBgOutOfVertex,mConeSize);
00291 (*ij).scaleEnergy (factorPU);
00292
00293 pOut->push_back(*ij);
00294 }
00295 iEvent.put(pOut);
00296
00297 }
00298
00299
00300
00301
00302 reco::TrackRefVector JetPlusTrackProducerAA::calculateBGtracksJet(reco::JPTJetCollection& fJets, std::vector <reco::TrackRef>& fTracks){
00303
00304 reco::TrackRefVector trBgOutOfVertex;
00305
00306 for (unsigned t = 0; t < fTracks.size(); ++t) {
00307
00308 int track_bg = 0;
00309
00310
00311
00312
00313
00314
00315
00316 const reco::Track* track = &*(fTracks[t]);
00317 double trackEta = track->eta();
00318 double trackPhi = track->phi();
00319
00320
00321
00322
00323
00324 for (unsigned j = 0; j < fJets.size(); ++j) {
00325
00326 const reco::Jet* jet = &(fJets[j]);
00327 double jetEta = jet->eta();
00328 double jetPhi = jet->phi();
00329
00330
00331
00332
00333 if(fabs(jetEta - trackEta) < mConeSize) {
00334 double dphiTrackJet = fabs(trackPhi - jetPhi);
00335 if(dphiTrackJet > M_PI) dphiTrackJet = 2*M_PI - dphiTrackJet;
00336
00337 if(dphiTrackJet < mConeSize)
00338 {
00339 track_bg = 1;
00340
00341
00342
00343 }
00344 }
00345 }
00346
00347 if( track_bg == 0 ) trBgOutOfVertex.push_back (fTracks[t]);
00348
00349
00350
00351
00352
00353 }
00354
00355 return trBgOutOfVertex;
00356 }
00357
00358 void
00359 JetPlusTrackProducerAA::beginJob()
00360 {
00361 }
00362
00363
00364 void
00365 JetPlusTrackProducerAA::endJob() {
00366 }
00367
00368
00369