CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_0/src/RecoEgamma/Examples/plugins/MCElectronAnalyzer.cc

Go to the documentation of this file.
00001 #include <iostream>
00002 //
00003 #include "RecoEgamma/Examples/plugins/MCElectronAnalyzer.h"
00004 #include "RecoEgamma/EgammaMCTools/interface/ElectronMCTruthFinder.h"
00005 #include "RecoEgamma/EgammaMCTools/interface/ElectronMCTruth.h"
00006 //
00007 #include "SimDataFormats/GeneratorProducts/interface/HepMCProduct.h"
00008 #include "SimDataFormats/Track/interface/SimTrack.h"
00009 #include "SimDataFormats/Track/interface/SimTrackContainer.h"
00010 #include "SimDataFormats/Vertex/interface/SimVertex.h"
00011 #include "SimDataFormats/Vertex/interface/SimVertexContainer.h"
00012 //
00013 #include "SimDataFormats/CrossingFrame/interface/CrossingFrame.h"
00014 #include "SimDataFormats/CrossingFrame/interface/MixCollection.h"
00015 #include "SimDataFormats/TrackingAnalysis/interface/TrackingParticleFwd.h"
00016 #include "SimDataFormats/TrackingAnalysis/interface/TrackingVertexContainer.h"
00017 //
00018 #include "DataFormats/Common/interface/Handle.h"
00019 //
00020 #include "FWCore/Framework/interface/Event.h"
00021 #include "FWCore/Framework/interface/EventSetup.h"
00022 #include "FWCore/Framework/interface/ESHandle.h"
00023 #include "FWCore/Framework/interface/MakerMacros.h"
00024 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00025 #include "FWCore/Utilities/interface/Exception.h"
00026 //
00027 #include "MagneticField/Records/interface/IdealMagneticFieldRecord.h"
00028 #include "MagneticField/Engine/interface/MagneticField.h"
00029 
00030 //
00031 #include "TFile.h"
00032 #include "TH1.h"
00033 #include "TH2.h"
00034 #include "TTree.h"
00035 #include "TVector3.h"
00036 #include "TProfile.h"
00037 //
00038 
00039 
00040 
00041 using namespace std;
00042 
00043 
00044 MCElectronAnalyzer::MCElectronAnalyzer( const edm::ParameterSet& pset )
00045    : fOutputFileName_( pset.getUntrackedParameter<string>("HistOutFile",std::string("TestConversions.root")) ),
00046      fOutputFile_(0)
00047 {
00048 
00049 
00050 }
00051 
00052 
00053 
00054 MCElectronAnalyzer::~MCElectronAnalyzer() {
00055 
00056 
00057   delete theElectronMCTruthFinder_;
00058 
00059 }
00060 
00061 
00062 void MCElectronAnalyzer::beginJob()
00063 {
00064 
00065 
00066   nEvt_=0;
00067 
00068   theElectronMCTruthFinder_ = new ElectronMCTruthFinder();
00069 
00070 
00071   fOutputFile_   = new TFile( fOutputFileName_.c_str(), "RECREATE" ) ;
00072 
00074   h_MCEleE_ = new TH1F("MCEleE","MC ele energy",100,0.,200.);
00075   h_MCElePhi_ = new TH1F("MCElePhi","MC ele phi",40,-3.14, 3.14);
00076   h_MCEleEta_ = new TH1F("MCEleEta","MC ele eta",40,-3., 3.);
00077   h_BremFrac_ = new TH1F("bremFrac","brem frac ", 100, 0., 1.);
00078   h_BremEnergy_ = new TH1F("BremE","Brem energy",100,0.,200.);
00079 
00080   p_BremVsR_ = new TProfile("BremVsR", " Mean Brem energy vs R ", 48, 0., 120.);
00081   p_BremVsEta_ = new TProfile("BremVsEta", " Mean Brem energy vs Eta ", 50, -2.5, 2.5);
00082 
00083 
00084 
00085 
00086   return ;
00087 }
00088 
00089 
00090 float MCElectronAnalyzer::etaTransformation(  float EtaParticle , float Zvertex)  {
00091 
00092 //---Definitions
00093         const float PI    = 3.1415927;
00094         //const float TWOPI = 2.0*PI;
00095 
00096 //---Definitions for ECAL
00097         const float R_ECAL           = 136.5;
00098         const float Z_Endcap         = 328.0;
00099         const float etaBarrelEndcap  = 1.479;
00100 
00101 //---ETA correction
00102 
00103         float Theta = 0.0  ;
00104         float ZEcal = R_ECAL*sinh(EtaParticle)+Zvertex;
00105 
00106         if(ZEcal != 0.0) Theta = atan(R_ECAL/ZEcal);
00107         if(Theta<0.0) Theta = Theta+PI ;
00108         float ETA = - log(tan(0.5*Theta));
00109 
00110         if( std::abs(ETA) > etaBarrelEndcap )
00111           {
00112            float Zend = Z_Endcap ;
00113            if(EtaParticle<0.0 )  Zend = -Zend ;
00114            float Zlen = Zend - Zvertex ;
00115            float RR = Zlen/sinh(EtaParticle);
00116            Theta = atan(RR/Zend);
00117            if(Theta<0.0) Theta = Theta+PI ;
00118            ETA = - log(tan(0.5*Theta));
00119           }
00120 //---Return the result
00121         return ETA;
00122 //---end
00123 }
00124 
00125 float MCElectronAnalyzer::phiNormalization(float & phi)
00126 {
00127 //---Definitions
00128  const float PI    = 3.1415927;
00129  const float TWOPI = 2.0*PI;
00130 
00131 
00132  if(phi >  PI) {phi = phi - TWOPI;}
00133  if(phi < -PI) {phi = phi + TWOPI;}
00134 
00135  //  cout << " Float_t PHInormalization out " << PHI << endl;
00136  return phi;
00137 
00138 }
00139 
00140 
00141 
00142 
00143 
00144 void MCElectronAnalyzer::analyze( const edm::Event& e, const edm::EventSetup& )
00145 {
00146 
00147 
00148   using namespace edm;
00149   //const float etaPhiDistance=0.01;
00150   // Fiducial region
00151   //const float TRK_BARL =0.9;
00152   //const float BARL = 1.4442; // DAQ TDR p.290
00153   //const float END_LO = 1.566;
00154   //const float END_HI = 2.5;
00155   // Electron mass
00156   //const Float_t mElec= 0.000511;
00157 
00158 
00159   nEvt_++;
00160   LogInfo("MCElectronAnalyzer") << "MCElectronAnalyzer Analyzing event number: " << e.id() << " Global Counter " << nEvt_ <<"\n";
00161   //  LogDebug("MCElectronAnalyzer") << "MCElectronAnalyzer Analyzing event number: "  << e.id() << " Global Counter " << nEvt_ <<"\n";
00162   std::cout << "MCElectronAnalyzer Analyzing event number: "  << e.id() << " Global Counter " << nEvt_ <<"\n";
00163 
00165   std::cout  << " MCElectronAnalyzer Looking for MC truth " << "\n";
00166 
00167   //get simtrack info
00168   std::vector<SimTrack> theSimTracks;
00169   std::vector<SimVertex> theSimVertices;
00170 
00171   edm::Handle<SimTrackContainer> SimTk;
00172   edm::Handle<SimVertexContainer> SimVtx;
00173   e.getByLabel("g4SimHits",SimTk);
00174   e.getByLabel("g4SimHits",SimVtx);
00175 
00176   theSimTracks.insert(theSimTracks.end(),SimTk->begin(),SimTk->end());
00177   theSimVertices.insert(theSimVertices.end(),SimVtx->begin(),SimVtx->end());
00178   std::cout << " MCElectronAnalyzer This Event has " <<  theSimTracks.size() << " sim tracks " << std::endl;
00179   std::cout << " MCElectronAnalyzer This Event has " <<  theSimVertices.size() << " sim vertices " << std::endl;
00180   if (  ! theSimTracks.size() ) std::cout << " Event number " << e.id() << " has NO sim tracks " << std::endl;
00181 
00182 
00183   std::vector<ElectronMCTruth> MCElectronctrons=theElectronMCTruthFinder_->find (theSimTracks,  theSimVertices);
00184   std::cout << " MCElectronAnalyzer MCElectronctrons size " <<  MCElectronctrons.size() << std::endl;
00185 
00186   for ( std::vector<ElectronMCTruth>::const_iterator iEl=MCElectronctrons.begin(); iEl !=MCElectronctrons.end(); ++iEl ){
00187 
00188     h_MCEleE_->Fill  ( (*iEl).fourMomentum().e() );
00189     h_MCEleEta_->Fill  ( (*iEl).fourMomentum().pseudoRapidity() );
00190     h_MCElePhi_->Fill  ( (*iEl).fourMomentum().phi() );
00191 
00192     float totBrem=0 ;
00193     unsigned int iBrem ;
00194     for ( iBrem=0; iBrem < (*iEl).bremVertices().size(); ++iBrem ) {
00195       float rBrem= (*iEl).bremVertices()[iBrem].perp();
00196       float etaBrem= (*iEl).bremVertices()[iBrem].eta();
00197       if ( rBrem < 120 ) {
00198         totBrem +=  (*iEl).bremMomentum()[iBrem].e();
00199         p_BremVsR_ ->Fill ( rBrem, (*iEl).bremMomentum()[iBrem].e() );
00200         p_BremVsEta_ ->Fill ( etaBrem, (*iEl).bremMomentum()[iBrem].e() );
00201 
00202       }
00203     }
00204 
00205 
00206     h_BremFrac_->Fill( totBrem/(*iEl).fourMomentum().e() );
00207     h_BremEnergy_->Fill (  totBrem  );
00208 
00209 
00210 
00211 
00212   }
00213 
00214 
00215 
00216 
00217 }
00218 
00219 
00220 
00221 
00222 void MCElectronAnalyzer::endJob()
00223 {
00224 
00225 
00226 
00227 
00228    fOutputFile_->Write() ;
00229    fOutputFile_->Close() ;
00230 
00231    edm::LogInfo("MCElectronAnalyzer") << "Analyzed " << nEvt_  << "\n";
00232    std::cout  << "MCElectronAnalyzer::endJob Analyzed " << nEvt_ << " events " << "\n";
00233 
00234    return ;
00235 }
00236