CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PythiaFilterIsolatedTrack.cc
Go to the documentation of this file.
6 
7 #include "CLHEP/Random/RandomEngine.h"
8 
9 #include <iostream>
10 #include<list>
11 #include<vector>
12 #include<cmath>
13 
14 //using namespace edm;
15 //using namespace std;
16 
17 std::pair<double,double> PythiaFilterIsolatedTrack::GetEtaPhiAtEcal(double etaIP, double phiIP, double pT, int charge, double vtxZ)
18 {
19  double deltaPhi;
20  double etaEC=100;
21  double phiEC=100;
22  double Rcurv=pT*33.3*100/38; //r(m)=pT(GeV)*33.3/B(kG)
23  double ecDist=317; //distance to ECAL andcap from IP (cm), 317 - ecal (not preshower), preshower -300
24  double ecRad=129; //radius of ECAL barrel (cm)
25  double theta=2*atan(exp(-etaIP));
26  double zNew;
27  if (theta>0.5*acos(-1)) theta=acos(-1)-theta;
28  if (fabs(etaIP)<1.479) {
29  deltaPhi=-charge*asin(0.5*ecRad/Rcurv);
30  double alpha1=2*asin(0.5*ecRad/Rcurv);
31  double z=ecRad/tan(theta);
32  if (etaIP>0) zNew=z*(Rcurv*alpha1)/ecRad+vtxZ; //new z-coordinate of track
33  else zNew=-z*(Rcurv*alpha1)/ecRad+vtxZ; //new z-coordinate of track
34  double zAbs=fabs(zNew);
35  if (zAbs<ecDist) {
36  etaEC=-log(tan(0.5*atan(ecRad/zAbs)));
37  deltaPhi=-charge*asin(0.5*ecRad/Rcurv);
38  }
39  if (zAbs>ecDist) {
40  zAbs=(fabs(etaIP)/etaIP)*ecDist;
41  double Zflight=fabs(zAbs-vtxZ);
42  double alpha=(Zflight*ecRad)/(z*Rcurv);
43  double Rec=2*Rcurv*sin(alpha/2);
44  deltaPhi=-charge*alpha/2;
45  etaEC=-log(tan(0.5*atan(Rec/ecDist)));
46  }
47  } else {
48  zNew=(fabs(etaIP)/etaIP)*ecDist;
49  double Zflight=fabs(zNew-vtxZ);
50  double Rvirt=fabs(Zflight*tan(theta));
51  double Rec=2*Rcurv*sin(Rvirt/(2*Rcurv));
52  deltaPhi=-(charge)*(Rvirt/(2*Rcurv));
53  etaEC=-log(tan(0.5*atan(Rec/ecDist)));
54  }
55 
56  if (zNew<0) etaEC=-etaEC;
57  phiEC=phiIP+deltaPhi;
58 
59  if (phiEC<-acos(-1)) phiEC=2*acos(-1)+phiEC;
60  if (phiEC>acos(-1)) phiEC=-2*acos(-1)+phiEC;
61 
62  std::pair<double,double> retVal(etaEC,phiEC);
63  return retVal;
64 }
65 
66 double PythiaFilterIsolatedTrack::getDistInCM(double eta1, double phi1, double eta2, double phi2)
67 {
68  double dR, Rec;
69  if (fabs(eta1)<1.479) Rec=129;
70  else Rec=317;
71  double ce1=cosh(eta1);
72  double ce2=cosh(eta2);
73  double te1=tanh(eta1);
74  double te2=tanh(eta2);
75 
76  double z=cos(phi1-phi2)/ce1/ce2+te1*te2;
77  if(z!=0) dR=fabs(Rec*ce1*sqrt(1./z/z-1.));
78  else dR=999999.;
79  return dR;
80 }
81 
83  token_(consumes<edm::HepMCProduct>(iConfig.getUntrackedParameter("ModuleLabel",edm::InputTag("generator","unsmeared")))),
84  MaxSeedEta_(iConfig.getUntrackedParameter<double>("MaxSeedEta", 2.3)),
85  MinSeedMom_(iConfig.getUntrackedParameter<double>("MinSeedMom", 20.)),
86  MinIsolTrackMom_(iConfig.getUntrackedParameter<double>("MinIsolTrackMom",2.0)),
87  IsolCone_(iConfig.getUntrackedParameter<double>("IsolCone", 40.0)),
88  PixelEfficiency_(iConfig.getUntrackedParameter<double>("PixelEfficiency", 0.8))
89 {
90 
91  // check if the random number generator service was configured
92  if(!rng_.isAvailable()) {
93  throw cms::Exception("Configuration") << "PythiaFilterIsolatedTrack requires the RandomNumberGeneratorService\n";
94  }
95 }
96 
97 
99 {
100 }
101 
102 
103 // ------------ method called to produce the data ------------
105 
106  CLHEP::HepRandomEngine& engine = rng_->getEngine(iEvent.streamID());
107 
109  iSetup.getData( pdt );
110 
112  iEvent.getByToken(token_, evt);
113 
114  const HepMC::GenEvent* myGenEvent = evt->GetEvent();
115 
116  // all of the stable, charged particles with momentum>MinIsolTrackMom_ and |eta|<MaxSeedEta_+0.5
117  std::vector<const HepMC::GenParticle *> chargedParticles;
118 
119  // all of the stable, charged particles with momentum>MinSeedMom_ and |eta|<MaxSeedEta_
120  std::vector<const HepMC::GenParticle *> seeds;
121 
122  // fill the vector of charged particles and seeds in the event
123  for(HepMC::GenEvent::particle_const_iterator iter=myGenEvent->particles_begin(); iter!=myGenEvent->particles_end(); ++iter) {
124  const HepMC::GenParticle *p=*iter;
125  int charge3 = pdt->particle(p->pdg_id())->ID().threeCharge();
126  int status = p->status();
127  double momentum = p->momentum().rho();
128  double abseta = fabs(p->momentum().eta());
129 
130  // only consider stable, charged particles
131  if(abs(charge3)==3 && status==1 && momentum>MinIsolTrackMom_ && abseta<MaxSeedEta_+0.5) {
132  chargedParticles.push_back(p);
133  if(momentum>MinSeedMom_ && abseta<MaxSeedEta_) {
134  seeds.push_back(p);
135  }
136  }
137  }
138 
139  // loop over all the seeds and see if any of them are isolated
140  for(std::vector<const HepMC::GenParticle *>::const_iterator it1=seeds.begin(); it1!=seeds.end(); ++it1) {
141  const HepMC::GenParticle *p1=*it1;
142 
143  std::pair<double,double> EtaPhi1=GetEtaPhiAtEcal(p1->momentum().eta(),
144  p1->momentum().phi(),
145  p1->momentum().perp(),
146  (pdt->particle(p1->pdg_id()))->ID().threeCharge()/3,
147  0.0);
148 
149  // loop over all of the other charged particles in the event, and see if any are close by
150  bool failsIso=false;
151  for(std::vector<const HepMC::GenParticle *>::const_iterator it2=chargedParticles.begin(); it2!=chargedParticles.end(); ++it2) {
152  const HepMC::GenParticle *p2=*it2;
153 
154  // don't consider the seed particle among the other charge particles
155  if(p1==p2) continue;
156 
157  std::pair<double,double> EtaPhi2=GetEtaPhiAtEcal(p2->momentum().eta(),
158  p2->momentum().phi(),
159  p2->momentum().perp(),
160  (pdt->particle(p2->pdg_id()))->ID().threeCharge()/3,
161  0.0);
162 
163  // find out how far apart the particles are
164  // if the seed fails the isolation requirement, try a different seed
165  // occasionally allow a seed to pass to isolation requirement
166  if(getDistInCM(EtaPhi1.first, EtaPhi1.second, EtaPhi2.first, EtaPhi2.second) < IsolCone_ &&
167  engine.flat() < PixelEfficiency_) {
168  failsIso=true;
169  break;
170  }
171  }
172 
173  if(!failsIso) return true;
174 
175  } //loop over seeds
176 
177  return false;
178 }
float alpha
Definition: AMPTWrapper.h:95
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:464
uint32_t ID
Definition: Definitions.h:26
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
Geom::Theta< T > theta() const
void getData(T &iHolder) const
Definition: EventSetup.h:79
int iEvent
Definition: GenABIO.cc:230
T sqrt(T t)
Definition: SSEVec.h:48
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
bool isAvailable() const
Definition: Service.h:46
virtual bool filter(edm::Event &, const edm::EventSetup &)
edm::EDGetTokenT< edm::HepMCProduct > token_
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double p2[4]
Definition: TauolaWrapper.h:90
static double getDistInCM(double eta1, double phi1, double eta2, double phi2)
PythiaFilterIsolatedTrack(const edm::ParameterSet &)
virtual CLHEP::HepRandomEngine & getEngine(StreamID const &) const =0
Use this engine in event methods.
edm::Service< edm::RandomNumberGenerator > rng_
double p1[4]
Definition: TauolaWrapper.h:89
StreamID streamID() const
Definition: Event.h:79
tuple status
Definition: ntuplemaker.py:245
static std::pair< double, double > GetEtaPhiAtEcal(double etaIP, double phiIP, double pT, int charge, double vtxZ)
tuple log
Definition: cmsBatch.py:341