CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_4_5_patch3/src/Calibration/HcalAlCaRecoProducers/src/AlCaIsoTracksProducer.cc

Go to the documentation of this file.
00001 #include "FWCore/Framework/interface/ESHandle.h"
00002 #include "FWCore/Framework/interface/EventSetup.h"
00003 #include "Calibration/HcalAlCaRecoProducers/interface/AlCaIsoTracksProducer.h"
00004 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00005 #include "Geometry/Records/interface/CaloGeometryRecord.h"
00006 #include "Geometry/CaloGeometry/interface/CaloGeometry.h"
00007 #include "Geometry/CaloGeometry/interface/CaloSubdetectorGeometry.h"
00008 #include "Geometry/CaloGeometry/interface/CaloCellGeometry.h" 
00009 #include "DataFormats/EcalDetId/interface/EcalSubdetector.h"
00010 #include "DataFormats/HcalDetId/interface/HcalSubdetector.h" 
00011 #include "DataFormats/GeometrySurface/interface/Cylinder.h"
00012 #include "DataFormats/GeometrySurface/interface/Plane.h" 
00013 #include "Geometry/CommonDetUnit/interface/GeomDetUnit.h"
00014 #include "MagneticField/Engine/interface/MagneticField.h" 
00015 #include "MagneticField/Records/interface/IdealMagneticFieldRecord.h"
00016 #include "TrackPropagation/SteppingHelixPropagator/interface/SteppingHelixPropagator.h"
00017 
00018 #include "DataFormats/HcalRecHit/interface/HcalRecHitCollections.h"
00019 #include "DataFormats/EcalRecHit/interface/EcalRecHitCollections.h"
00020 #include "FWCore/Utilities/interface/Exception.h"
00021 
00022 #include "Geometry/CaloTopology/interface/EcalPreshowerTopology.h"
00023 
00024 #include "DataFormats/TrackReco/interface/TrackExtra.h"
00025 #include "DataFormats/TrackReco/interface/TrackExtraFwd.h"
00026 #include "RecoTracker/TrackProducer/interface/TrackProducerBase.h"
00027 
00028 #include "DataFormats/HcalIsolatedTrack/interface/IsolatedPixelTrackCandidate.h"
00029 #include "DataFormats/HcalIsolatedTrack/interface/IsolatedPixelTrackCandidateFwd.h"
00030 
00031 #include "DataFormats/HLTReco/interface/TriggerEvent.h"
00032 
00033 #include "DataFormats/EcalDetId/interface/EBDetId.h"
00034 #include "DataFormats/EcalDetId/interface/EEDetId.h"
00035 
00036 #include "Math/GenVector/VectorUtil.h"
00037 #include "Math/GenVector/PxPyPzE4D.h"
00038 
00039 #include "DataFormats/L1Trigger/interface/L1JetParticleFwd.h"
00040 #include "DataFormats/L1Trigger/interface/L1JetParticle.h"
00041 
00042 #include "DataFormats/Provenance/interface/ProductID.h"
00043 
00044 #include "Calibration/HcalAlCaRecoProducers/plugins/ConeDefinition.h"
00045 
00046 #include <boost/regex.hpp> 
00047 
00048 double getDistInCM(double eta1, double phi1, double eta2, double phi2)
00049 {
00050   double dR, Rec;
00051   double theta1=2*atan(exp(-eta1));
00052   double theta2=2*atan(exp(-eta2));
00053   if (fabs(eta1)<1.479) Rec=129; //radius of ECAL barrel
00054   else Rec=tan(theta1)*317; //distance from IP to ECAL endcap
00055 
00056   //|vect| times tg of acos(scalar product)
00057   double angle=acos((sin(theta1)*sin(theta2)*(sin(phi1)*sin(phi2)+cos(phi1)*cos(phi2))+cos(theta1)*cos(theta2)));
00058   if (angle<acos(-1)/2)
00059     {
00060       dR=fabs((Rec/sin(theta1))*tan(angle));
00061       return dR;
00062     }
00063   else return 1000;
00064 }
00065 
00066 double getDist(double eta1, double phi1, double eta2, double phi2)
00067 {
00068   double dphi = fabs(phi1 - phi2); 
00069   if(dphi>acos(-1)) dphi = 2*acos(-1)-dphi;
00070   double dr = sqrt(dphi*dphi + std::pow(eta1-eta2,2));
00071   return dr;
00072 }
00073 
00074 bool checkHLTMatch(edm::Event& iEvent, edm::InputTag hltEventTag_, std::vector<std::string> hltFilterTag_, double eta, double phi, double hltMatchingCone_)
00075 {
00076   bool match =false;
00077   double minDDD=1000;
00078   
00079   edm::Handle<trigger::TriggerEvent> trEv;
00080   iEvent.getByLabel(hltEventTag_,trEv);
00081   const trigger::TriggerObjectCollection& TOCol(trEv->getObjects());
00082   
00083   trigger::Keys KEYS;
00084   const trigger::size_type nFilt(trEv->sizeFilters());
00085   for (trigger::size_type iFilt=0; iFilt!=nFilt; iFilt++) 
00086     {
00087       for (unsigned l=0; l<hltFilterTag_.size(); l++)
00088         {
00089           if ((trEv->filterTag(iFilt).label()).substr(0,27)==hltFilterTag_[l]) 
00090             {
00091               KEYS=trEv->filterKeys(iFilt);
00092             }
00093         }
00094     }
00095   trigger::size_type nReg=KEYS.size();
00096   for (trigger::size_type iReg=0; iReg<nReg; iReg++)
00097     {
00098       const trigger::TriggerObject& TObj(TOCol[KEYS[iReg]]);
00099       double dHit=getDist(TObj.eta(),TObj.phi(),eta,phi); 
00100       if (dHit<minDDD) minDDD=dHit;
00101     }
00102   if (minDDD>hltMatchingCone_) match=false;
00103   else match=true;
00104      
00105   return match;
00106 }
00107 
00108 std::pair<double,double> getL1triggerDirection(edm::Event& iEvent, edm::InputTag hltEventTag_, std::string l1FilterTag_)
00109 {
00110   edm::Handle<trigger::TriggerEvent> trEv;
00111   iEvent.getByLabel(hltEventTag_,trEv);
00112   const trigger::TriggerObjectCollection& TOCol(trEv->getObjects());
00113   
00114   trigger::Keys KEYS;
00115   const trigger::size_type nFilt(trEv->sizeFilters());
00116   for (trigger::size_type iFilt=0; iFilt!=nFilt; iFilt++)
00117     {
00118       if ((trEv->filterTag(iFilt).label()).substr(0,14)==l1FilterTag_) KEYS=trEv->filterKeys(iFilt); 
00119     }
00120   trigger::size_type nReg=KEYS.size();
00121   double etaTrig=-10000;
00122   double phiTrig=-10000;
00123   double ptMax=0;
00124   for (trigger::size_type iReg=0; iReg<nReg; iReg++)
00125     {
00126       const trigger::TriggerObject& TObj(TOCol[KEYS[iReg]]);
00127       if (TObj.pt()>ptMax)
00128         {
00129           etaTrig=TObj.eta();
00130           phiTrig=TObj.phi();
00131           ptMax=TObj.pt();
00132         }
00133     }
00134   return std::pair<double,double>(etaTrig,phiTrig);
00135 }
00136 
00137 
00138 AlCaIsoTracksProducer::AlCaIsoTracksProducer(const edm::ParameterSet& iConfig)
00139 { 
00140   
00141   m_inputTrackLabel_ = iConfig.getParameter<edm::InputTag>("InputTracksLabel");
00142 
00143   hoLabel_ = iConfig.getParameter<edm::InputTag>("HOInput");
00144 
00145   ecalLabels_=iConfig.getParameter<std::vector<edm::InputTag> >("ECALInputs");
00146 
00147   hbheLabel_= iConfig.getParameter<edm::InputTag>("HBHEInput");
00148   
00149   m_dvCut = iConfig.getParameter<double>("vtxCut");
00150   m_ddirCut = iConfig.getParameter<double>("RIsolAtHCALSurface");
00151   useConeCorr_=iConfig.getParameter<bool>("UseLowPtConeCorrection");
00152   m_pCut = iConfig.getParameter<double>("MinTrackP");
00153   m_ptCut = iConfig.getParameter<double>("MinTrackPt");
00154   m_ecalCut = iConfig.getUntrackedParameter<double>("NeutralIsolCut",8.);
00155 
00156   taECALCone_=iConfig.getUntrackedParameter<double>("TrackAssociatorECALCone",0.5);
00157   taHCALCone_=iConfig.getUntrackedParameter<double>("TrackAssociatorHCALCone",1.0);
00158 
00159   skipNeutrals_=iConfig.getUntrackedParameter<bool>("SkipNeutralIsoCheck",false);
00160 
00161   nHitsMinCore_=iConfig.getParameter<int>("MinNumberOfHitsCoreTrack");
00162   nHitsMinIso_=iConfig.getParameter<int>("MinNumberOfHitsInConeTracks");
00163 
00164   isolE_ = iConfig.getParameter<double>("MaxNearbyTrackEnergy");
00165   etaMax_= iConfig.getParameter<double>("MaxTrackEta");
00166   cluRad_ = iConfig.getParameter<double>("ECALClusterRadius");
00167   ringOutRad_ = iConfig.getParameter<double>("ECALRingOuterRadius");
00168   ringInnRad_=iConfig.getParameter<double>("ECALRingInnerRadius");  
00169 
00170   useECALCluMatrix_=iConfig.getParameter<bool>("ClusterECALasMatrix");
00171   matrixSize_=iConfig.getParameter<int>("ECALMatrixFullSize");
00172 
00173   checkHLTMatch_=iConfig.getParameter<bool>("CheckHLTMatch");
00174   hltEventTag_=iConfig.getParameter<edm::InputTag>("hltTriggerEventLabel");  
00175   hltFiltTag_=iConfig.getParameter<std::vector<std::string> >("hltL3FilterLabels");
00176   hltMatchingCone_=iConfig.getParameter<double>("hltMatchingCone");
00177   l1FilterTag_=iConfig.getParameter<std::string>("l1FilterLabel");
00178   l1jetVetoCone_=iConfig.getParameter<double>("l1JetVetoCone");
00179   
00181 //
00182 // Parameters for track associator   ===========================
00183 //  
00184   edm::ParameterSet parameters = iConfig.getParameter<edm::ParameterSet>("TrackAssociatorParameters");
00185   parameters_.loadParameters( parameters );
00186   trackAssociator_.useDefaultPropagator();
00187 // ===============================================================
00188 
00189   //create also IsolatedPixelTrackCandidateCollection which contains isolation info and reference to primary track
00190   produces<reco::IsolatedPixelTrackCandidateCollection>("HcalIsolatedTrackCollection");
00191 
00192   produces<reco::TrackCollection>("IsoTrackTracksCollection");
00193   produces<reco::TrackExtraCollection>("IsoTrackExtraTracksCollection");
00194 
00195   produces<EcalRecHitCollection>("IsoTrackEcalRecHitCollection");
00196   produces<EcalRecHitCollection>("IsoTrackPSEcalRecHitCollection");
00197   
00198   produces<HBHERecHitCollection>("IsoTrackHBHERecHitCollection");
00199   produces<HORecHitCollection>("IsoTrackHORecHitCollection");
00200 
00201 }
00202 
00203 
00204 AlCaIsoTracksProducer::~AlCaIsoTracksProducer() { }
00205 
00206 
00207 void AlCaIsoTracksProducer::produce(edm::Event& iEvent, const edm::EventSetup& iSetup)
00208 {
00209   std::auto_ptr<reco::IsolatedPixelTrackCandidateCollection> outputHcalIsoTrackColl(new reco::IsolatedPixelTrackCandidateCollection);
00210   std::auto_ptr<reco::TrackCollection> outputTColl(new reco::TrackCollection);
00211   std::auto_ptr<reco::TrackExtraCollection> outputExTColl(new reco::TrackExtraCollection);
00212   
00213   reco::TrackExtraRefProd rTrackExtras = iEvent.getRefBeforePut<reco::TrackExtraCollection>("IsoTrackExtraTracksCollection"); 
00214   std::auto_ptr<EcalRecHitCollection> outputEColl(new EcalRecHitCollection);
00215   std::auto_ptr<EcalRecHitCollection> outputESColl(new EcalRecHitCollection);
00216   std::auto_ptr<HBHERecHitCollection> outputHColl(new HBHERecHitCollection);
00217   std::auto_ptr<HORecHitCollection> outputHOColl(new HORecHitCollection);
00218   
00219   edm::ESHandle<CaloGeometry> pG;
00220   iSetup.get<CaloGeometryRecord>().get(pG);
00221   geo = pG.product();
00222   
00223   edm::Handle<reco::TrackCollection> trackCollection;
00224   iEvent.getByLabel(m_inputTrackLabel_,trackCollection);
00225   
00226   
00227   // temporary collection of EB+EE recHits
00228   std::auto_ptr<EcalRecHitCollection> tmpEcalRecHitCollection(new EcalRecHitCollection);
00229   
00230   std::vector<edm::InputTag>::const_iterator i;
00231   for (i=ecalLabels_.begin(); i!=ecalLabels_.end(); i++) 
00232     {
00233       edm::Handle<EcalRecHitCollection> ec;
00234       iEvent.getByLabel(*i,ec);
00235       for(EcalRecHitCollection::const_iterator recHit = (*ec).begin(); recHit != (*ec).end(); ++recHit)
00236         {
00237           tmpEcalRecHitCollection->push_back(*recHit);
00238         }
00239     }      
00240   
00241   edm::Handle<HBHERecHitCollection> hbheRHcol;
00242   iEvent.getByLabel(hbheLabel_, hbheRHcol);
00243 
00244   const reco::TrackCollection tC = *(trackCollection.product());
00245   
00246   int itrk=0;
00247   int nisotr=0;
00248   edm::Ref<reco::TrackExtraCollection>::key_type  idx = 0;
00249   
00250   //   Parameters for TrackDetAssociator ================================
00251   // For Low momentum tracks need to look for larger cone for search ====
00252   // ====================================================================
00253   
00254   parameters_.useEcal = true;
00255   parameters_.useHcal = true;
00256   parameters_.useCalo = false;
00257   parameters_.useMuon = false;
00258   parameters_.dREcal = taECALCone_;
00259   parameters_.dRHcal = taHCALCone_;
00260   
00262   
00264   std::vector<int> usedHitsHC;
00265   std::vector<int> usedHitsEC;
00267 
00268   // main loop over input tracks
00269   for (reco::TrackCollection::const_iterator track=tC.begin(); track!=tC.end(); track++) {
00270     bool noChargedTracks = true;
00271     int itrk1=0;
00272     itrk++;
00273     double px = track->px();
00274     double py = track->py();
00275     double pz = track->pz();
00276     double ptrack = sqrt(px*px+py*py+pz*pz);
00277     
00278     if (ptrack < m_pCut || track->pt() < m_ptCut ) continue; 
00279  
00280     if (track->hitPattern().numberOfValidHits() < nHitsMinCore_) continue;
00281 
00282     // check that track is not in the region of L1 jet
00283     double l1jDR=deltaR(track->eta(), track->phi(), getL1triggerDirection(iEvent,hltEventTag_,l1FilterTag_).first,getL1triggerDirection(iEvent,hltEventTag_,l1FilterTag_).second);
00284     if (l1jDR<l1jetVetoCone_) continue;
00286             
00287     TrackDetMatchInfo info = trackAssociator_.associate(iEvent, iSetup, trackAssociator_.getFreeTrajectoryState(iSetup, *track), parameters_);
00288     
00289     GlobalPoint gPointHcal(info.trkGlobPosAtHcal.x(),info.trkGlobPosAtHcal.y(),info.trkGlobPosAtHcal.z());
00290 
00291     GlobalVector trackMomAtHcal = info.trkMomAtHcal;
00292     
00293     double etaecal=info.trkGlobPosAtEcal.eta();
00294     double phiecal=info.trkGlobPosAtEcal.phi();
00295 
00296     if (etaecal==0&&phiecal==0) continue;    
00297 
00298     //check matching to HLT object (optional)
00299 
00300     if (checkHLTMatch_&&!checkHLTMatch(iEvent, hltEventTag_, hltFiltTag_, etaecal, phiecal,hltMatchingCone_)) continue;
00301     
00302     if (fabs(track->eta())>etaMax_) continue;
00303     
00304     // check charged isolation from all other tracks
00305     double maxPNearby=-10;
00306     double sumPNearby=0;
00307     for (reco::TrackCollection::const_iterator track1=tC.begin(); track1!=tC.end(); track1++)
00308       {
00309         itrk1++;
00310         if (track == track1) continue;
00311         if (track->hitPattern().numberOfValidHits() < nHitsMinIso_) continue;
00312         double ptrack1 = sqrt(track1->px()*track1->px()+track1->py()*track1->py()+track1->pz()*track1->pz());
00313         
00314         TrackDetMatchInfo info1 = trackAssociator_.associate(iEvent, iSetup, trackAssociator_.getFreeTrajectoryState(iSetup, *track1), parameters_);
00315         
00316         GlobalPoint gPointHcal1(info1.trkGlobPosAtHcal.x(),info1.trkGlobPosAtHcal.y(),info1.trkGlobPosAtHcal.z());
00317 
00318         double etaecal1=info1.trkGlobPosAtEcal.eta();
00319         double phiecal1=info1.trkGlobPosAtEcal.phi();
00320    
00321         if (etaecal1==0&&phiecal1==0) continue; 
00322 
00323         double hcDist=getDistInPlaneTrackDir(gPointHcal, trackMomAtHcal, gPointHcal1);
00324 
00325 //        double hcDist=getDistInCM(etaecal,phiecal,etaecal1,phiecal1);
00326 
00327         // increase required separation for low momentum tracks
00328         double factor=1.;
00329         double factor1=1.;
00330 
00331         if (useConeCorr_)
00332            {
00333              if(ptrack<10.)factor+=(10.-ptrack)/20.;
00334              if(ptrack1<10.)factor1+=(10.-ptrack1)/20.;
00335            }
00336         
00337         if( hcDist <  m_ddirCut*factor*factor1 ) 
00338           {
00339             //calculate maximum P and sum P near seed track
00340             if (track1->p()>maxPNearby)
00341               {
00342                 maxPNearby=track1->p();
00343               }
00344             sumPNearby+=track1->p();
00345             
00346             //apply loose isolation criteria
00347             if (track1->p()>isolE_) 
00348               {
00349                 noChargedTracks = false;
00350                 break;
00351               }
00353           }
00354         
00355       } //second track loop
00356 
00357     bool noNeutrals = false;
00358     
00359     // we have a good charge-isolated track, so check neutral isolation and write it out
00360         
00361     if(noChargedTracks) 
00362       {
00363         //find ecal cluster energy and write ecal recHits
00364         double ecClustR=0;
00365         double ecClustN=0;
00366         double ecOutRingR=0;
00367 
00368         //get index of ECAL crystal hit by track
00369         std::vector<const EcalRecHit*> crossedECids=info.crossedEcalRecHits;
00370         int etaIDcenter=-10000;
00371         int phiIDcenter=-10000;
00372         double enMax=0;
00373         for (unsigned int i=0; i<crossedECids.size(); i++)
00374           {
00375             if ((*crossedECids[i]).id().subdetId()==EcalEndcap)
00376               {
00377                 EEDetId did(crossedECids[i]->id());
00378                 if (crossedECids[i]->energy()>enMax)
00379                   {
00380                     enMax=crossedECids[i]->energy();
00381                     etaIDcenter=did.iy();
00382                     phiIDcenter=did.ix();
00383                   }
00384               }
00385             if ((*crossedECids[i]).id().subdetId()==EcalBarrel)
00386               {
00387                 EBDetId did(crossedECids[i]->id());
00388                 if (crossedECids[i]->energy()>enMax)
00389                   {
00390                     enMax=crossedECids[i]->energy();
00391                     etaIDcenter=did.ieta();
00392                     phiIDcenter=did.iphi();
00393                   }
00394               }
00395           }
00396         for (std::vector<EcalRecHit>::const_iterator ehit=tmpEcalRecHitCollection->begin(); ehit!=tmpEcalRecHitCollection->end(); ehit++) 
00397           {
00399             // R scheme of ECAL CLUSTERIZATION
00400             GlobalPoint posH = geo->getPosition((*ehit).detid());
00401             double phihit = posH.phi();
00402             double etahit = posH.eta();
00403             
00404             double dHit=deltaR(etaecal,phiecal,etahit,phihit);
00405             
00406             double dHitCM=getDistInPlaneTrackDir(gPointHcal, trackMomAtHcal, posH);
00407             
00408             if (dHitCM<cluRad_)
00409               {
00410                 ecClustR+=ehit->energy();
00411               }
00412             
00413             if (dHitCM>ringInnRad_&&dHitCM<ringOutRad_)
00414               {
00415                 ecOutRingR+=ehit->energy();
00416               }
00417 
00419             //NxN scheme & check whether hit was used or not, if not used push into usedHits
00420             bool hitIsUsed=false;
00421             int hitHashedIndex=-10000;
00422             if (ehit->id().subdetId()==EcalBarrel)
00423               {
00424                 EBDetId did(ehit->id());
00425                 hitHashedIndex=did.hashedIndex();
00426                 if (fabs(did.ieta()-etaIDcenter)<=matrixSize_/2&&fabs(did.iphi()-phiIDcenter)<=matrixSize_/2) ecClustN+=ehit->energy();
00427               }
00428             
00429             if (ehit->id().subdetId()==EcalEndcap)
00430               {
00431                 EEDetId did(ehit->id());
00432                 hitHashedIndex=did.hashedIndex();
00433                 if (fabs(did.iy()-etaIDcenter)<=matrixSize_/2&&fabs(did.ix()-phiIDcenter)<=matrixSize_/2) ecClustN+=ehit->energy();
00434               }
00435             for (uint32_t i=0; i<usedHitsEC.size(); i++)
00436               {
00437                 if (usedHitsEC[i]==hitHashedIndex) hitIsUsed=true;
00438               }
00439 
00440             if (hitIsUsed) continue; //skip if used
00441             usedHitsEC.push_back(hitHashedIndex);
00443             
00444             if(dHit<1.)  
00445               {
00446                 outputEColl->push_back(*ehit);
00447               }   
00448           }
00449 
00450         //check neutrals 
00451         if (ecOutRingR<m_ecalCut) noNeutrals=true;
00452         else noNeutrals=false;
00453 
00454         if (noNeutrals||skipNeutrals_)
00455         {
00456 
00457         // Take info on the track extras and keep it in the outercollection
00458 
00459         reco::TrackExtraRef myextra = (*track).extra();
00460         reco::TrackExtraRef teref= reco::TrackExtraRef ( rTrackExtras, idx ++ );
00461 
00462         outputTColl->push_back(*track);
00463         reco::Track & mytrack = outputTColl->back();
00464 
00465         mytrack.setExtra( teref );
00466         outputExTColl->push_back(*myextra);
00467         //      reco::TrackExtra & tx = outputExTColl->back();
00468 
00469         //Create IsolatedPixelTrackCandidate (will change naming in future release)
00470         reco::IsolatedPixelTrackCandidate newHITCandidate(reco::Candidate::LorentzVector(track->px(),track->py(),track->pz(),track->p()));
00471         newHITCandidate.SetSumPtPxl(sumPNearby);
00472         newHITCandidate.SetMaxPtPxl(maxPNearby);
00473 
00474         //set cluster energy deposition and ring energy deposition and push_back
00475         if (!useECALCluMatrix_) newHITCandidate.SetEnergyIn(ecClustR);
00476         else newHITCandidate.SetEnergyIn(ecClustN);
00477         newHITCandidate.SetEnergyOut(ecOutRingR);
00478         outputHcalIsoTrackColl->push_back(newHITCandidate);
00479 
00480         //save hcal recHits
00481         for (std::vector<HBHERecHit>::const_iterator hhit=hbheRHcol->begin(); hhit!=hbheRHcol->end(); hhit++) 
00482           {
00483             //check that this hit was not considered before and push it into usedHits
00484             bool hitIsUsed=false;
00485             int hitHashedIndex=hhit->id().hashed_index();
00486             for (uint32_t i=0; i<usedHitsHC.size(); i++)
00487               {
00488                 if (usedHitsHC[i]==hitHashedIndex) hitIsUsed=true;
00489               }
00490             if (hitIsUsed) continue;
00491             usedHitsHC.push_back(hitHashedIndex);
00493             
00494             GlobalPoint posH = geo->getPosition((*hhit).detid());
00495             double phihit = posH.phi();
00496             double etahit = posH.eta();
00497             
00498             double dHit=deltaR(etaecal,phiecal,etahit,phihit);
00499             
00500             if(dHit<1.)  
00501 
00502               {
00503                 outputHColl->push_back(*hhit);
00504               }
00505             }
00506           }
00507         
00508         nisotr++;
00509         
00510       } //if (noNeutrals....
00511   
00512   } // end of main track loop
00513   
00514   if(outputTColl->size() > 0)
00515     {
00516       //   Take HO collection
00517       edm::Handle<HORecHitCollection> ho;
00518       iEvent.getByLabel(hoLabel_,ho);
00519       
00520       const HORecHitCollection Hitho = *(ho.product());
00521       for(HORecHitCollection::const_iterator hoItr=Hitho.begin(); hoItr!=Hitho.end(); hoItr++)
00522         {
00523           outputHOColl->push_back(*hoItr);
00524         }
00525            
00526       // Take Preshower collection      
00527       
00528       // get the ps geometry
00529       edm::ESHandle<CaloGeometry> geoHandle;
00530       iSetup.get<CaloGeometryRecord>().get(geoHandle);
00531       
00532       // get the ps topology
00533       EcalPreshowerTopology psTopology(geoHandle);
00534       
00535       // process rechits
00536       edm::Handle<EcalRecHitCollection> pRecHits;
00537       iEvent.getByLabel("ecalPreshowerRecHit","EcalRecHitsES",pRecHits);
00538       const EcalRecHitCollection& psrechits = *(pRecHits.product());
00539 
00540       typedef EcalRecHitCollection::const_iterator IT;
00541       
00542       for(IT i=psrechits.begin(); i!=psrechits.end(); i++) 
00543         {
00544           outputESColl->push_back( *i );
00545         }
00546     }  
00547 
00548   iEvent.put( outputHcalIsoTrackColl, "HcalIsolatedTrackCollection");
00549   iEvent.put( outputTColl, "IsoTrackTracksCollection");
00550   iEvent.put( outputExTColl, "IsoTrackExtraTracksCollection");
00551   iEvent.put( outputEColl, "IsoTrackEcalRecHitCollection");
00552   iEvent.put( outputESColl, "IsoTrackPSEcalRecHitCollection");
00553   iEvent.put( outputHColl, "IsoTrackHBHERecHitCollection");
00554   iEvent.put( outputHOColl, "IsoTrackHORecHitCollection");
00555 }
00556 
00557 void AlCaIsoTracksProducer::endJob(void) {}
00558