CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
TrackParameterAnalyzer.cc
Go to the documentation of this file.
2 #include <string>
3 #include <vector>
7 
8 //
9 //
10 // constants, enums and typedefs
11 //
12 
13 //
14 // static data member definitions
15 //
16 
17 //
18 // constructors and destructor
19 //
21  simG4_( iConfig.getParameter<edm::InputTag>( "simG4" ) ) {
22  //now do whatever initialization is needed
23 
24  recoTrackProducer_ = iConfig.getUntrackedParameter<std::string>("recoTrackProducer");
25  // open output file to store histograms}
26  outputFile_ = iConfig.getUntrackedParameter<std::string>("outputFile");
27  TString tversion(edm::getReleaseVersion());
28  tversion = tversion.Remove(0,1);
29  tversion = tversion.Remove(tversion.Length()-1,tversion.Length());
30  outputFile_ = std::string(tversion)+"_"+outputFile_;
31 
32  rootFile_ = TFile::Open(outputFile_.c_str(),"RECREATE");
33  verbose_= iConfig.getUntrackedParameter<bool>("verbose", false);
34  simUnit_=1.0; // starting from CMSSW_1_2_x, I think
35  if ( (edm::getReleaseVersion()).find("CMSSW_1_1_",0)!=std::string::npos){
36  simUnit_=0.1; // for use in CMSSW_1_1_1 tutorial
37  }
38 }
39 
40 
42 {
43 
44  // do anything here that needs to be done at destruction time
45  // (e.g. close files, deallocate resources etc.)
46  delete rootFile_;
47 }
48 
49 
50 
51 //
52 // member functions
53 //
55  std::cout << " TrackParameterAnalyzer::beginJob conversion from sim units to rec units is " << simUnit_ << std::endl;
56 
57  rootFile_->cd();
58  h1_pull0_ = new TH1F("pull0","pull q/p",100,-25.,25.);
59  h1_pull1_ = new TH1F("pull1","pull lambda",100,-25.,25.);
60  h1_pull2_ = new TH1F("pull2","pull phi ",100,-25.,25.);
61  h1_pull3_ = new TH1F("pull3","pull dca ",100,-25.,25.);
62  h1_pull4_ = new TH1F("pull4","pull zdca ",100,-25.,25.);
63 
64  h1_res0_ = new TH1F("res0","res q/p",100,-0.1,0.1);
65  h1_res1_ = new TH1F("res1","res lambda",100,-0.1,0.1);
66  h1_res2_ = new TH1F("res2","res phi ",100,-0.1,0.1);
67  h1_res3_ = new TH1F("res3","res dca ",100,-0.1,0.1);
68  h1_res4_ = new TH1F("res4","res zdca ",100,-0.1,0.1);
69 
70  h1_Beff_ = new TH1F("Beff", "Beff",2000,-10.,10.);
71  h2_dvsphi_ = new TH2F("dvsphi","dvsphi",360,-M_PI,M_PI,100,-0.1,0.1);
72  h1_par0_ = new TH1F("par0","q/p",100,-0.1,0.1);
73  h1_par1_ = new TH1F("par1","lambda",100, -M_PI/2.,M_PI/2.);
74  h1_par2_ = new TH1F("par2","phi ",100,-M_PI,M_PI);
75  h1_par3_ = new TH1F("par3","dca ",100,-0.1,0.1);
76  h1_par4_ = new TH1F("par4","zdca ",1000,-10.,10.);
77 }
78 
79 
81  rootFile_->cd();
82  h1_pull0_->Write();
83  h1_pull1_->Write();
84  h1_pull2_->Write();
85  h1_pull3_->Write();
86  h1_pull4_->Write();
87 
88  h1_res0_->Write();
89  h1_res1_->Write();
90  h1_res2_->Write();
91  h1_res3_->Write();
92  h1_res4_->Write();
93 
94  h1_Beff_->Write();
95  h2_dvsphi_->Write();
96  h1_par0_->Write();
97  h1_par1_->Write();
98  h1_par2_->Write();
99  h1_par3_->Write();
100  h1_par4_->Write();
101 }
102 
103 // helper function
105  const ParameterVector &b){
106  double dtheta=a(1)-b(1);
107  double dphi =a(2)-b(2);
108  if (dphi>M_PI){ dphi-=M_2_PI; }else if(dphi<-M_PI){dphi+=M_2_PI;}
109  return ( (fabs(dtheta)<0.02) && (fabs(dphi)<0.04) );
110 }
111 
112 // ------------ method called to produce the data ------------
113 void
115 {
116  using namespace edm;
117  using CLHEP::HepLorentzVector;
118 
119  //edm::ESHandle<TransientTrackBuilder> theB;
120  //iSetup.get<TransientTrackRecord>().get("TransientTrackBuilder",theB);
121  //double fBfield=((*theB).field()->inTesla(GlobalPoint(0.,0.,0.))).z();
122  const double fBfield=3.8;
123 
125  iEvent.getByLabel( simG4_, simVtcs);
126  if(verbose_){
127  std::cout << "SimVertex " << simVtcs->size() << std::endl;
128  for(edm::SimVertexContainer::const_iterator v=simVtcs->begin();
129  v!=simVtcs->end(); ++v){
130  std::cout << "simvtx "
131  << std::setw(10) << std::setprecision(4)
132  << v->position().x() << " "
133  << v->position().y() << " "
134  << v->position().z() << " "
135  << v->parentIndex() << " "
136  << v->noParent() << " "
137  << std::endl;
138  }
139  }
140 
141  // get the simulated tracks, extract perigee parameters
143  iEvent.getByLabel( simG4_, simTrks);
144 
145  if(verbose_){std::cout << "simtrks " << simTrks->size() << std::endl;}
146  std::vector<ParameterVector > tsim;
147  for(edm::SimTrackContainer::const_iterator t=simTrks->begin();
148  t!=simTrks->end(); ++t){
149  if (t->noVertex()){
150  std::cout << "simtrk has no vertex" << std::endl;
151  return;
152  }else{
153  // get the vertex position
154  HepLorentzVector v((*simVtcs)[t->vertIndex()].position().x(),
155  (*simVtcs)[t->vertIndex()].position().y(),
156  (*simVtcs)[t->vertIndex()].position().z(),
157  (*simVtcs)[t->vertIndex()].position().e());
158  int pdgCode=t->type();
159 
160  if( pdgCode==-99 ){
161  // such entries cause crashes, no idea what they are
162  std::cout << "funny particle skipped , code=" << pdgCode << std::endl;
163  }else{
164  double Q=0; //double Q=HepPDT::theTable().getParticleData(pdgCode)->charge();
165  if ((pdgCode==11)||(pdgCode==13)||(pdgCode==15)||(pdgCode==-211)||(pdgCode==-2212)||(pdgCode==321)){Q=-1;}
166  else if((pdgCode==-11)||(pdgCode==-13)||(pdgCode==-15)||(pdgCode==211)||(pdgCode==2212)||(pdgCode==321)){Q=1;}
167  else {
168  std::cout << pdgCode << " " <<std::endl;
169  }
170  HepLorentzVector p(t->momentum().x(),t->momentum().y(),t->momentum().z(),t->momentum().e());
171  if(verbose_){
172  std::cout << "simtrk "
173  << " gen=" << std::setw( 4) << t->genpartIndex()
174  << " vtx=" << std::setw( 4) << t->vertIndex()
175  << " pdg=" << std::setw( 5) << t->type()
176  << " Q=" << std::setw( 3) << Q
177  << " pt=" << std::setw(6) << p.perp()
178  << " vx=" << std::setw(6) << v.x()
179  << " vy=" << std::setw(6) << v.y()
180  << " vz=" << std::setw(6) << v.z()
181  << std::endl;
182  }
183  if ( (Q != 0) && (p.perp()>0.1) ){
184  double x0=v.x()*simUnit_;
185  double y0=v.y()*simUnit_;
186  double z0=v.z()*simUnit_;
187  double kappa=-Q*0.002998*fBfield/p.perp();
188  double D0=x0*sin(p.phi())-y0*cos(p.phi())-0.5*kappa*(x0*x0+y0*y0);
189  double q=sqrt(1.-2.*kappa*D0);
190  double s0=(x0*cos(p.phi())+y0*sin(p.phi()))/q;
191  double s1;
192  if (fabs(kappa*s0)>0.001){
193  s1=asin(kappa*s0)/kappa;
194  }else{
195  double ks02=(kappa*s0)*(kappa*s0);
196  s1=s0*(1.+ks02/6.+3./40.*ks02*ks02+5./112.*pow(ks02,3));
197  }
198  ParameterVector par;
199  par[reco::TrackBase::i_qoverp] = Q/sqrt(p.perp2()+p.pz()*p.pz());
200  par[reco::TrackBase::i_lambda] = M_PI/2.-p.theta();
201  par[reco::TrackBase::i_phi] = p.phi()-asin(kappa*s0);
202  par[reco::TrackBase::i_dxy] = 2.*D0/(1.+q);
203  par[reco::TrackBase::i_dsz] = z0*sin(p.theta())-s1*cos(p.theta());
204  tsim.push_back(par);
205  }
206  }
207  }// has vertex
208  }//for loop
209 
210  // simtrack parameters are in now tsim
211  // loop over tracks and try to match them to simulated tracks
212 
213 
215  iEvent.getByLabel(recoTrackProducer_, recTracks);
216 
217  for(reco::TrackCollection::const_iterator t=recTracks->begin();
218  t!=recTracks->end(); ++t){
219  reco::TrackBase::ParameterVector p = t->parameters();
220  reco::TrackBase::CovarianceMatrix c = t->covariance();
221  if(verbose_){
222  std::cout << "reco pars= " << p << std::endl;
223  //std::cout << "z0=" << p(4) << std::endl;
224  }
225  for(std::vector<ParameterVector>::const_iterator s=tsim.begin();
226  s!=tsim.end(); ++s){
227  if (match(*s,p)){
228  //if(verbose_){ std::cout << "match found" << std::endl;}
229  h1_pull0_->Fill((p(0)-(*s)(0))/sqrt(c(0,0)));
230  h1_pull1_->Fill((p(1)-(*s)(1))/sqrt(c(1,1)));
231  h1_pull2_->Fill((p(2)-(*s)(2))/sqrt(c(2,2)));
232  h1_pull3_->Fill((p(3)-(*s)(3))/sqrt(c(3,3)));
233  h1_pull4_->Fill((p(4)-(*s)(4))/sqrt(c(4,4)));
234 
235  h1_res0_->Fill(p(0)-(*s)(0));
236  h1_res1_->Fill(p(1)-(*s)(1));
237  h1_res2_->Fill(p(2)-(*s)(2));
238  h1_res3_->Fill(p(3)-(*s)(3));
239  h1_res4_->Fill(p(4)-(*s)(4));
240 
241  h1_Beff_->Fill(p(0)/(*s)(0)*fBfield);
242  h2_dvsphi_->Fill(p(2),p(3));
243  h1_par0_->Fill(p(0));
244  h1_par1_->Fill(p(1));
245  h1_par2_->Fill(p(2));
246  h1_par3_->Fill(p(3));
247  h1_par4_->Fill(p(4));
248  }
249  }
250  }
251 
252 
253 
254 
255 }
256 
T getUntrackedParameter(std::string const &, T const &) const
Divides< arg, void > D0
Definition: Factorize.h:143
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
TrackParameterAnalyzer(const edm::ParameterSet &)
math::Vector< dimension >::type ParameterVector
parameter vector
Definition: TrackBase.h:70
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:7
bool match(const ParameterVector &a, const ParameterVector &b)
int iEvent
Definition: GenABIO.cc:243
T sqrt(T t)
Definition: SSEVec.h:46
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
virtual void analyze(const edm::Event &, const edm::EventSetup &)
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:356
std::string getReleaseVersion()
#define M_PI
Definition: BFit3D.cc:3
double b
Definition: hdecay.h:120
reco::TrackBase::ParameterVector ParameterVector
double a
Definition: hdecay.h:121
tuple cout
Definition: gather_cfg.py:121
mathSSE::Vec4< T > v
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
math::Error< dimension >::type CovarianceMatrix
5 parameter covariance matrix
Definition: TrackBase.h:72