CMS 3D CMS Logo

PhotonsWithConversionsAnalyzer Class Reference

#include <RecoEgamma/Examples/plugins/PhotonsWithConversionsAnalyzer.h>

Inheritance diagram for PhotonsWithConversionsAnalyzer:

edm::EDAnalyzer

List of all members.

Public Member Functions

virtual void analyze (const edm::Event &, const edm::EventSetup &)
virtual void beginJob (const edm::EventSetup &)
virtual void endJob ()
 PhotonsWithConversionsAnalyzer (const edm::ParameterSet &)
virtual ~PhotonsWithConversionsAnalyzer ()

Private Member Functions

float etaTransformation (float a, float b)

Private Attributes

TFile * fOutputFile_
std::string fOutputFileName_
TH2F * h2_tk_inPtVsR_
TH2F * h2_tk_nHitsVsR_
TH1F * h_deltaEta_
TH1F * h_deltaPhi_
TH1F * h_ErecoEMC_
TH1F * h_MCConvE_
TH1F * h_MCConvEta_
TH1F * h_MCConvPt_
TH1F * h_MCphoE_
TH1F * h_MCphoEta_
TH1F * h_MCphoPhi_
TH1F * h_phoE_
TH1F * h_phoEta_
TH1F * h_phoPhi_
TH1F * h_scE_
TH1F * h_scEt_
TH1F * h_scEta_
TH1F * h_scPhi_
std::string HepMCLabel
int nEvt_
int nMatched_
int nMCPho_
std::string photonCollection_
std::string photonCollectionProducer_
std::string SimHitLabel
std::string SimTkLabel
std::string SimVtxLabel
PhotonMCTruthFinderthePhotonMCTruthFinder_


Detailed Description

Definition at line 23 of file PhotonsWithConversionsAnalyzer.h.


Constructor & Destructor Documentation

PhotonsWithConversionsAnalyzer::PhotonsWithConversionsAnalyzer ( const edm::ParameterSet pset  )  [explicit]

Definition at line 51 of file PhotonsWithConversionsAnalyzer.cc.

References edm::ParameterSet::getParameter(), photonCollection_, and photonCollectionProducer_.

00052   {
00053 
00054   photonCollectionProducer_ = pset.getParameter<std::string>("phoProducer");
00055   photonCollection_ = pset.getParameter<std::string>("photonCollection");
00056   //
00057 
00058 
00059 }

PhotonsWithConversionsAnalyzer::~PhotonsWithConversionsAnalyzer (  )  [virtual]

Definition at line 63 of file PhotonsWithConversionsAnalyzer.cc.

References thePhotonMCTruthFinder_.

00063                                                                 {
00064 
00065 
00066   delete thePhotonMCTruthFinder_;
00067 
00068 }


Member Function Documentation

void PhotonsWithConversionsAnalyzer::analyze ( const edm::Event e,
const edm::EventSetup  
) [virtual]

Loop over recontructed photons

End loop over Reco particles

End loop over MC particles

Implements edm::EDAnalyzer.

Definition at line 162 of file PhotonsWithConversionsAnalyzer.cc.

References conversions_cfi::conversions, GenMuonPlsPt100GeV_cfg::cout, deltaPhi(), lat::endl(), etaTransformation(), PhotonMCTruthFinder::find(), edm::Event::getByLabel(), h2_tk_inPtVsR_, h2_tk_nHitsVsR_, h_deltaEta_, h_deltaPhi_, h_ErecoEMC_, h_MCConvEta_, h_MCphoE_, h_MCphoEta_, h_MCphoPhi_, h_phoE_, h_phoEta_, h_phoPhi_, h_scE_, h_scEt_, h_scEta_, h_scPhi_, i, edm::Event::id(), nEvt_, nMatched_, nMCPho_, photonCollection_, photonCollectionProducer_, Geom::pi(), funct::pow(), funct::sqrt(), thePhotonMCTruthFinder_, tracks, and Geom::twoPi().

00163 {
00164   
00165   
00166   using namespace edm;
00167   const float etaPhiDistance=0.01;
00168   // Fiducial region
00169   const float TRK_BARL =0.9;
00170   const float BARL = 1.4442; // DAQ TDR p.290
00171   const float END_LO = 1.566;
00172   const float END_HI = 2.5;
00173  // Electron mass
00174   const Float_t mElec= 0.000511;
00175 
00176 
00177   nEvt_++;  
00178   LogInfo("ConvertedPhotonAnalyzer") << "ConvertedPhotonAnalyzer Analyzing event number: " << e.id() << " Global Counter " << nEvt_ <<"\n";
00179   //  LogDebug("ConvertedPhotonAnalyzer") << "ConvertedPhotonAnalyzer Analyzing event number: "  << e.id() << " Global Counter " << nEvt_ <<"\n";
00180   std::cout << "ConvertedPhotonAnalyzer Analyzing event number: "  << e.id() << " Global Counter " << nEvt_ <<"\n";
00181  
00182   
00184   Handle<reco::PhotonCollection> photonHandle; 
00185   e.getByLabel(photonCollectionProducer_, photonCollection_ , photonHandle);
00186   const reco::PhotonCollection photonCollection = *(photonHandle.product());
00187   std::cout  << "ConvertedPhotonAnalyzer  Photons with conversions collection size " << photonCollection.size() << "\n";
00188 
00189 
00191   std::cout  << " ConvertedPhotonAnalyzer Looking for MC truth " << "\n";
00192   
00193   //get simtrack info
00194   std::vector<SimTrack> theSimTracks;
00195   std::vector<SimVertex> theSimVertices;
00196   
00197   edm::Handle<SimTrackContainer> SimTk;
00198   edm::Handle<SimVertexContainer> SimVtx;
00199   e.getByLabel("g4SimHits",SimTk);
00200   e.getByLabel("g4SimHits",SimVtx);
00201   
00202   theSimTracks.insert(theSimTracks.end(),SimTk->begin(),SimTk->end());
00203   theSimVertices.insert(theSimVertices.end(),SimVtx->begin(),SimVtx->end());
00204   
00205   std::vector<PhotonMCTruth> mcPhotons=thePhotonMCTruthFinder_->find (theSimTracks,  theSimVertices);  
00206   std::cout << " ConvertedPhotonAnalyzer mcPhotons size " <<  mcPhotons.size() << std::endl;
00207 
00208 
00209 
00210   // Loop over simulated photons
00211   int iDet=0;
00212   int iRadius=-1;
00213   int indPho=0;
00214 
00215   for ( std::vector<PhotonMCTruth>::const_iterator mcPho=mcPhotons.begin(); mcPho !=mcPhotons.end(); mcPho++) {
00216     float mcPhi= (*mcPho).fourMomentum().phi();
00217     float mcEta= (*mcPho).fourMomentum().pseudoRapidity();
00218     mcEta = etaTransformation(mcEta, (*mcPho).primaryVertex().z() ); 
00219         
00220     if ( (*mcPho).fourMomentum().et() < 20 ) continue;
00221     //    if ( ! (  fabs(mcEta) <= BARL || ( fabs(mcEta) >= END_LO && fabs(mcEta) <=END_HI ) ) ) {
00222     //     continue;
00223     //} // all ecal fiducial region
00224 
00225     
00226 
00227     h_MCphoE_->Fill(  (*mcPho).fourMomentum().e());
00228     h_MCphoEta_->Fill( (*mcPho).fourMomentum().eta());
00229     h_MCphoPhi_->Fill(  (*mcPho).fourMomentum().phi());
00230     
00231     
00232     
00233     // keep only visible conversions
00234     if (  (*mcPho).isAConversion() == 0 ) continue;
00235 
00236     
00237     nMCPho_++;
00238 
00239     h_MCConvEta_ ->Fill ( fabs( (*mcPho).fourMomentum().pseudoRapidity()) - 0.001);
00240     
00241     
00242     bool REJECTED;
00243 
00244 
00246     //std::cout   << " ConvertedPhotonAnalyzer  Starting loop over photon candidates " << "\n";
00247     for( reco::PhotonCollection::const_iterator  iPho = photonCollection.begin(); iPho != photonCollection.end(); iPho++) {
00248       REJECTED=false;
00249       
00250       //      std::cout  << " ConvertedPhotonAnalyzer Reco SC energy " << (*iPho).superCluster()->energy() <<  "\n";
00251 
00252       float phiClu=(*iPho).superCluster()->phi();
00253       float etaClu=(*iPho).superCluster()->eta();
00254       float deltaPhi = phiClu-mcPhi;
00255       float deltaEta = etaClu-mcEta;
00256       
00257       
00258       if ( deltaPhi > Geom::pi()  ) deltaPhi -= Geom::twoPi();
00259       if ( deltaPhi < -Geom::pi() ) deltaPhi += Geom::twoPi();
00260       deltaPhi=pow(deltaPhi,2);
00261       deltaEta=pow(deltaEta,2);
00262       float delta =  deltaPhi+deltaEta ; 
00263       if (  delta >= etaPhiDistance  )  REJECTED=true;
00264       
00265 
00266       //      if ( ! (  fabs(etaClu) <= BARL || ( fabs(etaClu) >= END_LO && fabs(etaClu) <=END_HI ) ) ) REJECTED=true;
00267       
00268       if ( REJECTED ) continue;
00269       std::cout << " MATCHED " << std::endl;
00270       nMatched_++;
00271 
00272 
00273       //      std::cout << " ConvertedPhotonAnalyzer Matching candidate " << std::endl;
00274       
00275       // std::cout << " ConvertedPhotonAnalyzer Photons isAconversion " << (*mcPho).isAConversion() << " mcMatchingPhoton energy " <<  (*mcPho).fourMomentum().e()  << " ConvertedPhotonAnalyzer conversion vertex R " <<  (*mcPho).vertex().perp() << " Z " <<  (*mcPho).vertex().z() <<  std::endl;
00276       
00277   
00278       h_ErecoEMC_->Fill(   (*iPho).superCluster()->energy()/(*mcPho).fourMomentum().e());
00279       h_deltaPhi_-> Fill ( (*iPho).superCluster()->position().phi()- mcPhi);
00280       h_deltaEta_-> Fill ( (*iPho).superCluster()->position().eta()- mcEta);
00281       
00282       h_scE_->Fill( (*iPho).superCluster()->energy() );
00283       h_scEt_->Fill( (*iPho).superCluster()->energy()/cosh( (*iPho).superCluster()->position().eta()) );
00284       h_scEta_->Fill( (*iPho).superCluster()->position().eta() );
00285       h_scPhi_->Fill( (*iPho).superCluster()->position().phi() );
00286 
00287 
00288       h_phoE_->Fill( (*iPho).energy() );
00289       h_phoEta_->Fill( (*iPho).eta() );
00290       h_phoPhi_->Fill( (*iPho).phi() );
00291       
00292       if ( !(*iPho).isConverted() ) continue;
00293       //   std::cout << " This photons has " << (*iPho).conversions().size() << " conversions candidates " << std::endl;      
00294       std::vector<reco::ConversionRef> conversions = (*iPho).conversions();
00295       
00296       for (unsigned int i=0; i<conversions.size(); i++) {
00297         //std::cout << " Conversion candidate Energy " << (*iPho).energy() << " number of tracks " << conversions[i]->nTracks() << std::endl;
00298         std::vector<reco::TrackRef> tracks = conversions[i]->tracks();
00299 
00300         for (unsigned int i=0; i<tracks.size(); i++) {
00301 
00302           //      std::cout  << " ConvertedPhotonAnalyzer Reco Track charge " <<  tracks[i]->charge() << "  Num of RecHits " << tracks[i]->recHitsSize() << " inner momentum " <<  sqrt ( tracks[i]->innerMomentum().Mag2() )  <<  "\n"; 
00303         
00304         
00305           h2_tk_nHitsVsR_ -> Fill (  (*mcPho).vertex().perp(), tracks[i]->recHitsSize()   );
00306           h2_tk_inPtVsR_->Fill  (  (*mcPho).vertex().perp(),  sqrt( tracks[i]->innerMomentum().Mag2() ) );
00307         }
00308         
00309       }
00310       
00311       
00312     }
00313     
00314   }
00315   
00316 
00317 }

void PhotonsWithConversionsAnalyzer::beginJob ( const edm::EventSetup setup  )  [virtual]

Reco - MC

Reimplemented from edm::EDAnalyzer.

Definition at line 71 of file PhotonsWithConversionsAnalyzer.cc.

References h2_tk_inPtVsR_, h2_tk_nHitsVsR_, h_deltaEta_, h_deltaPhi_, h_ErecoEMC_, h_MCConvE_, h_MCConvEta_, h_MCConvPt_, h_MCphoE_, h_MCphoEta_, h_MCphoPhi_, h_phoE_, h_phoEta_, h_phoPhi_, h_scE_, h_scEt_, h_scEta_, h_scPhi_, nEvt_, nMatched_, nMCPho_, and thePhotonMCTruthFinder_.

00072 {
00073 
00074 
00075   nEvt_=0;
00076   nMCPho_=0;
00077   nMatched_=0;
00078 
00079 
00080   thePhotonMCTruthFinder_ = new PhotonMCTruthFinder();
00081 
00082 
00083 
00084   edm::Service<TFileService> fs;
00085 
00086 
00088   h_ErecoEMC_  = fs->make<TH1F>("deltaE","    delta(reco-mc) energy",100,0.,2.);
00089   h_deltaPhi_ = fs->make<TH1F>("deltaPhi","  delta(reco-mc) phi",100,-0.1, 0.1);
00090   h_deltaEta_ = fs->make<TH1F>("deltaEta","  delta(reco-mc) eta",100,-0.05, 0.05);
00091   
00093   h_MCphoE_ = fs->make<TH1F>("MCphoE","MC photon energy",100,0.,100.);
00094   h_MCphoPhi_ = fs->make<TH1F>("MCphoPhi","MC photon phi",40,-3.14, 3.14);
00095   h_MCphoEta_ = fs->make<TH1F>("MCphoEta","MC photon eta",40,-3., 3.);
00096 
00097 
00099   h_MCConvE_ = fs->make<TH1F>("MCConvE","MC converted photon energy",100,0.,100.);
00100   h_MCConvPt_ = fs->make<TH1F>("MCConvPt","MC converted photon pt",100,0.,100.);
00101   h_MCConvEta_ = fs->make<TH1F>("MCConvEta","MC converted photon eta",50, 0., 2.5);
00102 
00104   h_scE_ = fs->make<TH1F>("scE","SC Energy ",100,0., 200.);
00105   h_scEt_ = fs->make<TH1F>("scEt","SC Et ",100,0., 200.);
00106   h_scEta_ = fs->make<TH1F>("scEta","SC Eta ",40,-3., 3.);
00107   h_scPhi_ = fs->make<TH1F>("scPhi","SC Phi ",40, -3.14, 3.14);
00108   //
00109   h_phoE_ = fs->make<TH1F>("phoE","Photon Energy ",100,0., 200.);
00110   h_phoEta_ = fs->make<TH1F>("phoEta","Photon Eta ",40,-3., 3.);
00111   h_phoPhi_ = fs->make<TH1F>("phoPhi","Photon  Phi ",40,  -3.14, 3.14);
00112 
00113   // Recontructed tracks from converted photon candidates
00114   h2_tk_nHitsVsR_ = fs->make<TH2F>("tknHitsVsR","Tracks Hits vs R  ", 12,0.,120.,20,0.5, 20.5);
00115   h2_tk_inPtVsR_  = fs->make<TH2F>("tkInPtvsR","Tracks inner Pt vs R  ", 12,0.,120.,100,0., 100.);
00116 
00117   
00118   return ;
00119 }

void PhotonsWithConversionsAnalyzer::endJob ( void   )  [virtual]

Reimplemented from edm::EDAnalyzer.

Definition at line 322 of file PhotonsWithConversionsAnalyzer.cc.

References GenMuonPlsPt100GeV_cfg::cout, and nEvt_.

00323 {
00324 
00325 
00326 
00327        
00328   //   fOutputFile_->Write() ;
00329   // fOutputFile_->Close() ;
00330   
00331    edm::LogInfo("ConvertedPhotonAnalyzer") << "Analyzed " << nEvt_  << "\n";
00332    // std::cout  << "::endJob Analyzed " << nEvt_ << " events " << " with total " << nPho_ << " Photons " << "\n";
00333    std::cout  << "ConvertedPhotonAnalyzer::endJob Analyzed " << nEvt_ << " events " << "\n";
00334    
00335    return ;
00336 }

float PhotonsWithConversionsAnalyzer::etaTransformation ( float  a,
float  b 
) [private]

Definition at line 122 of file PhotonsWithConversionsAnalyzer.cc.

References ETA, etaBarrelEndcap, funct::log(), PI, R_ECAL, funct::tan(), TWOPI, and Z_Endcap.

Referenced by analyze().

00122                                                                                             {
00123 
00124 //---Definitions
00125         const float PI    = 3.1415927;
00126         const float TWOPI = 2.0*PI;
00127 
00128 //---Definitions for ECAL
00129         const float R_ECAL           = 136.5;
00130         const float Z_Endcap         = 328.0;
00131         const float etaBarrelEndcap  = 1.479; 
00132    
00133 //---ETA correction
00134 
00135         float Theta = 0.0  ; 
00136         float ZEcal = R_ECAL*sinh(EtaParticle)+Zvertex;
00137 
00138         if(ZEcal != 0.0) Theta = atan(R_ECAL/ZEcal);
00139         if(Theta<0.0) Theta = Theta+PI ;
00140         float ETA = - log(tan(0.5*Theta));
00141          
00142         if( fabs(ETA) > etaBarrelEndcap )
00143           {
00144            float Zend = Z_Endcap ;
00145            if(EtaParticle<0.0 )  Zend = -Zend ;
00146            float Zlen = Zend - Zvertex ;
00147            float RR = Zlen/sinh(EtaParticle); 
00148            Theta = atan(RR/Zend);
00149            if(Theta<0.0) Theta = Theta+PI ;
00150            ETA = - log(tan(0.5*Theta));               
00151           } 
00152 //---Return the result
00153         return ETA;
00154 //---end
00155 }


Member Data Documentation

TFile* PhotonsWithConversionsAnalyzer::fOutputFile_ [private]

Definition at line 46 of file PhotonsWithConversionsAnalyzer.h.

std::string PhotonsWithConversionsAnalyzer::fOutputFileName_ [private]

Definition at line 45 of file PhotonsWithConversionsAnalyzer.h.

TH2F* PhotonsWithConversionsAnalyzer::h2_tk_inPtVsR_ [private]

Definition at line 93 of file PhotonsWithConversionsAnalyzer.h.

Referenced by analyze(), and beginJob().

TH2F* PhotonsWithConversionsAnalyzer::h2_tk_nHitsVsR_ [private]

Definition at line 91 of file PhotonsWithConversionsAnalyzer.h.

Referenced by analyze(), and beginJob().

TH1F* PhotonsWithConversionsAnalyzer::h_deltaEta_ [private]

Definition at line 64 of file PhotonsWithConversionsAnalyzer.h.

Referenced by analyze(), and beginJob().

TH1F* PhotonsWithConversionsAnalyzer::h_deltaPhi_ [private]

Definition at line 63 of file PhotonsWithConversionsAnalyzer.h.

Referenced by analyze(), and beginJob().

TH1F* PhotonsWithConversionsAnalyzer::h_ErecoEMC_ [private]

Definition at line 62 of file PhotonsWithConversionsAnalyzer.h.

Referenced by analyze(), and beginJob().

TH1F* PhotonsWithConversionsAnalyzer::h_MCConvE_ [private]

Definition at line 75 of file PhotonsWithConversionsAnalyzer.h.

Referenced by beginJob().

TH1F* PhotonsWithConversionsAnalyzer::h_MCConvEta_ [private]

Definition at line 77 of file PhotonsWithConversionsAnalyzer.h.

Referenced by analyze(), and beginJob().

TH1F* PhotonsWithConversionsAnalyzer::h_MCConvPt_ [private]

Definition at line 76 of file PhotonsWithConversionsAnalyzer.h.

Referenced by beginJob().

TH1F* PhotonsWithConversionsAnalyzer::h_MCphoE_ [private]

Definition at line 68 of file PhotonsWithConversionsAnalyzer.h.

Referenced by analyze(), and beginJob().

TH1F* PhotonsWithConversionsAnalyzer::h_MCphoEta_ [private]

Definition at line 70 of file PhotonsWithConversionsAnalyzer.h.

Referenced by analyze(), and beginJob().

TH1F* PhotonsWithConversionsAnalyzer::h_MCphoPhi_ [private]

Definition at line 69 of file PhotonsWithConversionsAnalyzer.h.

Referenced by analyze(), and beginJob().

TH1F* PhotonsWithConversionsAnalyzer::h_phoE_ [private]

Definition at line 86 of file PhotonsWithConversionsAnalyzer.h.

Referenced by analyze(), and beginJob().

TH1F* PhotonsWithConversionsAnalyzer::h_phoEta_ [private]

Definition at line 87 of file PhotonsWithConversionsAnalyzer.h.

Referenced by analyze(), and beginJob().

TH1F* PhotonsWithConversionsAnalyzer::h_phoPhi_ [private]

Definition at line 88 of file PhotonsWithConversionsAnalyzer.h.

Referenced by analyze(), and beginJob().

TH1F* PhotonsWithConversionsAnalyzer::h_scE_ [private]

Definition at line 81 of file PhotonsWithConversionsAnalyzer.h.

Referenced by analyze(), and beginJob().

TH1F* PhotonsWithConversionsAnalyzer::h_scEt_ [private]

Definition at line 82 of file PhotonsWithConversionsAnalyzer.h.

Referenced by analyze(), and beginJob().

TH1F* PhotonsWithConversionsAnalyzer::h_scEta_ [private]

Definition at line 83 of file PhotonsWithConversionsAnalyzer.h.

Referenced by analyze(), and beginJob().

TH1F* PhotonsWithConversionsAnalyzer::h_scPhi_ [private]

Definition at line 84 of file PhotonsWithConversionsAnalyzer.h.

Referenced by analyze(), and beginJob().

std::string PhotonsWithConversionsAnalyzer::HepMCLabel [private]

Definition at line 53 of file PhotonsWithConversionsAnalyzer.h.

int PhotonsWithConversionsAnalyzer::nEvt_ [private]

Definition at line 49 of file PhotonsWithConversionsAnalyzer.h.

Referenced by analyze(), beginJob(), and endJob().

int PhotonsWithConversionsAnalyzer::nMatched_ [private]

Definition at line 51 of file PhotonsWithConversionsAnalyzer.h.

Referenced by analyze(), and beginJob().

int PhotonsWithConversionsAnalyzer::nMCPho_ [private]

Definition at line 50 of file PhotonsWithConversionsAnalyzer.h.

Referenced by analyze(), and beginJob().

std::string PhotonsWithConversionsAnalyzer::photonCollection_ [private]

Definition at line 60 of file PhotonsWithConversionsAnalyzer.h.

Referenced by analyze(), and PhotonsWithConversionsAnalyzer().

std::string PhotonsWithConversionsAnalyzer::photonCollectionProducer_ [private]

Definition at line 59 of file PhotonsWithConversionsAnalyzer.h.

Referenced by analyze(), and PhotonsWithConversionsAnalyzer().

std::string PhotonsWithConversionsAnalyzer::SimHitLabel [private]

Definition at line 56 of file PhotonsWithConversionsAnalyzer.h.

std::string PhotonsWithConversionsAnalyzer::SimTkLabel [private]

Definition at line 54 of file PhotonsWithConversionsAnalyzer.h.

std::string PhotonsWithConversionsAnalyzer::SimVtxLabel [private]

Definition at line 55 of file PhotonsWithConversionsAnalyzer.h.

PhotonMCTruthFinder* PhotonsWithConversionsAnalyzer::thePhotonMCTruthFinder_ [private]

Definition at line 43 of file PhotonsWithConversionsAnalyzer.h.

Referenced by analyze(), beginJob(), and ~PhotonsWithConversionsAnalyzer().


The documentation for this class was generated from the following files:
Generated on Tue Jun 9 18:29:55 2009 for CMSSW by  doxygen 1.5.4