CMS 3D CMS Logo

/data/doxygen/doxygen-1.7.3/gen/CMSSW_4_2_8/src/GeneratorInterface/GenFilters/src/doubleEMEnrichingFilterAlgo.cc

Go to the documentation of this file.
00001 #include "GeneratorInterface/GenFilters/interface/doubleEMEnrichingFilterAlgo.h"
00002 
00003 #include "FWCore/Framework/interface/ESHandle.h"
00004 
00005 #include "DataFormats/Math/interface/deltaR.h"
00006 #include "DataFormats/Math/interface/deltaPhi.h"
00007 #include "DataFormats/GeometrySurface/interface/Cylinder.h"
00008 #include "DataFormats/GeometrySurface/interface/Plane.h"
00009 #include "DataFormats/GeometryVector/interface/GlobalPoint.h"
00010 #include "DataFormats/GeometryVector/interface/GlobalVector.h"
00011 
00012 #include "MagneticField/Engine/interface/MagneticField.h"
00013 #include "MagneticField/Records/interface/IdealMagneticFieldRecord.h"
00014 
00015 #include "TrackingTools/GeomPropagators/interface/AnalyticalPropagator.h"
00016 #include "TrackingTools/TrajectoryState/interface/FreeTrajectoryState.h"
00017 #include "TrackingTools/TrajectoryState/interface/TrajectoryStateOnSurface.h"
00018 
00019 #include "CLHEP/Vector/LorentzVector.h"
00020 
00021 
00022 using namespace edm;
00023 using namespace std;
00024 
00025 
00026 doubleEMEnrichingFilterAlgo::doubleEMEnrichingFilterAlgo(const edm::ParameterSet& iConfig) { 
00027 
00028   //set constants
00029   FILTER_TKISOCUT_=2;
00030   FILTER_CALOISOCUT_=2;
00031   FILTER_ETA_MIN_=0;
00032   FILTER_ETA_MAX_=2.5;
00033   ECALBARRELMAXETA_=1.479;
00034   ECALBARRELRADIUS_=129.0;
00035   ECALENDCAPZ_=304.5;
00036 
00037   // from bctoe
00038   eTThreshold_=(float)iConfig.getParameter<double>("eTThreshold");
00039   
00040   isoGenParETMin_=(float) iConfig.getParameter<double>("isoGenParETMin");
00041   isoGenParConeSize_=(float) iConfig.getParameter<double>("isoGenParConeSize");
00042   clusterThreshold_=(float) iConfig.getParameter<double>("clusterThreshold");
00043   seedThreshold_=(float) iConfig.getParameter<double>("seedThreshold");
00044   isoConeSize_=(float) iConfig.getParameter<double>("isoConeSize");
00045   hOverEMax_=(float) iConfig.getParameter<double>("hOverEMax");
00046   tkIsoMax_=(float) iConfig.getParameter<double>("tkIsoMax");
00047   caloIsoMax_=(float) iConfig.getParameter<double>("caloIsoMax");
00048   requireTrackMatch_=iConfig.getParameter<bool>("requireTrackMatch");
00049   genParSource_=iConfig.getParameter<edm::InputTag>("genParSource");
00050 
00051 }
00052 
00053 doubleEMEnrichingFilterAlgo::~doubleEMEnrichingFilterAlgo() {
00054 }
00055 
00056 
00057 bool doubleEMEnrichingFilterAlgo::filter(const edm::Event& iEvent, const edm::EventSetup& iSetup)  {
00058 
00059   
00060   //----------
00061 
00062   Handle<reco::GenParticleCollection> genParsHandle;
00063   iEvent.getByLabel(genParSource_,genParsHandle);
00064   reco::GenParticleCollection genPars=*genParsHandle;
00065 
00066   bool result = false;
00067   int BCtoEgood=0;
00068   int PhotElecgood=0;
00069   int IsoGenPargood=0;
00070 
00071   // cleaning
00072   sel1seeds.clear();
00073   sel2seeds.clear();
00074   selBCtoEseeds.clear();
00075 
00076   //bending of traj. of charged particles under influence of B-field
00077   std::vector<reco::GenParticle> genParsCurved=applyBFieldCurv(genPars,iSetup);
00078 
00079   PhotElecgood=filterPhotonElectronSeed(clusterThreshold_,
00080                                         seedThreshold_,
00081                                         isoConeSize_,
00082                                         hOverEMax_,
00083                                         tkIsoMax_,
00084                                         caloIsoMax_,
00085                                         requireTrackMatch_,
00086                                         genPars,
00087                                         genParsCurved);
00088   
00089   IsoGenPargood=filterIsoGenPar(isoGenParETMin_,isoGenParConeSize_,genPars,genParsCurved);
00090 
00091 
00093   for (uint32_t ig=0;ig<genParsCurved.size();ig++) {
00094     reco::GenParticle gp=genParsCurved.at(ig);
00095     if (gp.status()==1 && abs(gp.pdgId())==11 && gp.et()>eTThreshold_ && fabs(gp.eta())<FILTER_ETA_MAX_) {
00096       if (hasBCAncestors(gp)) {
00097         BCtoEgood++;
00098         selBCtoEseeds.push_back(gp);
00099       }
00100     }
00101   }
00102 
00103   // we want 2 different em candidates
00104   if ( PhotElecgood> 1) { 
00105     result = true;
00106   }
00107   else if ( IsoGenPargood> 1) {
00108     result = true;  
00109   }
00110   else if ( BCtoEgood > 1) {
00111     result = true;
00112   }
00113   else if ( PhotElecgood==1 && IsoGenPargood==1) {
00114     if ( (sel1seeds.at(0).eta()!=sel2seeds.at(0).eta()) && (sel1seeds.at(0).phi()!=sel2seeds.at(0).phi()) && (sel1seeds.at(0).et()!=sel2seeds.at(0).et()) ) result = true;
00115   }
00116   else if ( PhotElecgood==1 && BCtoEgood==1) {
00117     if ( (sel1seeds.at(0).eta()!=selBCtoEseeds.at(0).eta()) && (sel1seeds.at(0).phi()!=selBCtoEseeds.at(0).phi()) && (sel1seeds.at(0).et()!=selBCtoEseeds.at(0).et()) ) result = true;
00118   }
00119   else if ( BCtoEgood==1 && IsoGenPargood==1) {
00120     if ( (selBCtoEseeds.at(0).eta()!=sel2seeds.at(0).eta()) && (selBCtoEseeds.at(0).phi()!=sel2seeds.at(0).phi()) && (selBCtoEseeds.at(0).et()!=sel2seeds.at(0).et()) ) result = true;
00121   }
00122   
00123   return result;
00124         
00125 }
00126 
00127 
00128 
00129 //filter that uses clustering around photon and electron seeds
00130 //only electrons, photons, charged pions, and charged kaons are clustered
00131 //additional requirements:
00132 
00133 
00134 //seed threshold, total threshold, and cone size/shape are specified separately for the barrel and the endcap
00135 int doubleEMEnrichingFilterAlgo::filterPhotonElectronSeed(float clusterthreshold,
00136                                                      float seedthreshold,
00137                                                      float isoConeSize,
00138                                                      float hOverEMax,
00139                                                      float tkIsoMax,
00140                                                      float caloIsoMax,
00141                                                      bool requiretrackmatch,
00142                                                      const std::vector<reco::GenParticle> &genPars,
00143                                                      const std::vector<reco::GenParticle> &genParsCurved) {
00144 
00145   float conesizeendcap=15;
00146 
00147   int retval=0;
00148   
00149   vector<reco::GenParticle> seeds;
00150   //find electron and photon seeds - must have E>seedthreshold GeV
00151   for (uint32_t is=0;is<genParsCurved.size();is++) {
00152     reco::GenParticle gp=genParsCurved.at(is);
00153     if (gp.status()!=1 || fabs(gp.eta()) > FILTER_ETA_MAX_ || fabs(gp.eta()) < FILTER_ETA_MIN_) continue;
00154     int absid=abs(gp.pdgId());
00155     if (absid!=11 && absid!=22) continue;
00156     if (gp.et()>seedthreshold) seeds.push_back(gp);
00157   }
00158   
00159   bool matchtrack=false;
00160   
00161   //for every seed, try to cluster stable particles about it in cone
00162   for (uint32_t is=0;is<seeds.size();is++) {
00163     float eTInCone=0;//eT associated to the electron cluster
00164     float tkIsoET=0;//tracker isolation energy
00165     float caloIsoET=0;//calorimeter isolation energy
00166     float hadET=0;//isolation energy from heavy hadrons that goes in the same area as the "electron" - so contributes to H/E
00167     bool isBarrel=fabs(seeds.at(is).eta())<ECALBARRELMAXETA_;
00168     for (uint32_t ig=0;ig<genParsCurved.size();ig++) {
00169       reco::GenParticle gp=genParsCurved.at(ig);
00170       reco::GenParticle gpUnCurv=genPars.at(ig);//for tk isolation, p at vertex
00171       if (gp.status()!=1) continue;
00172       int gpabsid=abs(gp.pdgId());
00173       if (gp.et()<1) continue;//ignore very soft particles
00174       //BARREL
00175       if (isBarrel) {
00176         float dr=deltaR(seeds.at(is),gp);
00177         float dphi=deltaPhi(seeds.at(is).phi(),gp.phi());
00178         float deta=fabs(seeds.at(is).eta()-gp.eta());
00179         if (deta<0.03 && dphi<0.2) {
00180           if (gpabsid==22 || gpabsid==11 || gpabsid==211 || gpabsid==321) {
00181             //contributes to electron
00182             eTInCone+=gp.et();
00183             //check for a matched track with at least 5 GeV
00184             if ((gpabsid==11 || gpabsid==211 || gpabsid==321) && gp.et()>5) matchtrack=true;
00185           } else {
00186             //contributes to H/E
00187             hadET+=gp.et();
00188           }
00189         } else {
00190           float drUnCurv=deltaR(seeds.at(is),gpUnCurv);
00191           if ((gp.charge()==0 && dr<isoConeSize && gpabsid!=22) ||
00192               (gp.charge()!=0 && drUnCurv<isoConeSize)) {      
00193             //contributes to calo isolation energy
00194             caloIsoET+=gp.et();
00195           }
00196           if (gp.charge()!=0 && drUnCurv<isoConeSize) {
00197             //contributes to track isolation energy
00198             tkIsoET+=gp.et();
00199           }
00200         }
00201         //ENDCAP
00202       } else {
00203         float drxy=deltaRxyAtEE(seeds.at(is),gp);
00204         float dr=deltaR(seeds.at(is),gp);//the isolation is done in dR
00205         if (drxy<conesizeendcap) {
00206           if (gpabsid==22 || gpabsid==11 || gpabsid==211 || gpabsid==321) {
00207             //contributes to electron
00208             eTInCone+=gp.et();
00209             //check for a matched track with at least 5 GeV
00210             if ((gpabsid==11 || gpabsid==211 || gpabsid==321) && gp.et()>5) matchtrack=true;
00211           } else {
00212             //contributes to H/E
00213             hadET+=gp.et();
00214           }
00215         } else {
00216           float drUnCurv=deltaR(seeds.at(is),gpUnCurv);
00217           if ((gp.charge()==0 && dr<isoConeSize && gpabsid!=22) ||
00218               (gp.charge()!=0 && drUnCurv<isoConeSize)) {
00219             //contributes to calo isolation energy
00220             caloIsoET+=gp.et();
00221           }
00222           if (gp.charge()!=0 && drUnCurv<isoConeSize) {
00223             //contributes to track isolation energy
00224             tkIsoET+=gp.et();
00225           }
00226         }
00227       }
00228     }
00229 
00230     if (eTInCone>clusterthreshold && (!requiretrackmatch || matchtrack)) {
00231       //       cout <<"isoET: "<<isoET<<endl;
00232       if (hadET/eTInCone<hOverEMax && tkIsoET<tkIsoMax && caloIsoET<caloIsoMax) { 
00233         retval=retval+1;
00234         sel1seeds.push_back(seeds[is]);
00235         // break;
00236       }
00237     }
00238   }
00239   
00240   return retval;
00241 }
00242 
00243 
00244 
00245 
00246 //make new genparticles vector taking into account the bending of charged particles in the b field
00247 //only stable-final-state (status==1) particles, with ET>=1 GeV, have their trajectories bent
00248 std::vector<reco::GenParticle> doubleEMEnrichingFilterAlgo::applyBFieldCurv(const std::vector<reco::GenParticle> &genPars, const edm::EventSetup& iSetup) {
00249   
00250   
00251   vector<reco::GenParticle> curvedPars;
00252 
00253   edm::ESHandle<MagneticField> magField;
00254   iSetup.get<IdealMagneticFieldRecord>().get(magField);
00255 
00256   Cylinder::CylinderPointer theBarrel=Cylinder::build(Surface::PositionType(0,0,0),Surface::RotationType(),ECALBARRELRADIUS_);
00257   Plane::PlanePointer endCapPlus=Plane::build(Surface::PositionType(0,0,ECALENDCAPZ_),Surface::RotationType());
00258   Plane::PlanePointer endCapMinus=Plane::build(Surface::PositionType(0,0,-1*ECALENDCAPZ_),Surface::RotationType());
00259 
00260   AnalyticalPropagator propagator(&(*magField),alongMomentum);  
00261 
00262   for (uint32_t ig=0;ig<genPars.size();ig++) {
00263     reco::GenParticle gp=genPars.at(ig);
00264     //don't bend trajectories of neutral particles, unstable particles, particles with < 1 GeV
00265     //particles with < ~0.9 GeV don't reach the barrel
00266     //so just put them as-is into the new vector
00267     if (gp.charge()==0 || gp.status()!=1 || gp.et()<1) {
00268       curvedPars.push_back(gp);
00269       continue;
00270     }
00271     GlobalPoint vertex(gp.vx(),gp.vy(),gp.vz());
00272     GlobalVector gvect(gp.px(),gp.py(),gp.pz());
00273     FreeTrajectoryState fts(vertex,gvect,gp.charge(),&(*magField));
00274     TrajectoryStateOnSurface impact;
00275     //choose to propagate to barrel, +Z endcap, or -Z endcap, according to eta
00276     if (fabs(gp.eta())<ECALBARRELMAXETA_) {
00277       impact=propagator.propagate(fts,*theBarrel);
00278     } else if (gp.eta()>0) {
00279       impact=propagator.propagate(fts,*endCapPlus);
00280     } else {
00281       impact=propagator.propagate(fts,*endCapMinus);
00282     }
00283     //in case the particle doesn't reach the barrel/endcap, just put it as-is into the new vector
00284     //it should reach though.
00285     if (!impact.isValid()) {
00286       curvedPars.push_back(gp);
00287       continue;
00288     }
00289     math::XYZTLorentzVector newp4;
00290 
00291     //the magnitude of p doesn't change, only the direction
00292     //NB I do get some small change in magnitude by the following,
00293     //I think it is a numerical precision issue
00294     float et=gp.et();
00295     float phinew=impact.globalPosition().phi();
00296     float pxnew=et*cos(phinew);
00297     float pynew=et*sin(phinew);
00298     newp4.SetPx(pxnew);
00299     newp4.SetPy(pynew);
00300     newp4.SetPz(gp.pz());
00301     newp4.SetE(gp.energy());
00302     reco::GenParticle gpnew=gp;
00303     gpnew.setP4(newp4);
00304     curvedPars.push_back(gpnew);
00305   }
00306   return curvedPars;
00307 
00308 
00309 }
00310 
00311 //calculate the difference in the xy-plane positions of gp1 and gp1 at the ECAL endcap
00312 //if they go in different z directions returns 9999.
00313 float doubleEMEnrichingFilterAlgo::deltaRxyAtEE(const reco::GenParticle &gp1, const reco::GenParticle &gp2) {
00314   
00315   if (gp1.pz()*gp2.pz() < 0) return 9999.;
00316   
00317   float rxy1=ECALENDCAPZ_*tan(gp1.theta());
00318   float x1=cos(gp1.phi())*rxy1;
00319   float y1=sin(gp1.phi())*rxy1;
00320 
00321   float rxy2=ECALENDCAPZ_*tan(gp2.theta());
00322   float x2=cos(gp2.phi())*rxy2;
00323   float y2=sin(gp2.phi())*rxy2;
00324   
00325   float dxy=sqrt((x2-x1)*(x2-x1)+(y2-y1)*(y2-y1));
00326   return dxy;
00327 
00328 
00329 }
00330 
00331 
00332 //filter looking for isolated charged pions, charged kaons, and electrons.
00333 //isolation done in cone of given size, looking at charged particles and neutral hadrons
00334 //photons aren't counted in the isolation requirements
00335 
00336 //need to have both the the curved and uncurved genpar collections
00337 //because tracker iso has to be treated differently than calo iso
00338 int doubleEMEnrichingFilterAlgo::filterIsoGenPar(float etmin, float conesize,const reco::GenParticleCollection &gph,
00339                                            const reco::GenParticleCollection &gphCurved
00340                                            ) {
00341 
00342   int passed = 0;
00343   for (uint32_t ip=0;ip<gph.size();ip++) {
00344 
00345     reco::GenParticle gp=gph.at(ip);
00346     reco::GenParticle gpCurved=gphCurved.at(ip);    
00347     int gpabsid=abs(gp.pdgId());
00348     //find potential candidates
00349     if (gp.et()<=etmin || gp.status()!=1) continue;
00350     if (gpabsid!=11 && gpabsid != 211 && gpabsid!= 321) continue;
00351     if (fabs(gp.eta()) < FILTER_ETA_MIN_) continue;
00352     if (fabs(gp.eta()) > FILTER_ETA_MAX_) continue;
00353     
00354     //check isolation
00355     float tkiso=0;
00356     float caloiso=0;
00357     for (uint32_t jp=0;jp<gph.size();jp++) {
00358       if (jp==ip) continue;
00359       reco::GenParticle pp=gph.at(jp);
00360       reco::GenParticle ppCurved=gphCurved.at(jp);
00361       if (pp.status() != 1) continue;
00362       float dr=deltaR(gp,pp);
00363       float drCurved=deltaR(gpCurved,ppCurved);
00364       if (abs(pp.charge())==1 && pp.et()>2 && dr<conesize) {
00365         tkiso+=pp.et();
00366       }
00367       if (pp.et()>2 && abs(pp.pdgId())!=22 && drCurved<conesize) {
00368         caloiso+=pp.energy();
00369       }
00370     }
00371     if (tkiso<FILTER_TKISOCUT_ && caloiso<FILTER_CALOISOCUT_) { 
00372       sel2seeds.push_back(gpCurved); 
00373       passed++;
00374     }
00375   }
00376   return passed;
00377 }
00378 
00379 //does this particle have an ancestor (mother, mother of mother, etc.) that is a b or c hadron?
00380 //attention: the GenParticle argument must have the motherRef correctly filled for this
00381 //to work.  That is, you had better have got it out of the genParticles collection
00382 bool doubleEMEnrichingFilterAlgo::hasBCAncestors(reco::GenParticle gp) {
00383 
00384   //stopping condition: this particle is a b or c hadron
00385   if (isBCHadron(gp)) return true;
00386   //stopping condition: this particle has no mothers
00387   if (gp.numberOfMothers()==0) return false;
00388   //otherwise continue
00389   bool retval=false;
00390   for (uint32_t im=0;im<gp.numberOfMothers();im++) {
00391     retval=retval || hasBCAncestors(*gp.motherRef(im));
00392   }
00393   return retval;
00394   
00395 }
00396 
00397 bool doubleEMEnrichingFilterAlgo::isBCHadron(reco::GenParticle gp) {
00398   return isBCMeson(gp) || isBCBaryon(gp);
00399 }
00400 
00401 bool doubleEMEnrichingFilterAlgo::isBCMeson(reco::GenParticle gp) {
00402   
00403   uint32_t pdgid=abs(gp.pdgId());
00404   uint32_t hundreds=pdgid%1000;
00405   if (hundreds>=400 && hundreds<600) {
00406     return true;
00407   } else {
00408     return false;
00409   }
00410 
00411 }
00412 
00413 bool doubleEMEnrichingFilterAlgo::isBCBaryon(reco::GenParticle gp) {
00414   
00415   uint32_t pdgid=abs(gp.pdgId());
00416   uint32_t thousands=pdgid%10000;
00417   if (thousands>=4000 && thousands <6000) {
00418     return true;
00419   } else {
00420     return false;
00421   }
00422 
00423 }