CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_4/src/Calibration/HcalAlCaRecoProducers/src/AlCaGammaJetProducer.cc

Go to the documentation of this file.
00001 #include "Calibration/HcalAlCaRecoProducers/interface/AlCaGammaJetProducer.h"
00002 #include "DataFormats/DetId/interface/DetId.h"
00003 #include "Geometry/Records/interface/CaloGeometryRecord.h"
00004 #include "DataFormats/GeometryVector/interface/GlobalPoint.h"
00005 #include "DataFormats/HcalRecHit/interface/HcalRecHitCollections.h"
00006 #include "DataFormats/EcalRecHit/interface/EcalRecHitCollections.h"
00007 #include "TrackingTools/TransientTrack/interface/TransientTrack.h"
00008 #include "DataFormats/JetReco/interface/CaloJetCollection.h"
00009 #include "DataFormats/EgammaReco/interface/SuperCluster.h"
00010 #include "DataFormats/EgammaReco/interface/SuperClusterFwd.h"
00011 
00012 using namespace edm;
00013 using namespace std;
00014 using namespace reco;
00015 
00016 //namespace cms
00017 //{
00018 
00019 AlCaGammaJetProducer::AlCaGammaJetProducer(const edm::ParameterSet& iConfig)
00020 {
00021    // Take input 
00022    
00023    hbheLabel_= iConfig.getParameter<edm::InputTag>("hbheInput");
00024    hoLabel_=iConfig.getParameter<edm::InputTag>("hoInput");
00025    hfLabel_=iConfig.getParameter<edm::InputTag>("hfInput");
00026    mInputCalo = iConfig.getParameter<std::vector<edm::InputTag> >("srcCalo");
00027    ecalLabels_=iConfig.getParameter<std::vector<edm::InputTag> >("ecalInputs");
00028    m_inputTrackLabel = iConfig.getUntrackedParameter<std::string>("inputTrackLabel","generalTracks");
00029    correctedIslandBarrelSuperClusterCollection_ = iConfig.getParameter<std::string>("correctedIslandBarrelSuperClusterCollection");
00030    correctedIslandBarrelSuperClusterProducer_   = iConfig.getParameter<std::string>("correctedIslandBarrelSuperClusterProducer");
00031    correctedIslandEndcapSuperClusterCollection_ = iConfig.getParameter<std::string>("correctedIslandEndcapSuperClusterCollection");
00032    correctedIslandEndcapSuperClusterProducer_   = iConfig.getParameter<std::string>("correctedIslandEndcapSuperClusterProducer");  
00033    allowMissingInputs_=iConfig.getUntrackedParameter<bool>("AllowMissingInputs",true);
00034    
00035     
00036    //register your products
00037    produces<reco::TrackCollection>("GammaJetTracksCollection");
00038    produces<EcalRecHitCollection>("GammaJetEcalRecHitCollection");
00039    produces<HBHERecHitCollection>("GammaJetHBHERecHitCollection");
00040    produces<HORecHitCollection>("GammaJetHORecHitCollection");
00041    produces<HFRecHitCollection>("GammaJetHFRecHitCollection");
00042    produces<CaloJetCollection>("GammaJetJetBackToBackCollection");
00043    produces<reco::SuperClusterCollection>("GammaJetGammaBackToBackCollection");
00044 
00045     
00046    
00047 }
00048 void AlCaGammaJetProducer::beginJob()
00049 {
00050 }
00051 
00052 AlCaGammaJetProducer::~AlCaGammaJetProducer()
00053 {
00054  
00055 
00056 }
00057 
00058 
00059 // ------------ method called to produce the data  ------------
00060 void
00061 AlCaGammaJetProducer::produce(edm::Event& iEvent, const edm::EventSetup& iSetup)
00062 {
00063    using namespace edm;
00064    using namespace std;
00065 
00066    edm::ESHandle<CaloGeometry> pG;
00067    iSetup.get<CaloGeometryRecord>().get(pG);
00068    geo = pG.product();
00069 
00070 // Produced collections 
00071   
00072   std::auto_ptr<reco::SuperClusterCollection> result (new reco::SuperClusterCollection); //Corrected jets
00073   std::auto_ptr<CaloJetCollection> resultjet (new CaloJetCollection); //Corrected jets
00074   std::auto_ptr<EcalRecHitCollection> miniEcalRecHitCollection(new EcalRecHitCollection);
00075   std::auto_ptr<HBHERecHitCollection> miniHBHERecHitCollection(new HBHERecHitCollection);
00076   std::auto_ptr<HORecHitCollection> miniHORecHitCollection(new HORecHitCollection);
00077   std::auto_ptr<HFRecHitCollection> miniHFRecHitCollection(new HFRecHitCollection);
00078   std::auto_ptr<reco::TrackCollection> outputTColl(new reco::TrackCollection);
00079    
00080    
00081 // Get  Corrected supercluster collection
00082   int nclusb = 0;
00083   int ncluse = 0;
00084   reco::SuperClusterCollection::const_iterator maxclusbarrel;
00085   reco::SuperClusterCollection::const_iterator maxclusendcap;
00086 
00087   double vetmax = -100.;
00088   Handle<reco::SuperClusterCollection> pCorrectedIslandBarrelSuperClusters;
00089   iEvent.getByLabel(correctedIslandBarrelSuperClusterProducer_, 
00090                     correctedIslandBarrelSuperClusterCollection_, 
00091                     pCorrectedIslandBarrelSuperClusters);  
00092   if (!pCorrectedIslandBarrelSuperClusters.isValid()) {
00093     // can't find it!
00094     if (!allowMissingInputs_) {
00095       *pCorrectedIslandBarrelSuperClusters;  // will throw the proper exception
00096     }
00097   } else {
00098     const reco::SuperClusterCollection* correctedIslandBarrelSuperClusters = pCorrectedIslandBarrelSuperClusters.product();
00099     
00100     // loop over the super clusters and find the highest
00101     maxclusbarrel = correctedIslandBarrelSuperClusters->begin();
00102     for(reco::SuperClusterCollection::const_iterator aClus = correctedIslandBarrelSuperClusters->begin();
00103         aClus != correctedIslandBarrelSuperClusters->end(); aClus++) {
00104       double vet = aClus->energy()/cosh(aClus->eta());
00105           cout<<" Barrel supercluster " << vet <<" energy "<<aClus->energy()<<" eta "<<aClus->eta()<<endl;
00106       if(vet>20.) {
00107         if(vet > vetmax)
00108           {
00109             vetmax = vet;
00110             maxclusbarrel = aClus;
00111             nclusb = 1;
00112           }
00113       }
00114     }
00115     
00116   }
00117 
00118   Handle<reco::SuperClusterCollection> pCorrectedIslandEndcapSuperClusters;
00119   iEvent.getByLabel(correctedIslandEndcapSuperClusterProducer_, 
00120                     correctedIslandEndcapSuperClusterCollection_, 
00121                     pCorrectedIslandEndcapSuperClusters);  
00122   if (!pCorrectedIslandEndcapSuperClusters.isValid()) {
00123     // can't find it!
00124     if (!allowMissingInputs_) {
00125       *pCorrectedIslandEndcapSuperClusters;  // will throw the proper exception
00126     }
00127   } else {
00128     const reco::SuperClusterCollection* correctedIslandEndcapSuperClusters = pCorrectedIslandEndcapSuperClusters.product();
00129     
00130     // loop over the super clusters and find the highest
00131     maxclusendcap = correctedIslandEndcapSuperClusters->end();
00132     double vetmaxe = vetmax;
00133     for(reco::SuperClusterCollection::const_iterator aClus = correctedIslandEndcapSuperClusters->begin();
00134         aClus != correctedIslandEndcapSuperClusters->end(); aClus++) {
00135       double vet = aClus->energy()/cosh(aClus->eta());
00136       //   cout<<" Endcap supercluster " << vet <<" energy "<<aClus->energy()<<" eta "<<aClus->eta()<<endl;
00137       if(vet>20.) {
00138         if(vet > vetmaxe)
00139           {
00140             vetmaxe = vet;
00141             maxclusendcap = aClus;
00142             ncluse = 1;
00143           }
00144       }
00145     }
00146   }
00147    
00148   cout<<" Number of gammas "<<nclusb<<" "<<ncluse<<endl;  
00149   
00150   if( nclusb == 0 && ncluse == 0 ) {
00151    
00152   iEvent.put( outputTColl, "GammaJetTracksCollection");
00153   iEvent.put( miniEcalRecHitCollection, "GammaJetEcalRecHitCollection");
00154   iEvent.put( miniHBHERecHitCollection, "GammaJetHBHERecHitCollection");
00155   iEvent.put( miniHORecHitCollection, "GammaJetHORecHitCollection");
00156   iEvent.put( miniHFRecHitCollection, "GammaJetHFRecHitCollection");
00157   iEvent.put( result, "GammaJetGammaBackToBackCollection");
00158   iEvent.put( resultjet, "GammaJetJetBackToBackCollection");
00159   
00160   return;
00161   }
00162 
00163   double phigamma = -100.;
00164   double etagamma = -100.;
00165   if(ncluse == 1)
00166   {
00167     result->push_back(*maxclusendcap);
00168     phigamma = (*maxclusendcap).phi();
00169     etagamma = (*maxclusendcap).eta();
00170     
00171   } else
00172   {
00173     result->push_back(*maxclusbarrel);
00174     phigamma = (*maxclusbarrel).phi();
00175     etagamma = (*maxclusbarrel).eta();
00176   }
00177   
00178 //  cout<<" Size of egamma clusters "<<result->size()<<endl;
00179    
00180 //  
00181 // Jet Collection
00182 // Find jet in the angle ~ +- 170 degrees
00183 //
00184     
00185     double phijet =  -100.;
00186     double etajet = -100.;
00187     double phijet0 =  -100.;
00188     double etajet0 = -100.;
00189     
00190     std::vector<edm::InputTag>::const_iterator ic;
00191     for (ic=mInputCalo.begin(); ic!=mInputCalo.end(); ic++) {
00192 
00193       edm::Handle<reco::CaloJetCollection> jets;
00194       iEvent.getByLabel(*ic, jets);
00195       if (!jets.isValid()) {
00196         // can't find it!
00197         if (!allowMissingInputs_) {
00198           *jets;  // will throw the proper exception
00199         }
00200       } else {
00201        reco::CaloJetCollection::const_iterator jet = jets->begin ();
00202 
00203  //      cout<<" Size of jets "<<jets->size()<<endl;
00204 
00205        if( jets->size() == 0 ) continue;
00206 
00207          int iejet = 0;   
00208          int numjet = 0;
00209  
00210          for (; jet != jets->end (); jet++)
00211          {
00212            phijet0 = (*jet).phi();
00213            etajet0 = (*jet).eta();
00214 
00215 // Only 3 jets are kept
00216            numjet++;
00217            if(numjet > 3) break;
00218 
00219            cout<<" phi,eta "<< phigamma<<" "<< etagamma<<" "<<phijet0<<" "<<etajet0<<endl;
00220 
00221 // Find jet back to gamma
00222            
00223            double dphi = fabs(phigamma-phijet0); 
00224            if(dphi > 4.*atan(1.)) dphi = 8.*atan(1.) - dphi;
00225            double deta = fabs(etagamma-etajet0);
00226            double dr = sqrt(dphi*dphi+deta*deta);
00227            if(dr < 0.5 ) continue;
00228            resultjet->push_back ((*jet));
00229 
00230            dphi = dphi*180./(4.*atan(1.));
00231            if( fabs(dphi-180) < 30. )
00232            {
00233    //           cout<<" Jet is found "<<endl;
00234               iejet = 1;
00235               phijet = (*jet).phi();
00236               etajet = (*jet).eta();
00237            } // dphi
00238          } //jet collection
00239          if(iejet == 0) resultjet->clear(); 
00240       }
00241     } // Jet collection
00242      
00243      if( resultjet->size() == 0 ) {
00244 
00245       iEvent.put( outputTColl, "GammaJetTracksCollection");
00246       iEvent.put( miniEcalRecHitCollection, "GammaJetEcalRecHitCollection");
00247       iEvent.put( miniHBHERecHitCollection, "GammaJetHBHERecHitCollection");
00248       iEvent.put( miniHORecHitCollection, "GammaJetHORecHitCollection");
00249       iEvent.put( miniHFRecHitCollection, "GammaJetHFRecHitCollection");
00250       iEvent.put( result, "GammaJetGammaBackToBackCollection");
00251       iEvent.put( resultjet, "GammaJetJetBackToBackCollection");
00252 
00253      return;
00254      }
00255     // cout<<" Accepted event "<<resultjet->size()<<" PHI gamma "<<phigamma<<std::endl;
00256 //
00257 // Add Ecal, Hcal RecHits around Egamma caluster
00258 //
00259 
00260 // Load EcalRecHits
00261 
00262     
00263     std::vector<edm::InputTag>::const_iterator i;
00264     for (i=ecalLabels_.begin(); i!=ecalLabels_.end(); i++) {
00265       edm::Handle<EcalRecHitCollection> ec;
00266       iEvent.getByLabel(*i,ec);
00267       if (!ec.isValid()) {
00268         // can't find it!
00269         if (!allowMissingInputs_) { 
00270           cout<<" No ECAL input "<<endl;
00271           *ec;  // will throw the proper exception
00272         }
00273       } else {
00274         for(EcalRecHitCollection::const_iterator recHit = (*ec).begin();
00275             recHit != (*ec).end(); ++recHit)
00276           {
00277 // EcalBarrel = 1, EcalEndcap = 2
00278             GlobalPoint pos = geo->getPosition(recHit->detid());
00279             double phihit = pos.phi();
00280             double etahit = pos.eta();
00281             
00282             double dphi = fabs(phigamma - phihit); 
00283             if(dphi > 4.*atan(1.)) dphi = 8.*atan(1.) - dphi;
00284             double deta = fabs(etagamma - etahit); 
00285             double dr = sqrt(dphi*dphi + deta*deta);
00286             if(dr<1.)  miniEcalRecHitCollection->push_back(*recHit);
00287             
00288             dphi = fabs(phijet - phihit); 
00289             if(dphi > 4.*atan(1.)) dphi = 8.*atan(1.) - dphi;
00290             deta = fabs(etajet - etahit); 
00291             dr = sqrt(dphi*dphi + deta*deta);
00292             if(dr<1.4)  miniEcalRecHitCollection->push_back(*recHit);
00293             
00294           }
00295       }
00296     }
00297 
00298 //   cout<<" Ecal is done "<<endl; 
00299 
00300       edm::Handle<HBHERecHitCollection> hbhe;
00301       iEvent.getByLabel(hbheLabel_,hbhe);
00302       if (!hbhe.isValid()) {
00303         // can't find it!
00304         if (!allowMissingInputs_) {
00305           *hbhe;  // will throw the proper exception
00306         }
00307       } else {
00308         const HBHERecHitCollection Hithbhe = *(hbhe.product());
00309         for(HBHERecHitCollection::const_iterator hbheItr=Hithbhe.begin(); hbheItr!=Hithbhe.end(); hbheItr++)
00310           {
00311             GlobalPoint pos = geo->getPosition(hbheItr->detid());
00312             double phihit = pos.phi();
00313             double etahit = pos.eta();
00314             
00315             double dphi = fabs(phigamma - phihit); 
00316             if(dphi > 4.*atan(1.)) dphi = 8.*atan(1.) - dphi;
00317             double deta = fabs(etagamma - etahit); 
00318             double dr = sqrt(dphi*dphi + deta*deta);
00319             
00320             
00321             if(dr<1.)  miniHBHERecHitCollection->push_back(*hbheItr);
00322             dphi = fabs(phijet - phihit); 
00323             if(dphi > 4.*atan(1.)) dphi = 8.*atan(1.) - dphi;
00324             deta = fabs(etajet - etahit); 
00325             dr = sqrt(dphi*dphi + deta*deta);
00326             if(dr<1.4)  miniHBHERecHitCollection->push_back(*hbheItr);
00327           }
00328       }
00329 
00330 //   cout<<" HBHE is done "<<endl; 
00331         
00332       edm::Handle<HORecHitCollection> ho;
00333       iEvent.getByLabel(hoLabel_,ho);
00334       if (!ho.isValid()) {
00335         // can't find it!
00336         if (!allowMissingInputs_) {
00337           *ho;  // will throw the proper exception
00338         }
00339       } else {
00340         const HORecHitCollection Hitho = *(ho.product());
00341         for(HORecHitCollection::const_iterator hoItr=Hitho.begin(); hoItr!=Hitho.end(); hoItr++)
00342           {
00343             GlobalPoint pos = geo->getPosition(hoItr->detid());
00344             double phihit = pos.phi();
00345             double etahit = pos.eta();
00346             
00347             double dphi = fabs(phigamma - phihit); 
00348             if(dphi > 4.*atan(1.)) dphi = 8.*atan(1.) - dphi;
00349             double deta = fabs(etagamma - etahit); 
00350             double dr = sqrt(dphi*dphi + deta*deta);
00351             
00352             if(dr<1.)  miniHORecHitCollection->push_back(*hoItr);
00353             dphi = fabs(phijet - phihit); 
00354             if(dphi > 4.*atan(1.)) dphi = 8.*atan(1.) - dphi;
00355             deta = fabs(etajet - etahit); 
00356             dr = sqrt(dphi*dphi + deta*deta);
00357             if(dr<1.4)  miniHORecHitCollection->push_back(*hoItr);
00358             
00359           }
00360       }
00361 
00362  //  cout<<" HO is done "<<endl; 
00363    
00364       edm::Handle<HFRecHitCollection> hf;
00365       iEvent.getByLabel(hfLabel_,hf);
00366       if (!hf.isValid()) {
00367         // can't find it!
00368         if (!allowMissingInputs_) {
00369           *hf;  // will throw the proper exception
00370         }
00371       } else {
00372         const HFRecHitCollection Hithf = *(hf.product());
00373         //  cout<<" Size of HF collection "<<Hithf.size()<<endl;
00374         for(HFRecHitCollection::const_iterator hfItr=Hithf.begin(); hfItr!=Hithf.end(); hfItr++)
00375           {
00376             GlobalPoint pos = geo->getPosition(hfItr->detid());
00377             
00378             double phihit = pos.phi();
00379             double etahit = pos.eta();
00380             
00381             double dphi = fabs(phigamma - phihit); 
00382             if(dphi > 4.*atan(1.)) dphi = 8.*atan(1.) - dphi;
00383             double deta = fabs(etagamma - etahit); 
00384             double dr = sqrt(dphi*dphi + deta*deta);
00385             
00386             if(dr<1.)  miniHFRecHitCollection->push_back(*hfItr);
00387             dphi = fabs(phijet - phihit); 
00388             if(dphi > 4.*atan(1.)) dphi = 8.*atan(1.) - dphi;
00389             deta = fabs(etajet - etahit); 
00390             dr = sqrt(dphi*dphi + deta*deta);
00391             
00392             if( dr < 1.4 )  miniHFRecHitCollection->push_back(*hfItr);
00393           }
00394       }
00395 
00396  //  cout<<" Size of mini HF collection "<<miniHFRecHitCollection->size()<<endl;
00397      
00398 
00399 // Track Collection   
00400       edm::Handle<reco::TrackCollection> trackCollection;
00401       iEvent.getByLabel(m_inputTrackLabel,trackCollection);
00402       if (!trackCollection.isValid()) {
00403         // can't find it!
00404         if (!allowMissingInputs_) {
00405           *trackCollection;  // will throw the proper exception
00406         }
00407       } else {
00408         const reco::TrackCollection tC = *(trackCollection.product());
00409   
00410         //Create empty output collections
00411 
00412 
00413         for (reco::TrackCollection::const_iterator track=tC.begin(); track!=tC.end(); track++)
00414           {
00415             
00416             double deta = track->momentum().eta() - etagamma; 
00417             
00418             double dphi = fabs(track->momentum().phi() - phigamma);
00419             
00420             if (dphi > atan(1.)*4.) dphi = 8.*atan(1.) - dphi;
00421             double ddir1 = sqrt(deta*deta+dphi*dphi);
00422             
00423             
00424             
00425             deta = track->momentum().eta() - etajet;
00426             dphi = fabs(track->momentum().phi() - phijet);
00427             if (dphi > atan(1.)*4.) dphi = 8.*atan(1.) - dphi;
00428             double ddir2 = sqrt(deta*deta+dphi*dphi);
00429             
00430             if( ddir1 < 1.4  || ddir2 < 1.4)      
00431               {
00432                 outputTColl->push_back(*track);
00433               } 
00434           }
00435       }
00436    
00437   //Put selected information in the event
00438   
00439   iEvent.put( outputTColl, "GammaJetTracksCollection");
00440   iEvent.put( miniEcalRecHitCollection, "GammaJetEcalRecHitCollection");
00441   iEvent.put( miniHBHERecHitCollection, "GammaJetHBHERecHitCollection");
00442   iEvent.put( miniHORecHitCollection, "GammaJetHORecHitCollection");
00443   iEvent.put( miniHFRecHitCollection, "GammaJetHFRecHitCollection");
00444   iEvent.put( result, "GammaJetGammaBackToBackCollection");
00445   iEvent.put( resultjet, "GammaJetJetBackToBackCollection");
00446 }
00447 //}