CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_1_8_patch12/src/HLTrigger/HLTanalyzers/src/HLTEgamma.cc

Go to the documentation of this file.
00001 #include <iostream>
00002 #include <sstream>
00003 #include <istream>
00004 #include <fstream>
00005 #include <iomanip>
00006 #include <string>
00007 #include <cmath>
00008 #include <functional>
00009 #include <cstdlib>
00010 #include <cstring>
00011 
00012 #include "DataFormats/EgammaReco/interface/ElectronSeed.h"
00013 #include "DataFormats/EgammaReco/interface/ElectronSeedFwd.h"
00014 #include "DataFormats/GsfTrackReco/interface/GsfTrack.h" 
00015 #include "HLTrigger/HLTanalyzers/interface/HLTEgamma.h"
00016 
00017 
00018 
00019 #include "RecoEgamma/EgammaTools/interface/ECALPositionCalculator.h"
00020 
00021 #include "FWCore/Framework/interface/EventSetup.h"
00022 #include "MagneticField/Records/interface/IdealMagneticFieldRecord.h"
00023 #include "MagneticField/Engine/interface/MagneticField.h"
00024 
00025 #include "DataFormats/EgammaReco/interface/SuperCluster.h"
00026 
00027 #include "HLTMessages.h"
00028 
00029 static const size_t kMaxEl     = 10000;
00030 static const size_t kMaxPhot   = 10000;
00031 static const size_t kMaxhPhot  =   500;
00032 static const size_t kMaxhEle   =   500;
00033 
00034 HLTEgamma::HLTEgamma() {
00035 }
00036 
00037 /*  Setup the analysis to put the branch-variables size_to the tree. */
00038 void HLTEgamma::setup(const edm::ParameterSet& pSet, TTree* HltTree)
00039 {
00040   elpt              = new float[kMaxEl];
00041   elphi             = new float[kMaxEl];
00042   eleta             = new float[kMaxEl];
00043   elet              = new float[kMaxEl];
00044   ele               = new float[kMaxEl];
00045   eleId             = new int[kMaxEl];// RL  + 2*RT + 4*L +  4*T 
00046   elIP              = new float[kMaxEl];  
00047   elNLostHits       = new int[kMaxEl];  
00048   elTrkChi2NDF      = new float[kMaxEl];        
00049   elTrkIsoR03       = new float[kMaxEl];  
00050   elECaloIsoR03     = new float[kMaxEl];  
00051   elHCaloIsoR03     = new float[kMaxEl];  
00052   elIsEcalDriven    = new bool[kMaxEl];  
00053   elFbrem           = new float[kMaxEl];      
00054   elmishits         = new int[kMaxEl]; 
00055   eldist            = new float[kMaxEl]; 
00056   eldcot            = new float[kMaxEl]; 
00057   eltrkiso          = new float[kMaxEl]; 
00058   elecaliso         = new float[kMaxEl]; 
00059   elhcaliso         = new float[kMaxEl]; 
00060   elsigmaietaieta   = new float[kMaxEl]; 
00061   eldeltaPhiIn      = new float[kMaxEl]; 
00062   eldeltaEtaIn      = new float[kMaxEl]; 
00063   elhOverE          = new float[kMaxEl]; 
00064   elscEt            = new float[kMaxEl]; 
00065   eld0corr          = new float[kMaxEl]; 
00066   elqGsfCtfScPixConsistent = new bool[kMaxEl]; 
00067   
00068   photonpt          = new float[kMaxPhot];
00069   photonphi         = new float[kMaxPhot];
00070   photoneta         = new float[kMaxPhot];
00071   photonet          = new float[kMaxPhot];
00072   photone           = new float[kMaxPhot];
00073 
00074   hphotet           = new float[kMaxhPhot];
00075   hphoteta          = new float[kMaxhPhot];
00076   hphotphi          = new float[kMaxhPhot];
00077   hphoteiso         = new float[kMaxhPhot];
00078   hphothiso         = new float[kMaxhPhot];
00079   hphottiso         = new float[kMaxhPhot];
00080   hphotl1iso        = new int[kMaxhPhot];
00081   hphotClusShap     = new float[kMaxhPhot];
00082   hphotR9           = new float[kMaxhPhot]; 
00083   hphothovereh      = new float[kMaxhPhot];
00084   hphotR9ID         = new float[kMaxhPhot];
00085 
00086   heleet            = new float[kMaxhEle];
00087   heleeta           = new float[kMaxhEle];
00088   helephi           = new float[kMaxhEle];
00089   heleE             = new float[kMaxhEle];
00090   helep             = new float[kMaxhEle];
00091   helehiso          = new float[kMaxhEle];
00092   heleeiso          = new float[kMaxhEle];
00093   heletiso          = new float[kMaxhEle];
00094   helel1iso         = new int[kMaxhEle];
00095   helePixelSeeds    = new int[kMaxhEle];
00096   heleNewSC         = new int[kMaxhEle];
00097   heleClusShap      = new float[kMaxhEle];
00098   heleDeta          = new float[kMaxhEle];
00099   heleDphi          = new float[kMaxhEle];
00100   heleR9            = new float[kMaxhEle]; 
00101   helehovereh       = new float[kMaxhEle];
00102   heleR9ID          = new float[kMaxhEle];
00103 
00104   hhfelept         = new float[kMaxhEle];
00105   hhfeleeta        = new float[kMaxhEle]; 
00106   hhfclustere9e25  = new float[kMaxhEle]; 
00107   hhfcluster2Dcut  = new float[kMaxhEle]; 
00108 
00109   nele        = 0;
00110   nphoton     = 0;
00111   nhltgam     = 0;
00112   nhltele     = 0;
00113   nhlthfele   = 0;
00114   nhlthfeclus = 0;
00115 
00116   // Egamma-specific branches of the tree
00117   HltTree->Branch("NrecoElec",          & nele,             "NrecoElec/I");
00118   HltTree->Branch("recoElecPt",         elpt,               "recoElecPt[NrecoElec]/F");
00119   HltTree->Branch("recoElecPhi",        elphi,              "recoElecPhi[NrecoElec]/F");
00120   HltTree->Branch("recoElecEta",        eleta,              "recoElecEta[NrecoElec]/F");
00121   HltTree->Branch("recoElecEt",         elet,               "recoElecEt[NrecoElec]/F");
00122   HltTree->Branch("recoElecE",          ele,                "recoElecE[NrecoElec]/F");
00123   HltTree->Branch("recoElecEleID",      eleId,              "recoElecEleID[NrecoElec]/I");
00124   HltTree->Branch("recoElecIP",           elIP,            "recoElecIP[NrecoElec]/F");  
00125   HltTree->Branch("recoElecNLostHits",    elNLostHits,     "recoElecNLostHits[NrecoElec]/I");  
00126   HltTree->Branch("recoElecChi2NDF",      elTrkChi2NDF,    "recoElecChi2NDF[NrecoElec]/F");  
00127   HltTree->Branch("recoElecTrkIsoR03",    elTrkIsoR03,     "recoElecTrkIsoR03[NrecoElec]/F");  
00128   HltTree->Branch("recoElecECaloIsoR03",  elECaloIsoR03,   "recoElecECaloIsoR03[NrecoElec]/F");  
00129   HltTree->Branch("recoElecHCaloIsoR03",  elHCaloIsoR03,   "recoElecHCaloIsoR03[NrecoElec]/F");  
00130   HltTree->Branch("recoElecIsEcalDriven", elIsEcalDriven,  "recoElecIsEcalDriven[NrecoElec]/O");        
00131   HltTree->Branch("recoElecFbrem",        elFbrem,         "recoElecFbrem[NrecoElec]/F");  
00132   HltTree->Branch("recoElecmishits",                  elmishits,                "recoElecmishits[NrecoElec]/I"); 
00133   HltTree->Branch("recoElecdist",                     eldist,                   "recoElecdist[NrecoElec]/F"); 
00134   HltTree->Branch("recoElecdcot",                     eldcot,                   "recoElecdcot[NrecoElec]/F"); 
00135   HltTree->Branch("recoElectrkiso",                   eltrkiso,                 "recoElectrkiso[NrecoElec]/F"); 
00136   HltTree->Branch("recoElececaliso",                  elecaliso,                "recoElececaliso[NrecoElec]/F"); 
00137   HltTree->Branch("recoElechcaliso",                  elhcaliso,                "recoElechcaliso[NrecoElec]/F"); 
00138   HltTree->Branch("recoElecsigmaietaieta",            elsigmaietaieta,          "recoElecsigmaietaieta[NrecoElec]/F"); 
00139   HltTree->Branch("recoElecdeltaPhiIn",               eldeltaPhiIn,             "recoElecdeltaPhiIn[NrecoElec]/F"); 
00140   HltTree->Branch("recoElecdeltaEtaIn",               eldeltaEtaIn,             "recoElecdeltaEtaIn[NrecoElec]/F"); 
00141   HltTree->Branch("recoElechOverE",                   elhOverE,                 "recoElechOverE[NrecoElec]/F"); 
00142   HltTree->Branch("recoElecscEt",                     elscEt,                   "recoElecscEt[NrecoElec]/F"); 
00143   HltTree->Branch("recoElecd0corr",                   eld0corr,                 "recoElecd0corr[NrecoElec]/F"); 
00144   HltTree->Branch("recoElecqGsfCtfScPixConsistent",   elqGsfCtfScPixConsistent, "recoElecqGsfCtfScPixConsistent[NrecoElec]/O");  
00145 
00146   HltTree->Branch("NrecoPhot",          &nphoton,           "NrecoPhot/I");
00147   HltTree->Branch("recoPhotPt",         photonpt,           "recoPhotPt[NrecoPhot]/F");
00148   HltTree->Branch("recoPhotPhi",        photonphi,          "recoPhotPhi[NrecoPhot]/F");
00149   HltTree->Branch("recoPhotEta",        photoneta,          "recoPhotEta[NrecoPhot]/F");
00150   HltTree->Branch("recoPhotEt",         photonet,           "recoPhotEt[NrecoPhot]/F");
00151   HltTree->Branch("recoPhotE",          photone,            "recoPhotE[NrecoPhot]/F");
00152 
00153   HltTree->Branch("NohPhot",            & nhltgam,          "NohPhot/I");
00154   HltTree->Branch("ohPhotEt",           hphotet,            "ohPhotEt[NohPhot]/F");
00155   HltTree->Branch("ohPhotEta",          hphoteta,           "ohPhotEta[NohPhot]/F");
00156   HltTree->Branch("ohPhotPhi",          hphotphi,           "ohPhotPhi[NohPhot]/F");
00157   HltTree->Branch("ohPhotEiso",         hphoteiso,          "ohPhotEiso[NohPhot]/F");
00158   HltTree->Branch("ohPhotHiso",         hphothiso,          "ohPhotHiso[NohPhot]/F");
00159   HltTree->Branch("ohPhotTiso",         hphottiso,          "ohPhotTiso[NohPhot]/F");
00160   HltTree->Branch("ohPhotL1iso",        hphotl1iso,         "ohPhotL1iso[NohPhot]/I");
00161   HltTree->Branch("ohPhotClusShap",     hphotClusShap,      "ohPhotClusShap[NohPhot]/F");
00162   HltTree->Branch("ohPhotR9",           hphotR9,            "ohPhotR9[NohPhot]/F");  
00163   HltTree->Branch("ohPhotHforHoverE",   hphothovereh,       "ohPhotHforHoverE[NohPhot]/F");   
00164   HltTree->Branch("ohPhotR9ID",         hphotR9ID,          "ohPhotR9ID[NohPhot]/F");
00165 
00166   HltTree->Branch("NohEle",             & nhltele,          "NohEle/I");
00167   HltTree->Branch("ohEleEt",            heleet,             "ohEleEt[NohEle]/F");
00168   HltTree->Branch("ohEleEta",           heleeta,            "ohEleEta[NohEle]/F");
00169   HltTree->Branch("ohElePhi",           helephi,            "ohElePhi[NohEle]/F");
00170   HltTree->Branch("ohEleE",             heleE,              "ohEleE[NohEle]/F");
00171   HltTree->Branch("ohEleP",             helep,              "ohEleP[NohEle]/F");
00172   HltTree->Branch("ohEleHiso",          helehiso,           "ohEleHiso[NohEle]/F");
00173   HltTree->Branch("ohEleTiso",          heletiso,           "ohEleTiso[NohEle]/F");
00174   HltTree->Branch("ohEleEiso",          heleeiso,           "ohEleEiso[NohEle]/F"); 
00175   HltTree->Branch("ohEleL1iso",         helel1iso,          "ohEleLiso[NohEle]/I");
00176   HltTree->Branch("ohElePixelSeeds",    helePixelSeeds,     "ohElePixelSeeds[NohEle]/I");
00177   HltTree->Branch("ohEleNewSC",         heleNewSC,          "ohEleNewSC[NohEle]/I");
00178   HltTree->Branch("ohEleClusShap",      heleClusShap,       "ohEleClusShap[NohEle]/F");
00179   HltTree->Branch("ohEleDeta",          heleDeta,           "ohEleDeta[NohEle]/F");
00180   HltTree->Branch("ohEleDphi",          heleDphi,           "ohEleDphi[NohEle]/F");
00181   HltTree->Branch("ohEleR9",            heleR9,             "ohEleR9[NohEle]/F");  
00182   HltTree->Branch("ohEleHforHoverE",    helehovereh,        "ohEleHforHoverE[NohEle]/F");    
00183   HltTree->Branch("ohEleR9ID",          heleR9ID,           "ohEleR9ID[NohEle]/F");
00184   HltTree->Branch("NohHFEle",           &nhlthfele ,        "NohHFEle/I"); 
00185   HltTree->Branch("ohHFElePt",          hhfelept,           "ohHFElePt[NohHFEle]/F");
00186   HltTree->Branch("ohHFEleEta",         hhfeleeta,          "ohHFElePt[NohHFEle]/F");  
00187   HltTree->Branch("NohHFECALClus",      &nhlthfeclus,       "NohHFECALClus/I"); 
00188   HltTree->Branch("ohHFElee9e25",       hhfclustere9e25,    "ohHFElePt[NohHFECALClus]/F");  
00189   HltTree->Branch("ohHFEle2Dcut",       hhfcluster2Dcut,    "ohHFElePt[NohHFECALClus]/F");  
00190  
00191 }
00192 
00193 void HLTEgamma::clear(void)
00194 {
00195   std::memset(elpt,             '\0', kMaxEl     * sizeof(float));
00196   std::memset(elphi,            '\0', kMaxEl     * sizeof(float));
00197   std::memset(eleta,            '\0', kMaxEl     * sizeof(float));
00198   std::memset(elet,             '\0', kMaxEl     * sizeof(float));
00199   std::memset(ele,              '\0', kMaxEl     * sizeof(float));
00200   std::memset(ele,              '\0', kMaxEl     * sizeof(int));
00201   std::memset(elIP,             '\0', kMaxEl     * sizeof(float));  
00202   std::memset(elNLostHits,      '\0', kMaxEl     * sizeof(int));  
00203   std::memset(elTrkChi2NDF,     '\0', kMaxEl     * sizeof(float));     
00204   std::memset(elTrkIsoR03,      '\0', kMaxEl     * sizeof(float));  
00205   std::memset(elECaloIsoR03,    '\0', kMaxEl     * sizeof(float));  
00206   std::memset(elHCaloIsoR03,    '\0', kMaxEl     * sizeof(float));  
00207   std::memset(elIsEcalDriven,   '\0', kMaxEl     * sizeof(bool));  
00208   std::memset(elFbrem,          '\0', kMaxEl     * sizeof(float));  
00209 
00210   std::memset(photonpt,         '\0', kMaxPhot   * sizeof(float));
00211   std::memset(photonphi,        '\0', kMaxPhot   * sizeof(float));
00212   std::memset(photoneta,        '\0', kMaxPhot   * sizeof(float));
00213   std::memset(photonet,         '\0', kMaxPhot   * sizeof(float));
00214   std::memset(photone,          '\0', kMaxPhot   * sizeof(float));
00215 
00216   std::memset(hphotet,          '\0', kMaxhPhot  * sizeof(float));
00217   std::memset(hphoteta,         '\0', kMaxhPhot  * sizeof(float));
00218   std::memset(hphotphi,         '\0', kMaxhPhot  * sizeof(float));
00219   std::memset(hphoteiso,        '\0', kMaxhPhot  * sizeof(float));
00220   std::memset(hphothiso,        '\0', kMaxhPhot  * sizeof(float));
00221   std::memset(hphottiso,        '\0', kMaxhPhot  * sizeof(float));
00222   std::memset(hphotl1iso,       '\0', kMaxhPhot  * sizeof(int));
00223   std::memset(hphotClusShap,    '\0', kMaxhPhot  * sizeof(float));
00224 
00225   std::memset(heleet,           '\0', kMaxhEle   * sizeof(float));
00226   std::memset(heleeta,          '\0', kMaxhEle   * sizeof(float));
00227   std::memset(helephi,          '\0', kMaxhEle   * sizeof(float));
00228   std::memset(heleE,            '\0', kMaxhEle   * sizeof(float));
00229   std::memset(helep,            '\0', kMaxhEle   * sizeof(float));
00230   std::memset(helehiso,         '\0', kMaxhEle   * sizeof(float));
00231   std::memset(heletiso,         '\0', kMaxhEle   * sizeof(float));
00232   std::memset(heleeiso,         '\0', kMaxhEle   * sizeof(float)); 
00233   std::memset(helehovereh,      '\0', kMaxhEle   * sizeof(float));
00234   std::memset(helel1iso,        '\0', kMaxhEle   * sizeof(int));
00235   std::memset(helePixelSeeds,   '\0', kMaxhEle   * sizeof(int));
00236   std::memset(heleNewSC,        '\0', kMaxhEle   * sizeof(int));
00237   std::memset(heleClusShap,     '\0', kMaxhEle  * sizeof(float));
00238   std::memset(heleDeta,         '\0', kMaxhEle  * sizeof(float));
00239   std::memset(heleDphi,         '\0', kMaxhEle  * sizeof(float));
00240 
00241   std::memset(hhfelept,         '\0', kMaxhEle  * sizeof(float));
00242   std::memset(hhfeleeta,        '\0', kMaxhEle  * sizeof(float));
00243   std::memset(hhfclustere9e25,  '\0', kMaxhEle  * sizeof(float));
00244   std::memset(hhfcluster2Dcut,  '\0', kMaxhEle  * sizeof(float));
00245 
00246 
00247   nele      = 0;
00248   nphoton   = 0;
00249   nhltgam   = 0;
00250   nhltele   = 0;
00251   nhlthfele   = 0; 
00252   nhlthfeclus = 0; 
00253 }
00254 
00255 /* **Analyze the event** */
00256 void HLTEgamma::analyze(const edm::Handle<reco::GsfElectronCollection>         & electrons,
00257                         const edm::Handle<reco::PhotonCollection>              & photons,
00258                         const edm::Handle<reco::ElectronCollection>            & electronIsoHandle,
00259                         const edm::Handle<reco::ElectronCollection>            & electronNonIsoHandle,
00260                         const edm::Handle<reco::ElectronIsolationMap>          & NonIsoTrackEleIsolMap,
00261                         const edm::Handle<reco::ElectronIsolationMap>          & TrackEleIsolMap,
00262                         const edm::Handle<reco::ElectronSeedCollection>        & L1IsoPixelSeedsMap,
00263                         const edm::Handle<reco::ElectronSeedCollection>        & L1NonIsoPixelSeedsMap,
00264                         const edm::Handle<reco::RecoEcalCandidateCollection>   & recoIsolecalcands,
00265                         const edm::Handle<reco::RecoEcalCandidateCollection>   & recoNonIsolecalcands,
00266                         const edm::Handle<reco::RecoEcalCandidateIsolationMap> & EcalIsolMap,
00267                         const edm::Handle<reco::RecoEcalCandidateIsolationMap> & EcalNonIsolMap,
00268                         const edm::Handle<reco::RecoEcalCandidateIsolationMap> & HcalEleIsolMap,
00269                         const edm::Handle<reco::RecoEcalCandidateIsolationMap> & HcalEleNonIsolMap,
00270                         const edm::Handle<reco::RecoEcalCandidateIsolationMap> & HcalIsolMap,
00271                         const edm::Handle<reco::RecoEcalCandidateIsolationMap> & HcalNonIsolMap,
00272                         const edm::Handle<reco::RecoEcalCandidateIsolationMap> & TrackIsolMap,
00273                         const edm::Handle<reco::RecoEcalCandidateIsolationMap> & TrackNonIsolMap,
00274                         EcalClusterLazyTools& lazyTools,
00275                         const edm::ESHandle<MagneticField>& theMagField,
00276                         reco::BeamSpot::Point & BSPosition, 
00277                         std::vector<edm::Handle<edm::ValueMap<float> > > & eIDValueMap,
00278                         const edm::Handle<reco::RecoEcalCandidateIsolationMap> & photonR9IsoMap,   
00279                         const edm::Handle<reco::RecoEcalCandidateIsolationMap> & photonR9NonIsoMap,   
00280                         const edm::Handle<reco::RecoEcalCandidateIsolationMap> & electronR9IsoMap,   
00281                         const edm::Handle<reco::RecoEcalCandidateIsolationMap> & electronR9NonIsoMap,   
00282                         const edm::Handle<reco::RecoEcalCandidateIsolationMap> & photonHoverEHIsoMap,      
00283                         const edm::Handle<reco::RecoEcalCandidateIsolationMap> & photonHoverEHNonIsoMap,       
00284                         const edm::Handle<reco::RecoEcalCandidateIsolationMap> & photonR9IDIsoMap,
00285                         const edm::Handle<reco::RecoEcalCandidateIsolationMap> & photonR9IDNonIsoMap,
00286                         const edm::Handle<reco::RecoEcalCandidateIsolationMap> & electronR9IDIsoMap,
00287                         const edm::Handle<reco::RecoEcalCandidateIsolationMap> & electronR9IDNonIsoMap,
00288                         const edm::Handle<reco::SuperClusterCollection>        & electronHFECALClusters,  
00289                         const edm::Handle<reco::RecoEcalCandidateCollection>   & electronHFElectrons,   
00290                         TTree* HltTree)
00291 {
00292   // reset the tree variables
00293   clear();
00294 
00295   if (electrons.isValid()) {  
00296     reco::GsfElectronCollection myelectrons( electrons->begin(), electrons->end() );  
00297     nele = myelectrons.size();  
00298     std::sort(myelectrons.begin(), myelectrons.end(), EtGreater());  
00299     int iel = 0;  
00300     for (reco::GsfElectronCollection::const_iterator i = myelectrons.begin(); i != myelectrons.end(); i++) {  
00301       elpt[iel]  = i->pt();  
00302       elphi[iel] = i->phi();  
00303       eleta[iel] = i->eta();  
00304       elet[iel]  = i->et();  
00305       ele[iel]   = i->energy();  
00306         
00307       if(i->gsfTrack().isNonnull()){  
00308         elNLostHits[iel]   = i->gsfTrack()->trackerExpectedHitsInner().numberOfLostHits();  
00309         elIP[iel]          = i->gsfTrack()->dxy(BSPosition);    
00310         elTrkChi2NDF[iel]  = i->gsfTrack()->normalizedChi2();  
00311       }  
00312       else {  
00313         elNLostHits[iel]  = -99.;  
00314         elIP[iel]         = -99.;  
00315         elTrkChi2NDF[iel] = -99.;  
00316       }  
00317         
00318       elTrkIsoR03[iel]     = i->dr03TkSumPt();  
00319       elECaloIsoR03[iel]   = i->dr03EcalRecHitSumEt();  
00320       elHCaloIsoR03[iel]   = i->dr03HcalTowerSumEt();  
00321       elIsEcalDriven[iel]  = i->ecalDrivenSeed();  
00322       elFbrem[iel]         = i->fbrem();  
00323       elscEt[iel] = i->superCluster()->energy()*sin((2*atan(exp(-i->superCluster()->eta())))); 
00324       elhOverE[iel] = i->hadronicOverEm(); 
00325       elsigmaietaieta[iel] = i->sigmaIetaIeta(); 
00326       eldeltaPhiIn[iel] = i->deltaPhiSuperClusterTrackAtVtx(); 
00327       eldeltaEtaIn[iel] = i->deltaEtaSuperClusterTrackAtVtx(); 
00328       elmishits[iel] = i->gsfTrack()->trackerExpectedHitsInner().numberOfHits(); 
00329       eltrkiso[iel] = i->dr03TkSumPt(); 
00330       elecaliso[iel] = i->dr03EcalRecHitSumEt(); 
00331       elhcaliso[iel] = i->dr03HcalTowerSumEt(); 
00332       eld0corr[iel]= - (i->gsfTrack()->dxy(BSPosition)); 
00333       elqGsfCtfScPixConsistent[iel]=i->isGsfCtfScPixChargeConsistent();; 
00334  
00335       // conversion info will be available after 3_10_X 
00336       eldist[iel] = 0;// fabs(i->convDist()); 
00337       eldcot[iel] = 0; //fabs(i->convDcot()); 
00338         
00339       iel++;  
00340     }  
00341   } else {  
00342     nele = 0;  
00343   }
00344 
00345   if (photons.isValid()) {
00346     reco::PhotonCollection myphotons(* photons);
00347     nphoton = myphotons.size();
00348     std::sort(myphotons.begin(), myphotons.end(), EtGreater());
00349     int ipho = 0;
00350     for (reco::PhotonCollection::const_iterator i = myphotons.begin(); i!= myphotons.end(); i++) {
00351       photonpt[ipho] = i->pt();
00352       photonphi[ipho] = i->phi();
00353       photoneta[ipho] = i->eta();
00354       photonet[ipho] = i->et();
00355       photone[ipho] = i->energy();
00356       ipho++;
00357     }
00358   } else {
00359     nphoton = 0;
00360   }
00361 
00363 
00365   std::vector<OpenHLTPhoton> theHLTPhotons;
00366   MakeL1IsolatedPhotons(
00367       theHLTPhotons,
00368       recoIsolecalcands,
00369       EcalIsolMap,
00370       HcalIsolMap,
00371       TrackIsolMap,
00372       photonR9IsoMap,
00373       photonHoverEHIsoMap,
00374       photonR9IDIsoMap,
00375       lazyTools);
00376   MakeL1NonIsolatedPhotons(
00377       theHLTPhotons,
00378       recoNonIsolecalcands,
00379       EcalNonIsolMap,
00380       HcalNonIsolMap,
00381       TrackNonIsolMap,
00382       photonR9NonIsoMap,
00383       photonHoverEHNonIsoMap,
00384       photonR9IDNonIsoMap,
00385       lazyTools);
00386 
00387   std::sort(theHLTPhotons.begin(), theHLTPhotons.end(), EtGreater());
00388   nhltgam = theHLTPhotons.size();
00389 
00390   for (int u = 0; u < nhltgam; u++) {
00391     hphotet[u]    = theHLTPhotons[u].Et;
00392     hphoteta[u]   = theHLTPhotons[u].eta;
00393     hphotphi[u]   = theHLTPhotons[u].phi;
00394     hphoteiso[u]  = theHLTPhotons[u].ecalIsol;
00395     hphothiso[u]  = theHLTPhotons[u].hcalIsol;
00396     hphottiso[u]  = theHLTPhotons[u].trackIsol;
00397     hphotl1iso[u] = theHLTPhotons[u].L1Isolated;
00398     hphotClusShap[u] = theHLTPhotons[u].clusterShape;
00399     hphothovereh[u] = theHLTPhotons[u].hovereh; 
00400     hphotR9[u] = theHLTPhotons[u].r9;
00401     hphotR9ID[u] = theHLTPhotons[u].r9ID;
00402    }
00403 
00405   std::vector<OpenHLTElectron> theHLTElectrons;
00406   MakeL1IsolatedElectrons(
00407       theHLTElectrons,
00408       electronIsoHandle,
00409       recoIsolecalcands,
00410       HcalEleIsolMap,
00411       L1IsoPixelSeedsMap,
00412       TrackEleIsolMap,
00413       electronR9IsoMap,
00414       photonHoverEHIsoMap, 
00415       EcalIsolMap, 
00416       electronR9IDIsoMap,
00417       lazyTools,
00418       theMagField,
00419       BSPosition);
00420   MakeL1NonIsolatedElectrons(
00421       theHLTElectrons,
00422       electronNonIsoHandle,
00423       recoNonIsolecalcands,
00424       HcalEleNonIsolMap,
00425       L1NonIsoPixelSeedsMap,
00426       NonIsoTrackEleIsolMap,
00427       electronR9NonIsoMap,  
00428       photonHoverEHNonIsoMap, 
00429       EcalNonIsolMap, 
00430       electronR9IDNonIsoMap,
00431       lazyTools,
00432       theMagField,
00433       BSPosition);
00434 
00435   std::sort(theHLTElectrons.begin(), theHLTElectrons.end(), EtGreater());
00436   nhltele = theHLTElectrons.size();
00437 
00438   for (int u = 0; u < nhltele; u++) {
00439     heleet[u]         = theHLTElectrons[u].Et;
00440     heleeta[u]        = theHLTElectrons[u].eta;
00441     helephi[u]        = theHLTElectrons[u].phi;
00442     heleE[u]          = theHLTElectrons[u].E;
00443     helep[u]          = theHLTElectrons[u].p;
00444     helehiso[u]       = theHLTElectrons[u].hcalIsol;
00445     helePixelSeeds[u] = theHLTElectrons[u].pixelSeeds;
00446     heletiso[u]       = theHLTElectrons[u].trackIsol;
00447     heleeiso[u]       = theHLTElectrons[u].ecalIsol; 
00448     helel1iso[u]      = theHLTElectrons[u].L1Isolated;
00449     heleNewSC[u]      = theHLTElectrons[u].newSC;
00450     heleClusShap[u]   = theHLTElectrons[u].clusterShape;
00451     heleDeta[u]       = theHLTElectrons[u].Deta;
00452     heleDphi[u]       = theHLTElectrons[u].Dphi;
00453     heleR9[u]         = theHLTElectrons[u].r9;  
00454     helehovereh[u]    = theHLTElectrons[u].hovereh;
00455     heleR9ID[u]       = theHLTElectrons[u].r9ID;
00456   }
00457 
00458   if(electronHFElectrons.isValid()) {
00459     for (reco::RecoEcalCandidateCollection::const_iterator hfelecand = electronHFElectrons->begin(); 
00460          hfelecand!= electronHFElectrons->end(); hfelecand++) { 
00461       hhfelept[nhlthfele] = hfelecand->pt();
00462       hhfeleeta[nhlthfele] = hfelecand->eta(); 
00463       nhlthfele++;
00464     }
00465   }
00466   if(electronHFECALClusters.isValid()) { 
00467     for(reco::SuperClusterCollection::const_iterator hfeleclus = electronHFECALClusters->begin();
00468         hfeleclus!= electronHFECALClusters->end(); hfeleclus++) {
00469       hhfclustere9e25[nhlthfeclus] = -999.;
00470       hhfcluster2Dcut[nhlthfeclus] = -999.;
00471       nhlthfeclus++;
00472     }
00473   }
00474 }
00475 
00476 void HLTEgamma::MakeL1IsolatedPhotons(
00477     std::vector<OpenHLTPhoton> & theHLTPhotons,
00478     const edm::Handle<reco::RecoEcalCandidateCollection>   & recoIsolecalcands,
00479     const edm::Handle<reco::RecoEcalCandidateIsolationMap> & EcalIsolMap,
00480     const edm::Handle<reco::RecoEcalCandidateIsolationMap> & HcalIsolMap,
00481     const edm::Handle<reco::RecoEcalCandidateIsolationMap> & TrackIsolMap,
00482     const edm::Handle<reco::RecoEcalCandidateIsolationMap> & photonR9IsoMap,    
00483     const edm::Handle<reco::RecoEcalCandidateIsolationMap> & photonHoverEHIsoMap,     
00484     const edm::Handle<reco::RecoEcalCandidateIsolationMap> & photonR9IDIsoMap,
00485     EcalClusterLazyTools& lazyTools )
00486 {
00487   // Iterator to the isolation-map
00488   reco::RecoEcalCandidateIsolationMap::const_iterator mapi;
00489 
00490   if (recoIsolecalcands.isValid()) {
00491     // loop over SuperCluster and fill the HLTPhotons
00492   
00493 
00494     for (reco::RecoEcalCandidateCollection::const_iterator recoecalcand = recoIsolecalcands->begin();
00495          recoecalcand!= recoIsolecalcands->end(); recoecalcand++) {
00496 
00497       OpenHLTPhoton pho;
00498       pho.ecalIsol   = -999;
00499       pho.hcalIsol   = -999;
00500       pho.trackIsol  = -999;
00501       pho.clusterShape = -999;
00502       pho.L1Isolated = true;
00503       pho.Et         = recoecalcand->et();
00504       pho.eta        = recoecalcand->eta();
00505       pho.phi        = recoecalcand->phi();
00506       pho.r9         = -999.;
00507       pho.hovereh    = -999.;
00508       pho.r9ID       = -999.;
00509 
00510       //Get the cluster shape
00511       //      std::vector<float> vCov = lazyTools.covariances( *(recoecalcand->superCluster()->seed()) );
00512       std::vector<float> vCov = lazyTools.localCovariances( *(recoecalcand->superCluster()->seed()) );
00513       double sigmaee = sqrt(vCov[0]);
00514       //      float EtaSC = fabs(recoecalcand->eta());
00515       //      if(EtaSC > 1.479 ) {//Endcap
00516       //        sigmaee = sigmaee - 0.02*(EtaSC - 2.3);
00517       //      }
00518       pho.clusterShape = sigmaee;
00519 
00520       // Method to get the reference to the candidate
00521       reco::RecoEcalCandidateRef ref = reco::RecoEcalCandidateRef(recoIsolecalcands, distance(recoIsolecalcands->begin(), recoecalcand));
00522 
00523       // First/Second member of the Map: Ref-to-Candidate(mapi)/Isolation(->val)
00524       // fill the ecal Isolation
00525       if (EcalIsolMap.isValid()) {
00526         mapi = (*EcalIsolMap).find(ref);
00527         if (mapi !=(*EcalIsolMap).end()) { pho.ecalIsol = mapi->val;}
00528       }
00529       // fill the hcal Isolation
00530       if (HcalIsolMap.isValid()) {
00531         mapi = (*HcalIsolMap).find(ref);
00532         if (mapi !=(*HcalIsolMap).end()) { pho.hcalIsol = mapi->val;}
00533       }
00534       // fill the track Isolation
00535       if (TrackIsolMap.isValid()) {
00536         mapi = (*TrackIsolMap).find(ref);
00537         if (mapi !=(*TrackIsolMap).end()) { pho.trackIsol = mapi->val;}
00538       }
00539       // fill the R9
00540       if (photonR9IsoMap.isValid()) {
00541         mapi = (*photonR9IsoMap).find(ref); 
00542         if (mapi !=(*photonR9IsoMap).end()) { pho.r9 = mapi->val;} 
00543       }
00544       // fill the H for H/E
00545       if (photonHoverEHIsoMap.isValid()) {
00546         mapi = (*photonHoverEHIsoMap).find(ref);  
00547         if (mapi !=(*photonHoverEHIsoMap).end()) { pho.hovereh = mapi->val;}
00548       }
00549       // fill the R9ID
00550       if (photonR9IDIsoMap.isValid()) {
00551         mapi = (*photonR9IDIsoMap).find(ref);
00552         if (mapi !=(*photonR9IDIsoMap).end()) { pho.r9ID = mapi->val;}
00553       }
00554 
00555       // store the photon into the vector
00556       theHLTPhotons.push_back(pho);
00557     }
00558   }
00559 }
00560 
00561 void HLTEgamma::MakeL1NonIsolatedPhotons(
00562     std::vector<OpenHLTPhoton> & theHLTPhotons,
00563     const edm::Handle<reco::RecoEcalCandidateCollection>   & recoNonIsolecalcands,
00564     const edm::Handle<reco::RecoEcalCandidateIsolationMap> & EcalNonIsolMap,
00565     const edm::Handle<reco::RecoEcalCandidateIsolationMap> & HcalNonIsolMap,
00566     const edm::Handle<reco::RecoEcalCandidateIsolationMap> & TrackNonIsolMap,
00567     const edm::Handle<reco::RecoEcalCandidateIsolationMap> & photonR9NonIsoMap,     
00568     const edm::Handle<reco::RecoEcalCandidateIsolationMap> & photonHoverEHNonIsoMap,      
00569     const edm::Handle<reco::RecoEcalCandidateIsolationMap> & photonR9IDNonIsoMap,
00570     EcalClusterLazyTools& lazyTools )
00571 {
00572   reco::RecoEcalCandidateIsolationMap::const_iterator mapi;
00573 
00574   if (recoNonIsolecalcands.isValid()) {
00575     for (reco::RecoEcalCandidateCollection::const_iterator recoecalcand = recoNonIsolecalcands->begin();
00576          recoecalcand!= recoNonIsolecalcands->end(); recoecalcand++) {
00577       // loop over SuperCluster and fill the HLTPhotons
00578       OpenHLTPhoton pho;
00579       pho.ecalIsol   = -999;
00580       pho.hcalIsol   = -999;
00581       pho.trackIsol  = -999;
00582       pho.clusterShape = -999;
00583       pho.L1Isolated = false;
00584       pho.Et         = recoecalcand->et();
00585       pho.eta        = recoecalcand->eta();
00586       pho.phi        = recoecalcand->phi();
00587       pho.r9         = -999; 
00588       pho.hovereh    = -999.;
00589       pho.r9ID       = -999.;
00590 
00591       //Get the cluster shape
00592       //      std::vector<float> vCov = lazyTools.covariances( *(recoecalcand->superCluster()->seed()) );
00593       std::vector<float> vCov = lazyTools.localCovariances( *(recoecalcand->superCluster()->seed()) );
00594       double sigmaee = sqrt(vCov[0]);
00595       //      float EtaSC = fabs(recoecalcand->eta());
00596       //      if(EtaSC > 1.479 ) {//Endcap
00597       //        sigmaee = sigmaee - 0.02*(EtaSC - 2.3);
00598       //      }
00599       pho.clusterShape = sigmaee;
00600 
00601       reco::RecoEcalCandidateRef ref = reco::RecoEcalCandidateRef(recoNonIsolecalcands, distance(recoNonIsolecalcands->begin(), recoecalcand));
00602 
00603       // fill the ecal Isolation
00604       if (EcalNonIsolMap.isValid()) {
00605         mapi = (*EcalNonIsolMap).find(ref);
00606         if (mapi !=(*EcalNonIsolMap).end()) { pho.ecalIsol = mapi->val;}
00607       }
00608       // fill the hcal Isolation
00609       if (HcalNonIsolMap.isValid()) {
00610         mapi = (*HcalNonIsolMap).find(ref);
00611         if (mapi !=(*HcalNonIsolMap).end()) { pho.hcalIsol = mapi->val;}
00612       }
00613       // fill the track Isolation
00614       if (TrackNonIsolMap.isValid()) {
00615         mapi = (*TrackNonIsolMap).find(ref);
00616         if (mapi !=(*TrackNonIsolMap).end()) { pho.trackIsol = mapi->val;}
00617       }
00618       // fill the R9 
00619       if (photonR9NonIsoMap.isValid()) { 
00620         mapi = (*photonR9NonIsoMap).find(ref);  
00621         if (mapi !=(*photonR9NonIsoMap).end()) { pho.r9 = mapi->val;}  
00622       } 
00623       // fill the H for H/E 
00624       if (photonHoverEHNonIsoMap.isValid()) { 
00625         mapi = (*photonHoverEHNonIsoMap).find(ref);   
00626         if (mapi !=(*photonHoverEHNonIsoMap).end()) { pho.hovereh = mapi->val;} 
00627       } 
00628       // fill the R9ID
00629       if (photonR9IDNonIsoMap.isValid()) {
00630         mapi = (*photonR9IDNonIsoMap).find(ref);
00631         if (mapi !=(*photonR9IDNonIsoMap).end()) { pho.r9ID = mapi->val;}
00632       }
00633 
00634       // store the photon into the vector
00635       theHLTPhotons.push_back(pho);
00636     }
00637   }
00638 }
00639 
00640 void HLTEgamma::MakeL1IsolatedElectrons(
00641     std::vector<OpenHLTElectron> & theHLTElectrons,
00642     const edm::Handle<reco::ElectronCollection>            & electronIsoHandle,
00643     const edm::Handle<reco::RecoEcalCandidateCollection>   & recoIsolecalcands,
00644     const edm::Handle<reco::RecoEcalCandidateIsolationMap> & HcalEleIsolMap,
00645     const edm::Handle<reco::ElectronSeedCollection>        & L1IsoPixelSeedsMap,
00646     const edm::Handle<reco::ElectronIsolationMap>          & TrackEleIsolMap,
00647     const edm::Handle<reco::RecoEcalCandidateIsolationMap> & electronR9IsoMap,     
00648     const edm::Handle<reco::RecoEcalCandidateIsolationMap> & photonHoverEHIsoMap,      
00649     const edm::Handle<reco::RecoEcalCandidateIsolationMap> & EcalIsolMap, 
00650     const edm::Handle<reco::RecoEcalCandidateIsolationMap> & electronR9IDIsoMap,
00651     EcalClusterLazyTools& lazyTools,
00652     const edm::ESHandle<MagneticField>& theMagField,
00653     reco::BeamSpot::Point & BSPosition )
00654 {
00655   // if there are electrons, then the isolation maps and the SC should be in the event; if not it is an error
00656   if (recoIsolecalcands.isValid()) {
00657     for (reco::RecoEcalCandidateCollection::const_iterator recoecalcand = recoIsolecalcands->begin();
00658          recoecalcand!= recoIsolecalcands->end(); recoecalcand++) {
00659       // get the ref to the SC:
00660       reco::RecoEcalCandidateRef ref = reco::RecoEcalCandidateRef(recoIsolecalcands, distance(recoIsolecalcands->begin(), recoecalcand));
00661       reco::SuperClusterRef recrSC = ref->superCluster();
00662       //reco::SuperClusterRef recrSC = recoecalcand->superCluster();
00663 
00664       OpenHLTElectron ele;
00665       ele.hcalIsol   = -999;
00666       ele.trackIsol  = -999;
00667       ele.ecalIsol   = -999;
00668       ele.L1Isolated = true;
00669       ele.p          = -999;
00670       ele.pixelSeeds = -999;
00671       ele.newSC      = true;
00672       ele.clusterShape = -999;
00673       ele.Dphi = 700; 
00674       ele.Deta = 700;
00675       ele.hovereh = -999;
00676       ele.Et         = recoecalcand->et();
00677       ele.eta        = recoecalcand->eta();
00678       ele.phi        = recoecalcand->phi();
00679       ele.E          = recrSC->energy();
00680       //Get the cluster shape
00681       //      std::vector<float> vCov = lazyTools.covariances( *(recrSC->seed()) );
00682       std::vector<float> vCov = lazyTools.localCovariances( *(recrSC->seed()) );
00683       double sigmaee = sqrt(vCov[0]);
00684       //      float EtaSC = fabs(recoecalcand->eta());
00685       //      if(EtaSC > 1.479 ) {//Endcap
00686       //        sigmaee = sigmaee - 0.02*(EtaSC - 2.3);
00687       //      }
00688       ele.clusterShape = sigmaee;
00689       ele.r9 = -999.;
00690       ele.r9ID = -999.;
00691 
00692       // fill the ecal Isolation 
00693       if (EcalIsolMap.isValid()) { 
00694         reco::RecoEcalCandidateIsolationMap::const_iterator mapi = (*EcalIsolMap).find(ref); 
00695         if (mapi !=(*EcalIsolMap).end()) { ele.ecalIsol = mapi->val;} 
00696       } 
00697       // fill the hcal Isolation
00698       if (HcalEleIsolMap.isValid()) {
00699         //reco::RecoEcalCandidateIsolationMap::const_iterator mapi = (*HcalEleIsolMap).find( reco::RecoEcalCandidateRef(recoIsolecalcands, distance(recoIsolecalcands->begin(), recoecalcand)) );
00700         reco::RecoEcalCandidateIsolationMap::const_iterator mapi = (*HcalEleIsolMap).find( ref );
00701         if (mapi !=(*HcalEleIsolMap).end()) { ele.hcalIsol = mapi->val; }
00702       }
00703       // fill the R9   
00704       if (electronR9IsoMap.isValid()) {   
00705         reco::RecoEcalCandidateIsolationMap::const_iterator mapi = (*electronR9IsoMap).find( ref );   
00706         if (mapi !=(*electronR9IsoMap).end()) { ele.r9 = mapi->val; }   
00707       }   
00708       // fill the H for H/E 
00709       if (photonHoverEHIsoMap.isValid()) { 
00710         reco::RecoEcalCandidateIsolationMap::const_iterator mapi = (*photonHoverEHIsoMap).find(ref);   
00711         if (mapi !=(*photonHoverEHIsoMap).end()) { ele.hovereh = mapi->val;} 
00712       } 
00713       // fill the R9ID
00714       if (electronR9IDIsoMap.isValid()) {
00715         reco::RecoEcalCandidateIsolationMap::const_iterator mapi = (*electronR9IDIsoMap).find( ref );
00716         if (mapi !=(*electronR9IDIsoMap).end()) { ele.r9ID = mapi->val; }
00717       }
00718 
00719       // look if the SC has associated pixelSeeds
00720       int nmatch = 0;
00721 
00722       if (L1IsoPixelSeedsMap.isValid()) {
00723         for (reco::ElectronSeedCollection::const_iterator it = L1IsoPixelSeedsMap->begin();
00724              it != L1IsoPixelSeedsMap->end(); it++) {
00725           edm::RefToBase<reco::CaloCluster> caloCluster = it->caloCluster() ;
00726           reco::SuperClusterRef scRef = caloCluster.castTo<reco::SuperClusterRef>() ;
00727           if (&(*recrSC) ==  &(*scRef)) { nmatch++; }
00728         }
00729       }
00730 
00731       ele.pixelSeeds = nmatch;
00732 
00733       // look if the SC was promoted to an electron:
00734       if (electronIsoHandle.isValid()) {
00735         bool FirstElectron = true;
00736         reco::ElectronRef electronref;
00737         for (reco::ElectronCollection::const_iterator iElectron = electronIsoHandle->begin();
00738              iElectron != electronIsoHandle->end(); iElectron++) {
00739           // 1) find the SC from the electron
00740           electronref = reco::ElectronRef(electronIsoHandle, iElectron - electronIsoHandle->begin());
00741           const reco::SuperClusterRef theClus = electronref->superCluster(); // SC from the electron;
00742           if (&(*recrSC) ==  &(*theClus)) {     // ref is the RecoEcalCandidateRef corresponding to the electron
00743             if (FirstElectron) {                // the first electron is stored in ele, keeping the ele.newSC = true
00744               FirstElectron = false;
00745               ele.p = electronref->track()->momentum().R();
00746               float deta=-100, dphi=-100;
00747               CalculateDetaDphi(theMagField,BSPosition , electronref , deta, dphi, false);
00748               ele.Dphi=dphi; ele.Deta=deta;
00749               // fill the track Isolation
00750               if (TrackEleIsolMap.isValid()) {
00751                 reco::ElectronIsolationMap::const_iterator mapTr = (*TrackEleIsolMap).find(electronref);
00752                 if (mapTr != (*TrackEleIsolMap).end()) { ele.trackIsol = mapTr->val; }
00753               }
00754             }
00755             else {
00756               // FirstElectron is false, i.e. the SC of this electron is common to another electron.
00757               // A new  OpenHLTElectron is inserted in the theHLTElectrons vector setting newSC = false
00758               OpenHLTElectron ele2;
00759               ele2.hcalIsol  = ele.hcalIsol;
00760               ele2.trackIsol = -999;
00761               ele2.Dphi = 700; 
00762               ele2.Deta = 700;
00763               ele2.Et  = ele.Et;
00764               ele2.eta = ele.eta;
00765               ele2.phi = ele.phi;
00766               ele2.E   = ele.E;
00767               ele2.L1Isolated = ele.L1Isolated;
00768               ele2.pixelSeeds = ele.pixelSeeds;
00769               ele2.clusterShape = ele.clusterShape;
00770               ele2.newSC = false;
00771               ele2.p = electronref->track()->momentum().R();
00772               ele2.r9 = ele.r9;
00773               ele2.hovereh = ele.hovereh;
00774               ele2.ecalIsol = ele.ecalIsol;
00775               ele2.r9ID = ele.r9ID;
00776               float deta=-100, dphi=-100;
00777               CalculateDetaDphi(theMagField,BSPosition , electronref , deta, dphi, false);
00778               ele2.Dphi=dphi; ele2.Deta=deta;
00779               // fill the track Isolation
00780               if (TrackEleIsolMap.isValid()) {
00781                 reco::ElectronIsolationMap::const_iterator mapTr = (*TrackEleIsolMap).find( electronref);
00782                 if (mapTr !=(*TrackEleIsolMap).end()) { ele2.trackIsol = mapTr->val;}
00783               }
00784               theHLTElectrons.push_back(ele2);
00785             }
00786           }
00787         } // end of loop over electrons
00788       } // end of if (electronIsoHandle) {
00789 
00790       //store the electron into the vector
00791       theHLTElectrons.push_back(ele);
00792     } // end of loop over ecalCandidates
00793   } // end of if (recoIsolecalcands) {
00794 }
00795 
00796 
00797 void HLTEgamma::MakeL1NonIsolatedElectrons(
00798     std::vector<OpenHLTElectron> & theHLTElectrons,
00799     const edm::Handle<reco::ElectronCollection>            & electronNonIsoHandle,
00800     const edm::Handle<reco::RecoEcalCandidateCollection>   & recoNonIsolecalcands,
00801     const edm::Handle<reco::RecoEcalCandidateIsolationMap> & HcalEleIsolMap,
00802     const edm::Handle<reco::ElectronSeedCollection>   & L1NonIsoPixelSeedsMap,
00803     const edm::Handle<reco::ElectronIsolationMap>          & TrackEleIsolMap,
00804     const edm::Handle<reco::RecoEcalCandidateIsolationMap> & electronR9NonIsoMap,     
00805     const edm::Handle<reco::RecoEcalCandidateIsolationMap> & photonHoverEHNonIsoMap,      
00806     const edm::Handle<reco::RecoEcalCandidateIsolationMap> & EcalNonIsolMap, 
00807     const edm::Handle<reco::RecoEcalCandidateIsolationMap> & electronR9IDNonIsoMap,
00808     EcalClusterLazyTools& lazyTools,
00809     const edm::ESHandle<MagneticField>& theMagField,
00810     reco::BeamSpot::Point & BSPosition  )
00811 {
00812   // if there are electrons, then the isolation maps and the SC should be in the event; if not it is an error
00813   if (recoNonIsolecalcands.isValid()) {
00814     for (reco::RecoEcalCandidateCollection::const_iterator recoecalcand = recoNonIsolecalcands->begin();
00815          recoecalcand!= recoNonIsolecalcands->end(); recoecalcand++) {
00816       //get the ref to the SC:
00817       reco::RecoEcalCandidateRef ref = reco::RecoEcalCandidateRef(recoNonIsolecalcands, distance(recoNonIsolecalcands->begin(), recoecalcand));
00818       reco::SuperClusterRef recrSC = ref->superCluster();
00819       //reco::SuperClusterRef recrSC = recoecalcand->superCluster();
00820 
00821       OpenHLTElectron ele;
00822       ele.hcalIsol   = -999;
00823       ele.trackIsol  = -999;
00824       ele.ecalIsol   = -999; 
00825       ele.L1Isolated = false;
00826       ele.p          = -999;
00827       ele.pixelSeeds = -999;
00828       ele.newSC      = true;
00829       ele.clusterShape = -999;
00830       ele.Dphi = 700; 
00831       ele.Deta = 700;
00832       ele.r9 = -999.;
00833       ele.r9ID = -999.;
00834       ele.hovereh = -999; 
00835       ele.Et         = recoecalcand->et();
00836       ele.eta        = recoecalcand->eta();
00837       ele.phi        = recoecalcand->phi();
00838       ele.E          = recrSC->energy();
00839       //Get the cluster shape
00840       //      std::vector<float> vCov = lazyTools.covariances( *(recrSC->seed()) );
00841       std::vector<float> vCov = lazyTools.localCovariances( *(recrSC->seed()) );
00842       double sigmaee = sqrt(vCov[0]);
00843       //      float EtaSC = fabs(recoecalcand->eta());
00844       //      if(EtaSC > 1.479 ) {//Endcap
00845       //        sigmaee = sigmaee - 0.02*(EtaSC - 2.3);
00846       //      }
00847       ele.clusterShape = sigmaee;
00848 
00849       // fill the ecal Isolation 
00850       if (EcalNonIsolMap.isValid()) { 
00851         reco::RecoEcalCandidateIsolationMap::const_iterator mapi = (*EcalNonIsolMap).find(ref); 
00852         if (mapi !=(*EcalNonIsolMap).end()) { ele.ecalIsol = mapi->val;} 
00853       } 
00854       // fill the hcal Isolation
00855       if (HcalEleIsolMap.isValid()) {
00856         // reco::RecoEcalCandidateIsolationMap::const_iterator mapi = (*HcalEleIsolMap).find( reco::RecoEcalCandidateRef(recoNonIsolecalcands, distance(recoNonIsolecalcands->begin(), recoecalcand)) );
00857         reco::RecoEcalCandidateIsolationMap::const_iterator mapi = (*HcalEleIsolMap).find( ref );
00858         if (mapi !=(*HcalEleIsolMap).end()) {ele.hcalIsol = mapi->val;}
00859       }
00860       // fill the R9    
00861       if (electronR9NonIsoMap.isValid()) {    
00862         reco::RecoEcalCandidateIsolationMap::const_iterator mapi = (*electronR9NonIsoMap).find( ref );    
00863         if (mapi !=(*electronR9NonIsoMap).end()) { ele.r9 = mapi->val; }    
00864       }    
00865       // fill the H for H/E 
00866       if (photonHoverEHNonIsoMap.isValid()) { 
00867         reco::RecoEcalCandidateIsolationMap::const_iterator mapi = (*photonHoverEHNonIsoMap).find(ref);   
00868         if (mapi !=(*photonHoverEHNonIsoMap).end()) { ele.hovereh = mapi->val;} 
00869       } 
00870       // fill the R9ID
00871       if (electronR9IDNonIsoMap.isValid()) {
00872         reco::RecoEcalCandidateIsolationMap::const_iterator mapi = (*electronR9IDNonIsoMap).find( ref );
00873         if (mapi !=(*electronR9IDNonIsoMap).end()) { ele.r9ID = mapi->val; }
00874       }
00875 
00876       // look if the SC has associated pixelSeeds
00877       int nmatch = 0;
00878 
00879       if (L1NonIsoPixelSeedsMap.isValid()) {
00880         for (reco::ElectronSeedCollection::const_iterator it = L1NonIsoPixelSeedsMap->begin();
00881              it != L1NonIsoPixelSeedsMap->end(); it++) {
00882           edm::RefToBase<reco::CaloCluster> caloCluster = it->caloCluster() ;
00883           reco::SuperClusterRef scRef = caloCluster.castTo<reco::SuperClusterRef>() ;
00884           if (&(*recrSC) == &(*scRef)) { nmatch++;}
00885         }
00886       }
00887 
00888       ele.pixelSeeds = nmatch;
00889 
00890       // look if the SC was promoted to an electron:
00891       if (electronNonIsoHandle.isValid()) {
00892         bool FirstElectron = true;
00893         reco::ElectronRef electronref;
00894         for (reco::ElectronCollection::const_iterator iElectron = electronNonIsoHandle->begin(); 
00895              iElectron != electronNonIsoHandle->end();iElectron++) {
00896           // 1) find the SC from the electron
00897           electronref = reco::ElectronRef(electronNonIsoHandle, iElectron - electronNonIsoHandle->begin());
00898           const reco::SuperClusterRef theClus = electronref->superCluster(); //SC from the electron;
00899           if (&(*recrSC) ==  &(*theClus)) { // ref is the RecoEcalCandidateRef corresponding to the electron
00900             if (FirstElectron) { //the first electron is stored in ele, keeping the ele.newSC = true
00901               FirstElectron = false;
00902               ele.p = electronref->track()->momentum().R();
00903               float deta=-100, dphi=-100;
00904               CalculateDetaDphi(theMagField,BSPosition , electronref , deta, dphi, false);
00905               ele.Dphi=dphi; ele.Deta=deta;
00906 
00907               // fill the track Isolation
00908               if (TrackEleIsolMap.isValid()) {
00909                 reco::ElectronIsolationMap::const_iterator mapTr = (*TrackEleIsolMap).find( electronref);
00910                 if (mapTr !=(*TrackEleIsolMap).end()) { ele.trackIsol = mapTr->val;}
00911               }
00912             } else {
00913               // FirstElectron is false, i.e. the SC of this electron is common to another electron.
00914               // A new OpenHLTElectron is inserted in the theHLTElectrons vector setting newSC = false
00915               OpenHLTElectron ele2;
00916               ele2.hcalIsol   = ele.hcalIsol;
00917               ele2.trackIsol  =-999;
00918               ele2.ecalIsol = ele.ecalIsol;
00919               ele2.Dphi = 700; 
00920               ele2.Deta = 700;
00921               ele2.Et         = ele.Et;
00922               ele2.eta        = ele.eta;
00923               ele2.phi        = ele.phi;
00924               ele2.E          = ele.E;
00925               ele2.L1Isolated = ele.L1Isolated;
00926               ele2.pixelSeeds = ele.pixelSeeds;
00927               ele2.clusterShape = ele.clusterShape;
00928               ele2.newSC      = false;
00929               ele2.p          = electronref->track()->momentum().R();
00930               ele2.r9         = ele.r9;
00931               ele2.hovereh = ele.hovereh; 
00932               ele2.r9ID       = ele.r9ID;
00933               float deta=-100, dphi=-100;
00934               CalculateDetaDphi(theMagField,BSPosition , electronref , deta, dphi, false);
00935               ele2.Dphi=dphi; ele2.Deta=deta;
00936 
00937               // fill the track Isolation
00938               if (TrackEleIsolMap.isValid()) {
00939                 reco::ElectronIsolationMap::const_iterator mapTr = (*TrackEleIsolMap).find( electronref);
00940                 if (mapTr !=(*TrackEleIsolMap).end()) { ele2.trackIsol = mapTr->val;}
00941               }
00942               theHLTElectrons.push_back(ele2);
00943             }
00944           }
00945         } // end of loop over electrons
00946       } // end of if (electronNonIsoHandle) {
00947 
00948       // store the electron into the vector
00949       theHLTElectrons.push_back(ele);
00950     } // end of loop over ecalCandidates
00951   } // end of if (recoNonIsolecalcands) {
00952 }
00953 
00954 void HLTEgamma::CalculateDetaDphi(const edm::ESHandle<MagneticField>& theMagField, 
00955                                   reco::BeamSpot::Point & BSPosition, 
00956                                   const reco::ElectronRef eleref, 
00957                                   float& deltaeta, 
00958                                   float& deltaphi, bool useTrackProjectionToEcal )
00959 {
00960   
00961   const reco::SuperClusterRef theClus = eleref->superCluster();
00962   math::XYZVector scv(theClus->x(), theClus->y(), theClus->z());
00963   
00964   const math::XYZVector trackMom =  eleref->track()->momentum();
00965 
00966     math::XYZPoint SCcorrPosition(theClus->x()-BSPosition.x(), theClus->y()-BSPosition.y() , theClus->z()-eleref->track()->vz() );
00967     deltaeta = SCcorrPosition.eta()-eleref->track()->eta();
00968     
00969     if(useTrackProjectionToEcal){
00970       ECALPositionCalculator posCalc;
00971       const math::XYZPoint vertex(BSPosition.x(),BSPosition.y(),eleref->track()->vz());
00972       
00973       float phi1= posCalc.ecalPhi(theMagField.product(),trackMom,vertex,1);
00974       float phi2= posCalc.ecalPhi(theMagField.product(),trackMom,vertex,-1);
00975       
00976       float deltaphi1=fabs( phi1 - theClus->position().phi() );
00977       if(deltaphi1>6.283185308) deltaphi1 -= 6.283185308;
00978       if(deltaphi1>3.141592654) deltaphi1 = 6.283185308-deltaphi1;
00979       
00980       float deltaphi2=fabs( phi2 - theClus->position().phi() );
00981       if(deltaphi2>6.283185308) deltaphi2 -= 6.283185308;
00982       if(deltaphi2>3.141592654) deltaphi2 = 6.283185308-deltaphi2;
00983       
00984       deltaphi = deltaphi1;
00985       if(deltaphi2<deltaphi1){ deltaphi = deltaphi2;}
00986     }
00987     else {
00988       deltaphi=fabs(eleref->track()->outerPosition().phi()-theClus->phi());
00989       if(deltaphi>6.283185308) deltaphi -= 6.283185308;
00990       if(deltaphi>3.141592654) deltaphi = 6.283185308-deltaphi;
00991     }
00992 
00993 }