CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_4_5_patch3/src/RecoHI/HiEvtPlaneAlgos/src/EvtPlaneProducer.cc

Go to the documentation of this file.
00001 // -*- C++ -*-
00002 //
00003 // Package:    EvtPlaneProducer
00004 // Class:      EvtPlaneProducer
00005 // 
00013 //
00014 // Original Author:  Sergey Petrushanko
00015 //         Created:  Fri Jul 11 10:05:00 2008
00016 // $Id: EvtPlaneProducer.cc,v 1.18 2011/10/07 09:41:29 yilmaz Exp $
00017 //
00018 //
00019 #define TRACKCOLLECTION 1
00020 //#define RECOCHARGEDCANDIDATECOLLECTION 1
00021 
00022 
00023 // system include files
00024 #include <memory>
00025 #include <iostream>
00026 #include <time.h>
00027 #include "TMath.h"
00028 
00029 // user include files
00030 #include "FWCore/Framework/interface/Frameworkfwd.h"
00031 #include "FWCore/Framework/interface/EDProducer.h"
00032 #include "FWCore/Framework/interface/Event.h"
00033 #include "FWCore/Framework/interface/EventSetup.h"
00034 
00035 #include "FWCore/Framework/interface/MakerMacros.h"
00036 #include "FWCore/Framework/interface/ESHandle.h"
00037 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00038 
00039 #include "DataFormats/Candidate/interface/Candidate.h"
00040 #include "DataFormats/HeavyIonEvent/interface/EvtPlane.h"
00041 #include "DataFormats/RecoCandidate/interface/RecoChargedCandidate.h"
00042 #include "DataFormats/RecoCandidate/interface/RecoChargedCandidateFwd.h"
00043 
00044 #include "DataFormats/HcalRecHit/interface/HcalRecHitCollections.h"
00045 #include "DataFormats/Common/interface/EDProduct.h"
00046 #include "DataFormats/Common/interface/Ref.h"
00047 
00048 #include "DataFormats/Common/interface/Handle.h"
00049 #include "DataFormats/CaloTowers/interface/CaloTowerCollection.h"
00050 
00051 #include "FWCore/ServiceRegistry/interface/Service.h"
00052 #include "CommonTools/UtilAlgos/interface/TFileService.h"
00053 #include "DataFormats/TrackReco/interface/TrackFwd.h"
00054 #include "DataFormats/TrackReco/interface/Track.h"
00055 #include "DataFormats/VertexReco/interface/Vertex.h"
00056 #include "DataFormats/VertexReco/interface/VertexFwd.h"
00057 
00058 #include <cstdlib>
00059 #include "RecoHI/HiEvtPlaneAlgos/interface/HiEvtPlaneList.h"
00060 
00061 using namespace std;
00062 using namespace hi;
00063 
00064 //
00065 // class decleration
00066 //
00067 
00068 class EvtPlaneProducer : public edm::EDProducer {
00069 public:
00070   explicit EvtPlaneProducer(const edm::ParameterSet&);
00071   ~EvtPlaneProducer();
00072   
00073 private:
00074   //edm::Service<TFileService> fs;
00075   class GenPlane {
00076   public: 
00077     GenPlane(string name,double etaminval1,double etamaxval1,double etaminval2,double etamaxval2,int orderval){
00078       epname=name;
00079       etamin1=etaminval1;
00080       etamax1=etamaxval1;
00081       etamin2=etaminval2;
00082       etamax2=etamaxval2;
00083       sumsin=0;
00084       sumcos=0;
00085       order = (double) orderval;
00086     }
00087     ~GenPlane(){;}
00088     void addParticle(double w, double s, double c, double eta) {
00089       if(w < 0.001) return;
00090       if((eta>=etamin1 && eta<etamax1) || 
00091          (etamin2!= etamax2 && eta>=etamin2 && eta<etamax2 )) {
00092         sumsin+=w*s;
00093         sumcos+=w*c;
00094         //For tracking, w=1 and mult is the track multiplicity.  
00095         //For calorimetors, w is an energy that needs to be subsequently 
00096         //converted if an  multiplicity if needed
00097         mult+=w;
00098       }
00099     }
00100     
00101     double getAngle(double &ang, double &sv, double &cv){
00102       ang = -10;
00103       sv = 0;
00104       cv = 0;
00105       sv = sumsin;
00106       cv = sumcos;
00107       double q = sv*sv+cv*cv;
00108       if(q>0) ang = atan2(sv,cv)/order;
00109       return ang;
00110     }
00111     void reset() {
00112       sumsin=0;
00113       sumcos=0;
00114       mult = 0;
00115     }
00116   private:
00117     string epname;
00118     double etamin1;
00119     double etamax1;
00120 
00121     double etamin2;
00122     double etamax2;
00123     double sumsin;
00124     double sumcos;
00125     int mult;
00126     double order;
00127   };
00128   
00129 
00130   GenPlane *rp[NumEPNames];
00131 
00132   virtual void beginJob() ;
00133   virtual void produce(edm::Event&, const edm::EventSetup&);
00134   virtual void endJob() ;
00135   
00136   // ----------member data ---------------------------
00137   edm::InputTag vtxCollection_;
00138   edm::InputTag caloCollection_;
00139   edm::InputTag trackCollection_;
00140   bool useECAL_;
00141   bool useHCAL_;
00142   bool useTrack_;
00143   bool useTrackPtWeight_;
00144   double minet_;
00145   double maxet_;
00146   double minpt_;
00147   double maxpt_;
00148   double minvtx_;
00149   double maxvtx_;
00150   double dzerr_;
00151   double chi2_;
00152 
00153   bool storeNames_;
00154 };
00155 
00156 //
00157 // constants, enums and typedefs
00158 //
00159 
00160 
00161 //
00162 // static data member definitions
00163 //
00164 
00165 //
00166 // constructors and destructor
00167 //
00168 EvtPlaneProducer::EvtPlaneProducer(const edm::ParameterSet& iConfig)
00169 {
00170   
00171   
00172   //register your products
00173   vtxCollection_  = iConfig.getParameter<edm::InputTag>("vtxCollection_");
00174   caloCollection_  = iConfig.getParameter<edm::InputTag>("caloCollection_");
00175   trackCollection_  = iConfig.getParameter<edm::InputTag>("trackCollection_");
00176   useECAL_ = iConfig.getUntrackedParameter<bool>("useECAL_",true);
00177   useHCAL_ = iConfig.getUntrackedParameter<bool>("useHCAL_",true);
00178   useTrack_ = iConfig.getUntrackedParameter<bool>("useTrack",true);
00179   useTrackPtWeight_ = iConfig.getUntrackedParameter<bool>("useTrackPtWeight_",true);
00180   minet_ = iConfig.getUntrackedParameter<double>("minet_",0.2);
00181   maxet_ = iConfig.getUntrackedParameter<double>("maxet_",500.);
00182   minpt_ = iConfig.getUntrackedParameter<double>("minpt_",0.3);
00183   maxpt_ = iConfig.getUntrackedParameter<double>("maxpt_",2.5);
00184   minvtx_ = iConfig.getUntrackedParameter<double>("minvtx_",-50.);
00185   maxvtx_ = iConfig.getUntrackedParameter<double>("maxvtx_",50.);
00186   dzerr_ = iConfig.getUntrackedParameter<double>("dzerr_",10.);
00187   chi2_  = iConfig.getUntrackedParameter<double>("chi2_",40.);
00188 
00189   storeNames_ = 1;
00190   produces<reco::EvtPlaneCollection>("recoLevel");
00191   for(int i = 0; i<NumEPNames; i++ ) {
00192     rp[i] = new GenPlane(EPNames[i].data(),EPEtaMin1[i],EPEtaMax1[i],EPEtaMin2[i],EPEtaMax2[i],EPOrder[i]);
00193   }
00194 }
00195 
00196 
00197 EvtPlaneProducer::~EvtPlaneProducer()
00198 {
00199   
00200   // do anything here that needs to be done at desctruction time
00201   // (e.g. close files, deallocate resources etc.)
00202   
00203 }
00204 
00205 
00206 //
00207 // member functions
00208 //
00209 
00210 // ------------ method called to produce the data  ------------
00211 void
00212 EvtPlaneProducer::produce(edm::Event& iEvent, const edm::EventSetup& iSetup)
00213 {
00214   using namespace edm;
00215   using namespace reco;
00216 
00217   int vs_sell;
00218   float vzr_sell;
00219   float vzErr_sell;
00220   //
00221   //Get Vertex
00222   //
00223   edm::Handle<reco::VertexCollection> vertexCollection3;
00224   iEvent.getByLabel(vtxCollection_,vertexCollection3);
00225   const reco::VertexCollection * vertices3 = vertexCollection3.product();
00226   vs_sell = vertices3->size();
00227   if(vs_sell>0) {
00228     vzr_sell = vertices3->begin()->z();
00229     vzErr_sell = vertices3->begin()->zError();
00230   } else
00231     vzr_sell = -999.9;
00232   //
00233   for(int i = 0; i<NumEPNames; i++) rp[i]->reset();
00234   if(vzr_sell>minvtx_ && vzr_sell<maxvtx_) {
00235     //calorimetry part
00236     
00237     double tower_eta, tower_phi;
00238     double tower_energy, tower_energy_e, tower_energy_h;
00239     double tower_energyet, tower_energyet_e, tower_energyet_h;
00240     double s1, s2, s11, s21,s13,s23,s14,s24,s15,s25,s16,s26;
00241     Handle<CaloTowerCollection> calotower;
00242     iEvent.getByLabel(caloCollection_,calotower);
00243     
00244     if(calotower.isValid()){
00245       
00246       for (CaloTowerCollection::const_iterator j = calotower->begin();j !=calotower->end(); j++) {   
00247         tower_eta        = j->eta();
00248         tower_phi        = j->phi();
00249         tower_energy_e   = j->emEnergy();
00250         tower_energy_h   = j->hadEnergy();
00251         tower_energy     = tower_energy_e + tower_energy_h;
00252         tower_energyet_e   = j->emEt();
00253         tower_energyet_h   = j->hadEt();
00254         tower_energyet     = tower_energyet_e + tower_energyet_h;
00255         
00256         s1 = sin(2.*tower_phi);
00257         s2 = cos(2.*tower_phi);    
00258         s11 = sin(tower_phi);
00259         s21 = cos(tower_phi);
00260         s13 = sin(3.*tower_phi);
00261         s23 = cos(3.*tower_phi);
00262         s14 = sin(4.*tower_phi);
00263         s24 = cos(4.*tower_phi);
00264         s15 = sin(5.*tower_phi);
00265         s25 = cos(5.*tower_phi);
00266         s16 = sin(6.*tower_phi);
00267         s26 = cos(6.*tower_phi);
00268         
00269         if(tower_energyet<minet_) continue;
00270         if(tower_energyet>maxet_) continue;;
00271 
00272         rp[etEcal    ]->addParticle (tower_energyet_e, s1,    s2,    tower_eta);
00273         rp[etEcalP   ]->addParticle (tower_energyet_e, s1,    s2,    tower_eta);
00274         rp[etEcalM   ]->addParticle (tower_energyet_e, s1,    s2,    tower_eta);
00275         rp[etHcal    ]->addParticle (tower_energyet_h, s1,    s2,    tower_eta);
00276         rp[etHcalP   ]->addParticle (tower_energyet_h, s1,    s2,    tower_eta);
00277         rp[etHcalM   ]->addParticle (tower_energyet_h, s1,    s2,    tower_eta);
00278         rp[etHF      ]->addParticle (tower_energyet,   s1,    s2,    tower_eta);
00279         rp[etHFp     ]->addParticle (tower_energyet,   s1,    s2,    tower_eta);
00280         rp[etHFm     ]->addParticle (tower_energyet,   s1,    s2,    tower_eta);
00281         rp[etCaloHFP ]->addParticle (tower_energyet,   s1,    s2,    tower_eta);
00282         rp[etCaloHFM ]->addParticle (tower_energyet,   s1,    s2,    tower_eta);
00283 
00284         rp[etHF3     ]->addParticle (tower_energyet,   s13,   s23,   tower_eta);
00285         rp[etHFp3    ]->addParticle (tower_energyet,   s13,   s23,   tower_eta);    
00286         rp[etHFm3    ]->addParticle (tower_energyet,   s13,   s23,   tower_eta);
00287 
00288         rp[etHF4     ]->addParticle (tower_energyet,   s14,   s24,   tower_eta);
00289         rp[etHFp4    ]->addParticle (tower_energyet,   s14,   s24,   tower_eta);    
00290         rp[etHFm4    ]->addParticle (tower_energyet,   s14,   s24,   tower_eta);
00291 
00292         rp[etHF5     ]->addParticle (tower_energyet,   s15,   s25,   tower_eta);
00293         rp[etHFp5    ]->addParticle (tower_energyet,   s15,   s25,   tower_eta);    
00294         rp[etHFm5    ]->addParticle (tower_energyet,   s15,   s25,   tower_eta);
00295 
00296         rp[etHF6     ]->addParticle (tower_energyet,   s16,   s26,   tower_eta);
00297         rp[etHFp6    ]->addParticle (tower_energyet,   s16,   s26,   tower_eta);    
00298         rp[etHFm6    ]->addParticle (tower_energyet,   s16,   s26,   tower_eta);
00299       } 
00300     }
00301     //Tracking part
00302     
00303     double track_eta;
00304     double track_phi;
00305     double track_pt;
00306     double track_charge;
00307 #ifdef TRACKCOLLECTION  
00308     Handle<reco::TrackCollection> tracks;
00309     iEvent.getByLabel(trackCollection_, tracks);
00310     
00311     if(tracks.isValid()){
00312       for(reco::TrackCollection::const_iterator j = tracks->begin(); j != tracks->end(); j++){
00313         
00314         //Find possible collections with command line: edmDumpEventContent *.root
00315 #endif
00316 #ifdef RECOCHARGEDCANDIDATECOLLECTION
00317         edm::Handle<reco::RecoChargedCandidateCollection> trackCollection;
00318         iEvent.getByLabel("allMergedPtSplit12Tracks",trackCollection);
00319         //      iEvent.getByLabel("hiGoodMergedTracks",trackCollection);
00320         if(trackCollection.isValid()){
00321           
00322           const reco::RecoChargedCandidateCollection * tracks = trackCollection.product();
00323           for(reco::RecoChargedCandidateCollection::const_iterator j = tracks->begin(); j != tracks->end(); j++){
00324 #endif  
00325         edm::Handle<reco::VertexCollection> vertex;
00326         iEvent.getByLabel(vtxCollection_, vertex);
00327         
00328 // find the vertex point and error
00329 
00330         math::XYZPoint vtxPoint(0.0,0.0,0.0);
00331         double vzErr =0.0, vxErr=0.0, vyErr=0.0;
00332         if(vertex->size()>0) {
00333           vtxPoint=vertex->begin()->position();
00334           vzErr=vertex->begin()->zError();
00335           vxErr=vertex->begin()->xError();
00336           vyErr=vertex->begin()->yError();
00337         }
00338         bool accepted = true;
00339         bool isPixel = false;
00340         // determine if the track is a pixel track
00341         if ( j->numberOfValidHits() < 7 ) isPixel = true;
00342         
00343         // determine the vertex significance 
00344         double d0=0.0, dz=0.0, d0sigma=0.0, dzsigma=0.0;
00345         d0 = -1.*j->dxy(vtxPoint);
00346         dz = j->dz(vtxPoint);
00347         d0sigma = sqrt(j->d0Error()*j->d0Error()+vxErr*vyErr);
00348         dzsigma = sqrt(j->dzError()*j->dzError()+vzErr*vzErr);
00349         
00350         // cuts for pixel tracks
00351         if( isPixel )
00352           {
00353             // dz significance cut 
00354             if ( fabs(dz/dzsigma) > dzerr_ ) accepted = false;
00355             
00356             // chi2/ndof cut 
00357             if ( j->normalizedChi2() > chi2_ ) accepted = false;
00358           }
00359         
00360         // cuts for full tracks
00361         if ( ! isPixel)
00362           {
00363             // dz and d0 significance cuts 
00364             if ( fabs(dz/dzsigma) > 3 ) accepted = false;
00365             if ( fabs(d0/d0sigma) > 3 ) accepted = false;
00366             
00367             // pt resolution cut
00368             if ( j->ptError()/j->pt() > 0.05 ) accepted = false;
00369             
00370             // number of valid hits cut
00371             if ( j->numberOfValidHits() < 12 ) accepted = false;
00372           }
00373         if( accepted ) {
00374           track_eta = j->eta();
00375           track_phi = j->phi();
00376           track_pt = j->pt();
00377           track_charge = j->charge();
00378           double s =sin(2*track_phi);
00379           double c =cos(2*track_phi);
00380           double s3 =sin(3*track_phi);
00381           double c3 =cos(3*track_phi);
00382           double s4 =sin(4*track_phi);
00383           double c4 =cos(4*track_phi);
00384           double s5 =sin(5*track_phi);
00385           double c5 =cos(5*track_phi);
00386           double s6 =sin(6*track_phi);
00387           double c6 =cos(6*track_phi);
00388           double w = 1;
00389           if(useTrackPtWeight_) {
00390             w = track_pt;
00391             if(w>2.5) w=2.0;   //v2 starts decreasing above ~2.5 GeV/c
00392           }
00393           if(track_pt<minpt_) continue;
00394           if(track_pt>maxpt_) continue;
00395           rp[EvtPlaneFromTracksMidEta]->addParticle(w,s,c,track_eta);
00396           rp[EvtPTracksPosEtaGap]->addParticle(w,s,c,track_eta);
00397           rp[EvtPTracksNegEtaGap]->addParticle(w,s,c,track_eta);
00398           rp[EPTracksMid3]->addParticle(w,s3,c3,track_eta);
00399           rp[EPTracksPos3]->addParticle(w,s3,c3,track_eta);
00400           rp[EPTracksNeg3]->addParticle(w,s3,c3,track_eta);
00401           rp[EPTracksMid4]->addParticle(w,s4,c4,track_eta);
00402           rp[EPTracksPos4]->addParticle(w,s4,c4,track_eta);
00403           rp[EPTracksNeg4]->addParticle(w,s4,c4,track_eta);
00404           rp[EPTracksMid5]->addParticle(w,s5,c5,track_eta);
00405           rp[EPTracksPos5]->addParticle(w,s5,c5,track_eta);
00406           rp[EPTracksNeg5]->addParticle(w,s5,c5,track_eta);
00407           rp[EPTracksMid6]->addParticle(w,s6,c6,track_eta);
00408           rp[EPTracksPos6]->addParticle(w,s6,c6,track_eta);
00409           rp[EPTracksNeg6]->addParticle(w,s6,c6,track_eta);
00410           
00411         }
00412         
00413       }
00414     }
00415       }
00416       std::auto_ptr<EvtPlaneCollection> evtplaneOutput(new EvtPlaneCollection);
00417       
00418   EvtPlane *ep[NumEPNames];
00419   
00420   double ang=-10;
00421   double sv = 0;
00422   double cv = 0;
00423 
00424   for(int i = 0; i<NumEPNames; i++) {
00425     rp[i]->getAngle(ang,sv,cv);
00426     if(storeNames_) ep[i] = new EvtPlane(ang,sv,cv,EPNames[i]);
00427     else ep[i] = new EvtPlane(ang,sv,cv,"");
00428   }
00429   if(useTrack_) {
00430     for(int i = 0; i<9; i++) {
00431       evtplaneOutput->push_back(*ep[i]);
00432     }  
00433   }
00434   for(int i = 9; i<NumEPNames; i++) {
00435     if(useECAL_ && !useHCAL_) {
00436       if(EPNames[i].rfind("Ecal")!=string::npos) {
00437         evtplaneOutput->push_back(*ep[i]);
00438       }
00439     } else if (useHCAL_ && !useECAL_) {
00440       if(EPNames[i].rfind("Hcal")!=string::npos) {
00441         evtplaneOutput->push_back(*ep[i]);
00442       }
00443     }else if (useECAL_ && useHCAL_) {
00444       evtplaneOutput->push_back(*ep[i]);
00445     }
00446   }
00447   
00448   iEvent.put(evtplaneOutput, "recoLevel");
00449   //  storeNames_ = 0;
00450 }
00451   
00452   // ------------ method called once each job just before starting event loop  ------------
00453 void 
00454 EvtPlaneProducer::beginJob()
00455 {
00456 }
00457 
00458 // ------------ method called once each job just after ending the event loop  ------------
00459 void 
00460 EvtPlaneProducer::endJob() {
00461 }
00462 
00463 //define this as a plug-in
00464 DEFINE_FWK_MODULE(EvtPlaneProducer);