CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_1_8_patch9/src/DQMOffline/Trigger/src/EgHLTOffHelper.cc

Go to the documentation of this file.
00001 #include "DQMOffline/Trigger/interface/EgHLTOffHelper.h"
00002 
00003 
00004 #include "DataFormats/EcalDetId/interface/EcalSubdetector.h"
00005 #include "DataFormats/DetId/interface/DetId.h"
00006 #include "DataFormats/JetReco/interface/CaloJetCollection.h"
00007 #include "FWCore/Common/interface/TriggerNames.h"
00008 
00009 #include "RecoEcal/EgammaCoreTools/interface/EcalClusterTools.h"
00010 
00011 #include "DataFormats/TrackReco/interface/Track.h"
00012 #include "DataFormats/GsfTrackReco/interface/GsfTrack.h"
00013 
00014 #include "RecoEgamma/EgammaHLTAlgos/interface/EgammaHLTTrackIsolation.h"
00015 #include "RecoEgamma/EgammaIsolationAlgos/interface/EgammaTowerIsolation.h"
00016 #include "RecoEgamma/EgammaIsolationAlgos/interface/EgammaRecHitIsolation.h"
00017 
00018 #include "DQMOffline/Trigger/interface/EgHLTTrigCodes.h"
00019 #include "DQMOffline/Trigger/interface/EgHLTTrigTools.h"
00020 #include "DQMOffline/Trigger/interface/EgHLTErrCodes.h"
00021 using namespace egHLT;
00022 
00023 OffHelper::~OffHelper()
00024 {
00025   if(hltEleTrkIsolAlgo_) delete hltEleTrkIsolAlgo_;
00026   if(hltPhoTrkIsolAlgo_) delete hltPhoTrkIsolAlgo_;
00027 }
00028 
00029 void OffHelper::setup(const edm::ParameterSet& conf)
00030 {
00031 
00032   ecalRecHitsEBTag_ = conf.getParameter<edm::InputTag>("BarrelRecHitCollection");
00033   ecalRecHitsEETag_ = conf.getParameter<edm::InputTag>("EndcapRecHitCollection");
00034   caloJetsTag_ = conf.getParameter<edm::InputTag>("CaloJetCollection");
00035   isolTrkTag_ = conf.getParameter<edm::InputTag>("IsolTrackCollection");
00036   hbheHitsTag_ = conf.getParameter<edm::InputTag>("HBHERecHitCollection");
00037   hfHitsTag_ = conf.getParameter<edm::InputTag>("HFRecHitCollection");
00038   electronsTag_ = conf.getParameter<edm::InputTag>("ElectronCollection");
00039   photonsTag_ = conf.getParameter<edm::InputTag>("PhotonCollection");
00040   triggerSummaryLabel_ = conf.getParameter<edm::InputTag>("triggerSummaryLabel");
00041   hltTag_ = conf.getParameter<std::string>("hltTag");
00042   beamSpotTag_ = conf.getParameter<edm::InputTag>("BeamSpotProducer");
00043   caloTowersTag_ = conf.getParameter<edm::InputTag>("CaloTowers");
00044   trigResultsTag_ = conf.getParameter<edm::InputTag>("TrigResults");
00045 
00046   eleCuts_.setup(conf.getParameter<edm::ParameterSet>("eleCuts"));
00047   eleLooseCuts_.setup(conf.getParameter<edm::ParameterSet>("eleLooseCuts"));
00048   phoCuts_.setup(conf.getParameter<edm::ParameterSet>("phoCuts"));
00049   phoLooseCuts_.setup(conf.getParameter<edm::ParameterSet>("phoLooseCuts"));
00050  
00051   //now we have the isolations completely configurable via python
00052   hltEMIsolOuterCone_ = conf.getParameter<double>("hltEMIsolOuterCone");
00053   hltEMIsolInnerConeEB_ = conf.getParameter<double>("hltEMIsolInnerConeEB");
00054   hltEMIsolEtaSliceEB_ = conf.getParameter<double>("hltEMIsolEtaSliceEB");
00055   hltEMIsolEtMinEB_ = conf.getParameter<double>("hltEMIsolEtMinEB");
00056   hltEMIsolEMinEB_ = conf.getParameter<double>("hltEMIsolEMinEB");
00057   hltEMIsolInnerConeEE_ = conf.getParameter<double>("hltEMIsolInnerConeEE");
00058   hltEMIsolEtaSliceEE_ = conf.getParameter<double>("hltEMIsolEtaSliceEE");
00059   hltEMIsolEtMinEE_ = conf.getParameter<double>("hltEMIsolEtMinEE");
00060   hltEMIsolEMinEE_ = conf.getParameter<double>("hltEMIsolEMinEE");
00061 
00062   hltPhoTrkIsolPtMin_ = conf.getParameter<double>("hltPhoTrkIsolPtMin");
00063   hltPhoTrkIsolOuterCone_ = conf.getParameter<double>("hltPhoTrkIsolOuterCone");
00064   hltPhoTrkIsolInnerCone_ = conf.getParameter<double>("hltPhoTrkIsolInnerCone"); 
00065   hltPhoTrkIsolZSpan_ = conf.getParameter<double>("hltPhoTrkIsolZSpan");
00066   hltPhoTrkIsolRSpan_ = conf.getParameter<double>("hltPhoTrkIsolZSpan");
00067   hltPhoTrkIsolCountTrks_ = conf.getParameter<bool>("hltPhoTrkIsolCountTrks");
00068 
00069   hltEleTrkIsolPtMin_ = conf.getParameter<double>("hltEleTrkIsolPtMin");
00070   hltEleTrkIsolOuterCone_ = conf.getParameter<double>("hltEleTrkIsolOuterCone");
00071   hltEleTrkIsolInnerCone_ = conf.getParameter<double>("hltEleTrkIsolInnerCone"); 
00072   hltEleTrkIsolZSpan_ = conf.getParameter<double>("hltEleTrkIsolZSpan");
00073   hltEleTrkIsolRSpan_ = conf.getParameter<double>("hltEleTrkIsolZSpan");
00074 
00075   hltHadIsolOuterCone_ = conf.getParameter<double>("hltHadIsolOuterCone");
00076   hltHadIsolInnerCone_ = conf.getParameter<double>("hltHadIsolInnerCone");
00077   hltHadIsolEtMin_ = conf.getParameter<double>("hltHadIsolEtMin");
00078   hltHadIsolDepth_ = conf.getParameter<int>("hltHadIsolDepth"); 
00079 
00080   calHLTHcalIsol_ = conf.getParameter<bool>("calHLTHcalIsol");
00081   calHLTEmIsol_ = conf.getParameter<bool>("calHLTEmIsol");
00082   calHLTEleTrkIsol_ = conf.getParameter<bool>("calHLTEleTrkIsol");
00083   calHLTPhoTrkIsol_ = conf.getParameter<bool>("calHLTPhoTrkIsol");
00084  
00085   trigCutParams_ = conf.getParameter<std::vector<edm::ParameterSet> >("triggerCuts"); //setupTriggers used to be in this function but had to be moved due to HLTConfigChanges (has to be called beginRun) so we have to save this for later.
00086 
00087   hltEleTrkIsolAlgo_ = new EgammaHLTTrackIsolation(hltEleTrkIsolPtMin_,hltEleTrkIsolOuterCone_,hltEleTrkIsolZSpan_,hltEleTrkIsolRSpan_,hltEleTrkIsolInnerCone_);
00088   hltPhoTrkIsolAlgo_ = new EgammaHLTTrackIsolation(hltPhoTrkIsolPtMin_,hltPhoTrkIsolOuterCone_,hltPhoTrkIsolZSpan_,hltPhoTrkIsolRSpan_,hltPhoTrkIsolInnerCone_);
00089                                                
00090 
00091 }
00092 
00093 //this code was taken out of OffHelper::setup due to HLTConfigProvider changes
00094 //it still assumes that this is called only once
00095 void OffHelper::setupTriggers(const HLTConfigProvider& hltConfig,const std::vector<std::string>& hltFiltersUsed)
00096 {
00097   hltFiltersUsed_ = hltFiltersUsed; //expensive but only do this once and faster ways could make things less clear
00098   //now work out how many objects are requires to pass filter for it to accept
00099   hltFiltersUsedWithNrCandsCut_.clear();
00100   for(size_t filterNr=0;filterNr<hltFiltersUsed_.size();filterNr++){
00101     hltFiltersUsedWithNrCandsCut_.push_back(std::make_pair(hltFiltersUsed_[filterNr],egHLT::trigTools::getMinNrObjsRequiredByFilter(hltFiltersUsed_[filterNr])));
00102   }
00103 
00104   //now loading the cuts for every trigger into our vector which stores them
00105   //only load cuts for triggers that are in hltFiltersUsed
00106   
00107   for(size_t trigNr=0;trigNr<trigCutParams_.size();trigNr++) {
00108     std::string trigName = trigCutParams_[trigNr].getParameter<std::string>("trigName");
00109     if(std::find(hltFiltersUsed_.begin(),hltFiltersUsed_.end(),trigName)!=hltFiltersUsed_.end()){ //perhaps I should sort hltFiltersUsed_....
00110       trigCuts_.push_back(std::make_pair(TrigCodes::getCode(trigName),OffEgSel(trigCutParams_[trigNr])));
00111     }
00112   }
00113   trigCutParams_.clear();//dont need it any more, get rid of it
00114 
00115   //to make my life difficult, the scaled l1 paths are special
00116   //and arent stored in trigger event
00117   //to I have to figure out the path, see if it passes
00118   //and then hunt down the l1 seed filter and use that to match to the pho/ele
00119   //matching on l1 seed filter is not enough as that will be passed for normal 
00120   //electron triggers even if pre-scale hasnt fired
00121   l1PreScaledFilters_.clear();
00122   l1PreScaledPaths_.clear();
00123   l1PreAndSeedFilters_.clear();
00124   for(size_t filterNr=0;filterNr<hltFiltersUsed_.size();filterNr++){
00125     if(hltFiltersUsed_[filterNr].find("hltPreL1")==0){ //l1 prescaled path
00126       l1PreScaledFilters_.push_back(hltFiltersUsed_[filterNr]);
00127     }
00128   }
00129 
00130   egHLT::trigTools::translateFiltersToPathNames(hltConfig,l1PreScaledFilters_,l1PreScaledPaths_);
00131   if(l1PreScaledPaths_.size()==l1PreScaledFilters_.size()){
00132     for(size_t pathNr=0;pathNr<l1PreScaledPaths_.size();pathNr++){
00133      
00134       std::string l1SeedFilter =egHLT::trigTools::getL1SeedFilterOfPath(hltConfig,l1PreScaledPaths_[pathNr]);
00135      
00136       l1PreAndSeedFilters_.push_back(std::make_pair(l1PreScaledFilters_[pathNr],l1SeedFilter));
00137     }
00138   }
00139 }
00140 
00141 int OffHelper::makeOffEvt(const edm::Event& edmEvent,const edm::EventSetup& setup,egHLT::OffEvt& offEvent)
00142 {
00143   offEvent.clear();
00144   int errCode=0; //excution stops as soon as an error is flagged
00145   if(errCode==0) errCode = getHandles(edmEvent,setup);
00146   if(errCode==0) errCode = fillOffEleVec(offEvent.eles());
00147   if(errCode==0) errCode = fillOffPhoVec(offEvent.phos());
00148   if(errCode==0) errCode =  setTrigInfo(edmEvent, offEvent);
00149   if(errCode==0) offEvent.setJets(recoJets_);
00150   return errCode;
00151 }
00152 
00153 int OffHelper::getHandles(const edm::Event& event,const edm::EventSetup& setup)
00154 {
00155   try { 
00156     setup.get<CaloGeometryRecord>().get(caloGeom_);
00157     setup.get<CaloTopologyRecord>().get(caloTopology_);
00158   }catch(...){
00159     return errCodes::Geom;
00160   }
00161   try {
00162     setup.get<IdealMagneticFieldRecord>().get(magField_);
00163   }catch(...){
00164     return errCodes::MagField;
00165   }
00166 
00167   //get objects
00168   if(!getHandle(event,triggerSummaryLabel_,trigEvt_)) return errCodes::TrigEvent; //must have this, otherwise skip event
00169   if(!getHandle(event,trigResultsTag_,trigResults_)) return errCodes::TrigEvent; //re using bit to minimise bug fix code changes
00170   if(!getHandle(event,electronsTag_,recoEles_)) return errCodes::OffEle; //need for electrons
00171   if(!getHandle(event,photonsTag_, recoPhos_)) return errCodes::OffPho; //need for photons
00172   if(!getHandle(event,caloJetsTag_,recoJets_)) return errCodes::OffJet; //need for electrons and photons
00173 
00174   //need for HLT isolations (rec hits also need for sigmaIPhiIPhi (ele/pho) and r9 pho)
00175   if(!getHandle(event,ecalRecHitsEBTag_,ebRecHits_)) return errCodes::EBRecHits;
00176   if(!getHandle(event,ecalRecHitsEETag_,eeRecHits_)) return errCodes::EERecHits;
00177   if(!getHandle(event,isolTrkTag_,isolTrks_)) return errCodes::IsolTrks;
00178   if(!getHandle(event,hbheHitsTag_, hbheHits_)) return errCodes::HBHERecHits; //I dont think we need hbhe rec-hits any more
00179   if(!getHandle(event,hfHitsTag_, hfHits_)) return errCodes::HFRecHits;//I dont think we need hf rec-hits any more
00180   if(!getHandle(event,beamSpotTag_,beamSpot_)) return errCodes::BeamSpot;
00181   if(!getHandle(event,caloTowersTag_,caloTowers_)) return errCodes::CaloTowers;
00182 
00183   
00184   return 0;
00185 }
00186 
00187 //this function coverts GsfElectrons to a format which is actually useful to me
00188 int OffHelper::fillOffEleVec(std::vector<OffEle>& egHLTOffEles)
00189 {
00190   egHLTOffEles.clear();
00191   egHLTOffEles.reserve(recoEles_->size());
00192   for(reco::GsfElectronCollection::const_iterator gsfIter=recoEles_->begin(); gsfIter!=recoEles_->end();++gsfIter){
00193     if(!gsfIter->ecalDrivenSeed()) continue; //avoid PF electrons (this is Eg HLT validation and HLT is ecal driven)
00194 
00195     OffEle::IsolData isolData;   
00196     fillIsolData(*gsfIter,isolData);
00197     
00198     OffEle::ClusShapeData clusShapeData;
00199     fillClusShapeData(*gsfIter,clusShapeData);
00200 
00201     OffEle::HLTData hltData;
00202     fillHLTData(*gsfIter,hltData);
00203 
00204     egHLTOffEles.push_back(OffEle(*gsfIter,clusShapeData,isolData,hltData));
00205     
00206     //now we would like to set the cut results
00207     OffEle& ele =  egHLTOffEles.back();
00208     ele.setCutCode(eleCuts_.getCutCode(ele));
00209     ele.setLooseCutCode(eleLooseCuts_.getCutCode(ele));
00210     
00211     std::vector<std::pair<TrigCodes::TrigBitSet,int> >trigCutsCutCodes;
00212     for(size_t i=0;i<trigCuts_.size();i++) trigCutsCutCodes.push_back(std::make_pair(trigCuts_[i].first,trigCuts_[i].second.getCutCode(ele)));
00213    
00214     ele.setTrigCutsCutCodes(trigCutsCutCodes);
00215   }//end loop over gsf electron collection
00216   return 0;
00217 }
00218 
00219 void OffHelper::fillIsolData(const reco::GsfElectron& ele,OffEle::IsolData& isolData)
00220 {
00221   EgammaTowerIsolation hcalIsolAlgo(hltHadIsolOuterCone_,hltHadIsolInnerCone_,hltHadIsolEtMin_,hltHadIsolDepth_,caloTowers_.product());
00222   EcalRecHitMetaCollection ebHits(*ebRecHits_);
00223   EcalRecHitMetaCollection eeHits(*eeRecHits_);
00224   EgammaRecHitIsolation ecalIsolAlgoEB(hltEMIsolOuterCone_,hltEMIsolInnerConeEB_,hltEMIsolEtaSliceEB_,
00225                                        hltEMIsolEtMinEB_,hltEMIsolEMinEB_,caloGeom_,&ebHits,DetId::Ecal);
00226   EgammaRecHitIsolation ecalIsolAlgoEE(hltEMIsolOuterCone_,hltEMIsolInnerConeEE_,hltEMIsolEtaSliceEE_,
00227                                        hltEMIsolEtMinEE_,hltEMIsolEMinEE_,caloGeom_,&eeHits,DetId::Ecal);
00228   
00229   isolData.ptTrks=ele.dr03TkSumPt();
00230   isolData.nrTrks=999; //no longer supported
00231   isolData.em= ele.dr03EcalRecHitSumEt();
00232   isolData.hadDepth1 =  ele.dr03HcalDepth1TowerSumEt();
00233   isolData.hadDepth2 =  ele.dr03HcalDepth2TowerSumEt();   
00234 
00235   //now time to do the HLT algos
00236   if(calHLTHcalIsol_) isolData.hltHad=hcalIsolAlgo.getTowerESum(&ele);
00237   else isolData.hltHad = 0.;
00238   if(calHLTEleTrkIsol_) isolData.hltTrksEle=hltEleTrkIsolAlgo_->electronPtSum(&(*(ele.gsfTrack())),isolTrks_.product());
00239   else isolData.hltTrksEle = 0.;
00240   if(calHLTPhoTrkIsol_){
00241     if(hltPhoTrkIsolCountTrks_) isolData.hltTrksPho=hltPhoTrkIsolAlgo_->photonTrackCount(&ele,isolTrks_.product(),false);
00242     else isolData.hltTrksPho=hltPhoTrkIsolAlgo_->photonPtSum(&ele,isolTrks_.product(),false);
00243   }
00244   else isolData.hltTrksPho = 0.;
00245   if(calHLTEmIsol_) isolData.hltEm = ecalIsolAlgoEB.getEtSum(&ele) + ecalIsolAlgoEE.getEtSum(&ele);
00246   else isolData.hltEm = 0.;
00247   
00248 }
00249 
00250 
00251 void OffHelper::fillClusShapeData(const reco::GsfElectron& ele,OffEle::ClusShapeData& clusShapeData)
00252 {
00253   clusShapeData.sigmaEtaEta = ele.sigmaEtaEta();
00254   clusShapeData.sigmaIEtaIEta = ele.sigmaIetaIeta();
00255   double e5x5 = ele.e5x5();
00256   if(e5x5!=0.){
00257     clusShapeData.e1x5Over5x5 = ele.e1x5()/e5x5;
00258     clusShapeData.e2x5MaxOver5x5 = ele.e2x5Max()/e5x5;
00259   }else{
00260     clusShapeData.e1x5Over5x5 = -1;
00261     clusShapeData.e2x5MaxOver5x5 = -1;
00262   }
00263   
00264   //want to calculate r9, sigmaPhiPhi and sigmaIPhiIPhi, have to do old fashioned way
00265   const reco::BasicCluster& seedClus = *(ele.superCluster()->seed());
00266   const DetId seedDetId = seedClus.hitsAndFractions()[0].first; //note this may not actually be the seed hit but it doesnt matter because all hits will be in the barrel OR endcap
00267   if(seedDetId.subdetId()==EcalBarrel){
00268     std::vector<float> stdCov = EcalClusterTools::covariances(seedClus,ebRecHits_.product(),caloTopology_.product(),caloGeom_.product());
00269     std::vector<float> crysCov = EcalClusterTools::localCovariances(seedClus,ebRecHits_.product(),caloTopology_.product());
00270     clusShapeData.sigmaPhiPhi = sqrt(stdCov[2]);
00271       clusShapeData.sigmaIPhiIPhi =  sqrt(crysCov[2]);
00272       if(ele.superCluster()->rawEnergy()!=0.){
00273         clusShapeData.r9 = EcalClusterTools::e3x3(seedClus,ebRecHits_.product(),caloTopology_.product()) / ele.superCluster()->rawEnergy();
00274       }else clusShapeData.r9 = -1.;
00275        
00276   }else{
00277     std::vector<float> stdCov = EcalClusterTools::covariances(seedClus,eeRecHits_.product(),caloTopology_.product(),caloGeom_.product()); 
00278     std::vector<float> crysCov = EcalClusterTools::localCovariances(seedClus,eeRecHits_.product(),caloTopology_.product());
00279     clusShapeData.sigmaPhiPhi = sqrt(stdCov[2]);
00280     clusShapeData.sigmaIPhiIPhi = sqrt(crysCov[2]);
00281     if(ele.superCluster()->rawEnergy()!=0.){
00282         clusShapeData.r9 = EcalClusterTools::e3x3(seedClus,eeRecHits_.product(),caloTopology_.product()) / ele.superCluster()->rawEnergy();
00283     }else clusShapeData.r9 = -1.;
00284   } 
00285 }
00286 
00287 //reco approximations of hlt quantities
00288 void OffHelper::fillHLTData(const reco::GsfElectron& ele,OffEle::HLTData& hltData)
00289 {
00290   if(ele.closestCtfTrackRef().isNonnull() && 
00291      ele.closestCtfTrackRef()->extra().isNonnull()){
00292     reco::TrackRef ctfTrack = ele.closestCtfTrackRef();
00293     reco::SuperClusterRef scClus = ele.superCluster();
00294 
00295     //dEta
00296     const reco::BeamSpot::Point& bsPos = beamSpot_->position();     
00297     math::XYZPoint scPosWRTVtx(scClus->x()-bsPos.x(), scClus->y()-bsPos.y() , scClus->z()-ctfTrack->vz());
00298     hltData.dEtaIn = fabs(scPosWRTVtx.eta()-ctfTrack->eta());
00299 
00300     //dPhi: lifted straight from hlt code
00301     float deltaPhi=fabs(ctfTrack->outerPosition().phi()-scClus->phi());
00302     if(deltaPhi>6.283185308) deltaPhi -= 6.283185308;
00303     if(deltaPhi>3.141592654) deltaPhi = 6.283185308-deltaPhi;
00304     hltData.dPhiIn = deltaPhi;
00305     
00306     //invEInvP
00307     if(ele.ecalEnergy()!=0 && ctfTrack->p()!=0) hltData.invEInvP= 1/ele.ecalEnergy() - 1/ctfTrack->p();
00308     else hltData.invEInvP = 0;
00309   }else{
00310     hltData.dEtaIn =999;
00311     hltData.dPhiIn =999;
00312     hltData.invEInvP = 999;
00313 
00314   }
00315 }
00316 
00317 //this function coverts Photons to a format which more useful to me
00318 int OffHelper::fillOffPhoVec(std::vector<OffPho>& egHLTOffPhos)
00319 {
00320   egHLTOffPhos.clear();
00321   egHLTOffPhos.reserve(recoPhos_->size());
00322   for(reco::PhotonCollection::const_iterator phoIter=recoPhos_->begin(); phoIter!=recoPhos_->end();++phoIter){
00323 
00324     OffPho::IsolData isolData;  
00325     OffPho::ClusShapeData clusShapeData;
00326   
00327     fillIsolData(*phoIter,isolData);
00328     fillClusShapeData(*phoIter,clusShapeData);
00329  
00330     egHLTOffPhos.push_back(OffPho(*phoIter,clusShapeData,isolData));
00331     OffPho& pho =  egHLTOffPhos.back();
00332     pho.setCutCode(phoCuts_.getCutCode(pho));
00333     pho.setLooseCutCode(phoLooseCuts_.getCutCode(pho));
00334 
00335     std::vector<std::pair<TrigCodes::TrigBitSet,int> >trigCutsCutCodes;
00336     for(size_t i=0;i<trigCuts_.size();i++) trigCutsCutCodes.push_back(std::make_pair(trigCuts_[i].first,trigCuts_[i].second.getCutCode(pho)));
00337     pho.setTrigCutsCutCodes(trigCutsCutCodes); 
00338 
00339 
00340   }//end loop over photon collection
00341   return 0;
00342 }
00343 
00344 
00345 void OffHelper::fillIsolData(const reco::Photon& pho,OffPho::IsolData& isolData)
00346 {
00347   EgammaTowerIsolation hcalIsolAlgo(hltHadIsolOuterCone_,hltHadIsolInnerCone_,hltHadIsolEtMin_,hltHadIsolDepth_,caloTowers_.product());
00348   EcalRecHitMetaCollection ebHits(*ebRecHits_);
00349   EcalRecHitMetaCollection eeHits(*ebRecHits_);
00350   EgammaRecHitIsolation ecalIsolAlgoEB(hltEMIsolOuterCone_,hltEMIsolInnerConeEB_,hltEMIsolEtaSliceEB_,
00351                                        hltEMIsolEtMinEB_,hltEMIsolEMinEB_,caloGeom_,&ebHits,DetId::Ecal);
00352   EgammaRecHitIsolation ecalIsolAlgoEE(hltEMIsolOuterCone_,hltEMIsolInnerConeEE_,hltEMIsolEtaSliceEE_,
00353                                        hltEMIsolEtMinEE_,hltEMIsolEMinEE_,caloGeom_,&eeHits,DetId::Ecal);
00354   
00355   isolData.nrTrks = pho.nTrkHollowConeDR03();
00356   isolData.ptTrks = pho.trkSumPtHollowConeDR03();
00357   isolData.em = pho.ecalRecHitSumEtConeDR03();
00358   isolData.had = pho.hcalTowerSumEtConeDR03();  
00359 
00360   //now calculate hlt algos
00361   if(calHLTHcalIsol_) isolData.hltHad=hcalIsolAlgo.getTowerESum(&pho);
00362   else isolData.hltHad = 0.;
00363   if(calHLTPhoTrkIsol_){
00364     if(hltPhoTrkIsolCountTrks_) isolData.hltTrks=hltPhoTrkIsolAlgo_->photonTrackCount(&pho,isolTrks_.product(),false);
00365     else isolData.hltTrks=hltPhoTrkIsolAlgo_->photonPtSum(&pho,isolTrks_.product(),false);
00366   }
00367   else isolData.hltTrks = 0.;
00368   if(calHLTEmIsol_) isolData.hltEm = ecalIsolAlgoEB.getEtSum(&pho) + ecalIsolAlgoEE.getEtSum(&pho);
00369   else isolData.hltEm = 0.;
00370   
00371 }
00372 
00373 void OffHelper::fillClusShapeData(const reco::Photon& pho,OffPho::ClusShapeData& clusShapeData)
00374 {
00375   clusShapeData.sigmaEtaEta = pho.sigmaEtaEta();
00376   clusShapeData.sigmaIEtaIEta =  pho.sigmaIetaIeta();
00377   double e5x5 =  pho.e5x5();
00378   if(e5x5!=0.){ //even though it is almost impossible for this to be 0., this code can never ever crash under any situation
00379     clusShapeData.e1x5Over5x5 =  pho.e1x5()/e5x5;
00380     clusShapeData.e2x5MaxOver5x5 = pho.e2x5()/e5x5;
00381   }else{
00382     clusShapeData.e1x5Over5x5 = -1;
00383     clusShapeData.e2x5MaxOver5x5 = -1;
00384   }
00385   clusShapeData.r9 = pho.r9();
00386   
00387   //sigmaPhiPhi and sigmaIPhiIPhi are not in object (and nor should they be) so have to get them old fashioned way
00388   //need to figure out if its in the barrel or endcap
00389   //get the first hit of the cluster and figure out if its barrel or endcap 
00390   const reco::BasicCluster& seedClus = *(pho.superCluster()->seed());
00391   const DetId seedDetId = seedClus.hitsAndFractions()[0].first; //note this may not actually be the seed hit but it doesnt matter because all hits will be in the barrel OR endcap (it is also incredably inefficient as it getHitsByDetId passes the vector by value not reference
00392   if(seedDetId.subdetId()==EcalBarrel){
00393     std::vector<float> stdCov = EcalClusterTools::covariances(seedClus,ebRecHits_.product(),caloTopology_.product(),caloGeom_.product());
00394     std::vector<float> crysCov = EcalClusterTools::localCovariances(seedClus,ebRecHits_.product(),caloTopology_.product());
00395     clusShapeData.sigmaPhiPhi = sqrt(stdCov[2]);
00396     clusShapeData.sigmaIPhiIPhi =  sqrt(crysCov[2]);
00397     }else{
00398     std::vector<float> stdCov = EcalClusterTools::covariances(seedClus,eeRecHits_.product(),caloTopology_.product(),caloGeom_.product()); 
00399     std::vector<float> crysCov = EcalClusterTools::localCovariances(seedClus,eeRecHits_.product(),caloTopology_.product());
00400     
00401     clusShapeData.sigmaPhiPhi = sqrt(stdCov[2]);
00402     clusShapeData.sigmaIPhiIPhi = sqrt(crysCov[2]); 
00403   }
00404 }  
00405 
00406 int OffHelper::setTrigInfo(const edm::Event & edmEvent, egHLT::OffEvt& offEvent)
00407 {
00408   TrigCodes::TrigBitSet evtTrigBits = trigTools::getFiltersPassed(hltFiltersUsedWithNrCandsCut_,trigEvt_.product(),hltTag_);
00409   //the l1 prescale paths dont have a filter with I can figure out if it passed or failed with so have to use TriggerResults
00410   if(l1PreScaledPaths_.size()==l1PreScaledFilters_.size()){ //check to ensure both vectors have same number of events incase of screw ups     
00411     const edm::TriggerNames & triggerNames = edmEvent.triggerNames(*trigResults_);
00412     for(size_t pathNr=0;pathNr<l1PreScaledPaths_.size();pathNr++){ //now we have to check the prescaled l1 trigger paths
00413       unsigned int pathIndex = triggerNames.triggerIndex(l1PreScaledPaths_[pathNr]);
00414       if(pathIndex<trigResults_->size() && trigResults_->accept(pathIndex)){
00415         evtTrigBits |=TrigCodes::getCode(l1PreScaledFilters_[pathNr]);
00416       }    
00417     }
00418   }
00419 
00420   offEvent.setEvtTrigBits(evtTrigBits);
00421 
00422   trigTools::setFiltersObjPasses(offEvent.eles(),hltFiltersUsed_,l1PreAndSeedFilters_,evtTrigBits,trigEvt_.product(),hltTag_);
00423   trigTools::setFiltersObjPasses(offEvent.phos(),hltFiltersUsed_,l1PreAndSeedFilters_,evtTrigBits,trigEvt_.product(),hltTag_); 
00424   return 0;
00425 }