Go to the documentation of this file.00001 #include "GeneratorInterface/GenFilters/interface/PythiaFilterIsolatedTrack.h"
00002 #include "SimDataFormats/GeneratorProducts/interface/HepMCProduct.h"
00003 #include "SimGeneral/HepPDTRecord/interface/ParticleDataTable.h"
00004 #include "FWCore/Framework/interface/EventSetup.h"
00005 #include "FWCore/Framework/interface/ESHandle.h"
00006 #include <iostream>
00007 #include<list>
00008 #include<vector>
00009 #include<cmath>
00010
00011
00012
00013
00014 std::pair<double,double> PythiaFilterIsolatedTrack::GetEtaPhiAtEcal(double etaIP, double phiIP, double pT, int charge, double vtxZ)
00015 {
00016 double deltaPhi;
00017 double etaEC=100;
00018 double phiEC=100;
00019 double Rcurv=pT*33.3*100/38;
00020 double ecDist=317;
00021 double ecRad=129;
00022 double theta=2*atan(exp(-etaIP));
00023 double zNew;
00024 if (theta>0.5*acos(-1)) theta=acos(-1)-theta;
00025 if (fabs(etaIP)<1.479) {
00026 deltaPhi=-charge*asin(0.5*ecRad/Rcurv);
00027 double alpha1=2*asin(0.5*ecRad/Rcurv);
00028 double z=ecRad/tan(theta);
00029 if (etaIP>0) zNew=z*(Rcurv*alpha1)/ecRad+vtxZ;
00030 else zNew=-z*(Rcurv*alpha1)/ecRad+vtxZ;
00031 double zAbs=fabs(zNew);
00032 if (zAbs<ecDist) {
00033 etaEC=-log(tan(0.5*atan(ecRad/zAbs)));
00034 deltaPhi=-charge*asin(0.5*ecRad/Rcurv);
00035 }
00036 if (zAbs>ecDist) {
00037 zAbs=(fabs(etaIP)/etaIP)*ecDist;
00038 double Zflight=fabs(zAbs-vtxZ);
00039 double alpha=(Zflight*ecRad)/(z*Rcurv);
00040 double Rec=2*Rcurv*sin(alpha/2);
00041 deltaPhi=-charge*alpha/2;
00042 etaEC=-log(tan(0.5*atan(Rec/ecDist)));
00043 }
00044 } else {
00045 zNew=(fabs(etaIP)/etaIP)*ecDist;
00046 double Zflight=fabs(zNew-vtxZ);
00047 double Rvirt=fabs(Zflight*tan(theta));
00048 double Rec=2*Rcurv*sin(Rvirt/(2*Rcurv));
00049 deltaPhi=-(charge)*(Rvirt/(2*Rcurv));
00050 etaEC=-log(tan(0.5*atan(Rec/ecDist)));
00051 }
00052
00053 if (zNew<0) etaEC=-etaEC;
00054 phiEC=phiIP+deltaPhi;
00055
00056 if (phiEC<-acos(-1)) phiEC=2*acos(-1)+phiEC;
00057 if (phiEC>acos(-1)) phiEC=-2*acos(-1)+phiEC;
00058
00059 std::pair<double,double> retVal(etaEC,phiEC);
00060 return retVal;
00061 }
00062
00063 double PythiaFilterIsolatedTrack::getDistInCM(double eta1, double phi1, double eta2, double phi2)
00064 {
00065 double dR, Rec;
00066 if (fabs(eta1)<1.479) Rec=129;
00067 else Rec=317;
00068 double ce1=cosh(eta1);
00069 double ce2=cosh(eta2);
00070 double te1=tanh(eta1);
00071 double te2=tanh(eta2);
00072
00073 double z=cos(phi1-phi2)/ce1/ce2+te1*te2;
00074 if(z!=0) dR=fabs(Rec*ce1*sqrt(1./z/z-1.));
00075 else dR=999999.;
00076 return dR;
00077 }
00078
00079 PythiaFilterIsolatedTrack::PythiaFilterIsolatedTrack(const edm::ParameterSet& iConfig) :
00080 ModuleLabel_(iConfig.getUntrackedParameter("ModuleLabel",std::string("generator"))),
00081 MaxSeedEta_(iConfig.getUntrackedParameter<double>("MaxSeedEta", 2.3)),
00082 MinSeedMom_(iConfig.getUntrackedParameter<double>("MinSeedMom", 20.)),
00083 MinIsolTrackMom_(iConfig.getUntrackedParameter<double>("MinIsolTrackMom",2.0)),
00084 IsolCone_(iConfig.getUntrackedParameter<double>("IsolCone", 40.0)),
00085 PixelEfficiency_(iConfig.getUntrackedParameter<double>("PixelEfficiency", 0.8))
00086 {
00087
00088
00089 if(!rng_.isAvailable()) {
00090 throw cms::Exception("Configuration") << "PythiaFilterIsolatedTrack requires the RandomNumberGeneratorService\n";
00091 }
00092 CLHEP::HepRandomEngine& engine = rng_->getEngine();
00093 flatDistribution_ = new CLHEP::RandFlat(engine, 0.0, 1.0);
00094 }
00095
00096
00097 PythiaFilterIsolatedTrack::~PythiaFilterIsolatedTrack()
00098 {
00099 delete flatDistribution_;
00100 }
00101
00102
00103
00104 bool PythiaFilterIsolatedTrack::filter(edm::Event& iEvent, const edm::EventSetup& iSetup){
00105
00106 edm::ESHandle<ParticleDataTable> pdt;
00107 iSetup.getData( pdt );
00108
00109 edm::Handle<edm::HepMCProduct> evt;
00110 iEvent.getByLabel(ModuleLabel_, evt);
00111
00112 const HepMC::GenEvent* myGenEvent = evt->GetEvent();
00113
00114
00115 std::vector<const HepMC::GenParticle *> chargedParticles;
00116
00117
00118 std::vector<const HepMC::GenParticle *> seeds;
00119
00120
00121 for(HepMC::GenEvent::particle_const_iterator iter=myGenEvent->particles_begin(); iter!=myGenEvent->particles_end(); ++iter) {
00122 const HepMC::GenParticle *p=*iter;
00123 int charge3 = pdt->particle(p->pdg_id())->ID().threeCharge();
00124 int status = p->status();
00125 double momentum = p->momentum().rho();
00126 double abseta = fabs(p->momentum().eta());
00127
00128
00129 if(abs(charge3)==3 && status==1 && momentum>MinIsolTrackMom_ && abseta<MaxSeedEta_+0.5) {
00130 chargedParticles.push_back(p);
00131 if(momentum>MinSeedMom_ && abseta<MaxSeedEta_) {
00132 seeds.push_back(p);
00133 }
00134 }
00135 }
00136
00137
00138 for(std::vector<const HepMC::GenParticle *>::const_iterator it1=seeds.begin(); it1!=seeds.end(); ++it1) {
00139 const HepMC::GenParticle *p1=*it1;
00140
00141 std::pair<double,double> EtaPhi1=GetEtaPhiAtEcal(p1->momentum().eta(),
00142 p1->momentum().phi(),
00143 p1->momentum().perp(),
00144 (pdt->particle(p1->pdg_id()))->ID().threeCharge()/3,
00145 0.0);
00146
00147
00148 bool failsIso=false;
00149 for(std::vector<const HepMC::GenParticle *>::const_iterator it2=chargedParticles.begin(); it2!=chargedParticles.end(); ++it2) {
00150 const HepMC::GenParticle *p2=*it2;
00151
00152
00153 if(p1==p2) continue;
00154
00155 std::pair<double,double> EtaPhi2=GetEtaPhiAtEcal(p2->momentum().eta(),
00156 p2->momentum().phi(),
00157 p2->momentum().perp(),
00158 (pdt->particle(p2->pdg_id()))->ID().threeCharge()/3,
00159 0.0);
00160
00161
00162
00163
00164 if(getDistInCM(EtaPhi1.first, EtaPhi1.second, EtaPhi2.first, EtaPhi2.second) < IsolCone_ &&
00165 flatDistribution_->fire() < PixelEfficiency_) {
00166 failsIso=true;
00167 break;
00168 }
00169 }
00170
00171 if(!failsIso) return true;
00172
00173 }
00174
00175 return false;
00176 }