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