CMS 3D CMS Logo

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>(
40  edm::InputTag(iConfig.getUntrackedParameter<std::string>("recoTrackProducer")))),
41  outputFile_(iConfig.getUntrackedParameter<std::string>("outputFile")),
42  simUnit_(1.0) // starting from CMSSW_1_2_x, I think
43  ,
44  verbose_(iConfig.getUntrackedParameter<bool>("verbose", false)) {
45  //now do whatever initialization is needed
46  // open output file to store histograms}
47  auto tversion = edm::getReleaseVersion();
48  tversion = tversion.erase(tversion.size() - 1, 1).erase(0, 1);
49  outputFile_ = tversion + "_" + outputFile_;
50  rootFile_ = TFile::Open(outputFile_.c_str(), "RECREATE");
51  if ((edm::getReleaseVersion()).find("CMSSW_1_1_", 0) != std::string::npos) {
52  simUnit_ = 0.1; // for use in CMSSW_1_1_1 tutorial
53  }
54 }
55 
57  // do anything here that needs to be done at destruction time
58  // (e.g. close files, deallocate resources etc.)
59  delete rootFile_;
60 }
61 
62 //
63 // member functions
64 //
66  std::cout << " TrackParameterAnalyzer::beginJob conversion from sim units to rec units is " << simUnit_ << std::endl;
67 
68  rootFile_->cd();
69  h1_pull0_ = new TH1F("pull0", "pull q/p", 100, -25., 25.);
70  h1_pull1_ = new TH1F("pull1", "pull lambda", 100, -25., 25.);
71  h1_pull2_ = new TH1F("pull2", "pull phi ", 100, -25., 25.);
72  h1_pull3_ = new TH1F("pull3", "pull dca ", 100, -25., 25.);
73  h1_pull4_ = new TH1F("pull4", "pull zdca ", 100, -25., 25.);
74 
75  h1_res0_ = new TH1F("res0", "res q/p", 100, -0.1, 0.1);
76  h1_res1_ = new TH1F("res1", "res lambda", 100, -0.1, 0.1);
77  h1_res2_ = new TH1F("res2", "res phi ", 100, -0.1, 0.1);
78  h1_res3_ = new TH1F("res3", "res dca ", 100, -0.1, 0.1);
79  h1_res4_ = new TH1F("res4", "res zdca ", 100, -0.1, 0.1);
80 
81  h1_Beff_ = new TH1F("Beff", "Beff", 2000, -10., 10.);
82  h2_dvsphi_ = new TH2F("dvsphi", "dvsphi", 360, -M_PI, M_PI, 100, -0.1, 0.1);
83  h1_par0_ = new TH1F("par0", "q/p", 100, -0.1, 0.1);
84  h1_par1_ = new TH1F("par1", "lambda", 100, -M_PI / 2., M_PI / 2.);
85  h1_par2_ = new TH1F("par2", "phi ", 100, -M_PI, M_PI);
86  h1_par3_ = new TH1F("par3", "dca ", 100, -0.1, 0.1);
87  h1_par4_ = new TH1F("par4", "zdca ", 1000, -10., 10.);
88 }
89 
91  rootFile_->cd();
92  h1_pull0_->Write();
93  h1_pull1_->Write();
94  h1_pull2_->Write();
95  h1_pull3_->Write();
96  h1_pull4_->Write();
97 
98  h1_res0_->Write();
99  h1_res1_->Write();
100  h1_res2_->Write();
101  h1_res3_->Write();
102  h1_res4_->Write();
103 
104  h1_Beff_->Write();
105  h2_dvsphi_->Write();
106  h1_par0_->Write();
107  h1_par1_->Write();
108  h1_par2_->Write();
109  h1_par3_->Write();
110  h1_par4_->Write();
111 }
112 
113 // helper function
115  double dtheta = a(1) - b(1);
116  double dphi = a(2) - b(2);
117  if (dphi > M_PI) {
118  dphi -= M_2_PI;
119  } else if (dphi < -M_PI) {
120  dphi += M_2_PI;
121  }
122  return ((fabs(dtheta) < 0.02) && (fabs(dphi) < 0.04));
123 }
124 
125 // ------------ method called to produce the data ------------
127  using CLHEP::HepLorentzVector;
128 
129  const double fBfield = 3.8;
130 
132  iEvent.getByToken(edmSimVertexContainerToken_, simVtcs);
133  if (verbose_) {
134  std::cout << "SimVertex " << simVtcs->size() << std::endl;
135  for (edm::SimVertexContainer::const_iterator v = simVtcs->begin(); v != simVtcs->end(); ++v) {
136  std::cout << "simvtx " << std::setw(10) << std::setprecision(4) << v->position().x() << " " << v->position().y()
137  << " " << v->position().z() << " " << v->parentIndex() << " " << v->noParent() << " " << std::endl;
138  }
139  }
140 
141  // get the simulated tracks, extract perigee parameters
143  iEvent.getByToken(edmSimTrackContainerToken_, simTrks);
144 
145  if (verbose_) {
146  std::cout << "simtrks " << simTrks->size() << std::endl;
147  }
148  std::vector<ParameterVector> tsim;
149  for (edm::SimTrackContainer::const_iterator t = simTrks->begin(); t != simTrks->end(); ++t) {
150  if (t->noVertex()) {
151  std::cout << "simtrk has no vertex" << std::endl;
152  return;
153  } else {
154  // get the vertex position
155  HepLorentzVector v((*simVtcs)[t->vertIndex()].position().x(),
156  (*simVtcs)[t->vertIndex()].position().y(),
157  (*simVtcs)[t->vertIndex()].position().z(),
158  (*simVtcs)[t->vertIndex()].position().e());
159  int pdgCode = t->type();
160 
161  if (pdgCode == -99) {
162  // such entries cause crashes, no idea what they are
163  std::cout << "funny particle skipped , code=" << pdgCode << std::endl;
164  } else {
165  double Q = 0;
166  if ((pdgCode == 11) || (pdgCode == 13) || (pdgCode == 15) || (pdgCode == -211) || (pdgCode == -2212) ||
167  (pdgCode == 321)) {
168  Q = -1;
169  } else if ((pdgCode == -11) || (pdgCode == -13) || (pdgCode == -15) || (pdgCode == 211) || (pdgCode == 2212) ||
170  (pdgCode == 321)) {
171  Q = 1;
172  } else {
173  std::cout << pdgCode << " " << std::endl;
174  }
175  HepLorentzVector p(t->momentum().x(), t->momentum().y(), t->momentum().z(), t->momentum().e());
176  if (verbose_) {
177  std::cout << "simtrk "
178  << " gen=" << std::setw(4) << t->genpartIndex() << " vtx=" << std::setw(4) << t->vertIndex()
179  << " pdg=" << std::setw(5) << t->type() << " Q=" << std::setw(3) << Q << " pt=" << std::setw(6)
180  << p.perp() << " vx=" << std::setw(6) << v.x() << " vy=" << std::setw(6) << v.y()
181  << " vz=" << std::setw(6) << v.z() << 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 
215 
216  for (reco::TrackCollection::const_iterator t = recTracks->begin(); t != recTracks->end(); ++t) {
217  reco::TrackBase::ParameterVector p = t->parameters();
218  reco::TrackBase::CovarianceMatrix c = t->covariance();
219  if (verbose_) {
220  std::cout << "reco pars= " << p << std::endl;
221  }
222  for (std::vector<ParameterVector>::const_iterator s = tsim.begin(); s != tsim.end(); ++s) {
223  if (match(*s, p)) {
224  h1_pull0_->Fill((p(0) - (*s)(0)) / sqrt(c(0, 0)));
225  h1_pull1_->Fill((p(1) - (*s)(1)) / sqrt(c(1, 1)));
226  h1_pull2_->Fill((p(2) - (*s)(2)) / sqrt(c(2, 2)));
227  h1_pull3_->Fill((p(3) - (*s)(3)) / sqrt(c(3, 3)));
228  h1_pull4_->Fill((p(4) - (*s)(4)) / sqrt(c(4, 4)));
229 
230  h1_res0_->Fill(p(0) - (*s)(0));
231  h1_res1_->Fill(p(1) - (*s)(1));
232  h1_res2_->Fill(p(2) - (*s)(2));
233  h1_res3_->Fill(p(3) - (*s)(3));
234  h1_res4_->Fill(p(4) - (*s)(4));
235 
236  h1_Beff_->Fill(p(0) / (*s)(0) * fBfield);
237  h2_dvsphi_->Fill(p(2), p(3));
238  h1_par0_->Fill(p(0));
239  h1_par1_->Fill(p(1));
240  h1_par2_->Fill(p(2));
241  h1_par3_->Fill(p(3));
242  h1_par4_->Fill(p(4));
243  }
244  }
245  }
246 }
TrackParameterAnalyzer::h1_pull2_
TH1 * h1_pull2_
Definition: TrackParameterAnalyzer.h:66
TrackParameterAnalyzer::edmSimTrackContainerToken_
edm::EDGetTokenT< edm::SimTrackContainer > edmSimTrackContainerToken_
Definition: TrackParameterAnalyzer.h:59
electrons_cff.bool
bool
Definition: electrons_cff.py:393
MessageLogger.h
funct::false
false
Definition: Factorize.h:29
TrackParameterAnalyzer::h1_res2_
TH1 * h1_res2_
Definition: TrackParameterAnalyzer.h:71
TrackParameterAnalyzer.h
funct::D0
Divides< arg, void > D0
Definition: Factorize.h:135
TrackParameterAnalyzer::~TrackParameterAnalyzer
~TrackParameterAnalyzer() override
Definition: TrackParameterAnalyzer.cc:56
reco::TrackBase::i_phi
Definition: TrackBase.h:83
edm
HLT enums.
Definition: AlignableModifier.h:19
TrackParameterAnalyzer::h1_pull4_
TH1 * h1_pull4_
Definition: TrackParameterAnalyzer.h:68
AlCaHLTBitMon_ParallelJobs.p
p
Definition: AlCaHLTBitMon_ParallelJobs.py:153
gather_cfg.cout
cout
Definition: gather_cfg.py:144
TrackParameterAnalyzer::beginJob
void beginJob() override
Definition: TrackParameterAnalyzer.cc:65
HLT_FULL_cff.InputTag
InputTag
Definition: HLT_FULL_cff.py:89287
TrackParameterAnalyzer::h1_par4_
TH1 * h1_par4_
Definition: TrackParameterAnalyzer.h:80
TrackParameterAnalyzer::h1_par2_
TH1 * h1_par2_
Definition: TrackParameterAnalyzer.h:78
reco::TrackBase::i_dxy
Definition: TrackBase.h:83
TrackParameterAnalyzer::analyze
void analyze(const edm::Event &, const edm::EventSetup &) override
Definition: TrackParameterAnalyzer.cc:126
reco
fixed size matrix
Definition: AlignmentAlgorithmBase.h:45
TrackParameterAnalyzer::h1_pull0_
TH1 * h1_pull0_
Definition: TrackParameterAnalyzer.h:64
findQualityFiles.v
v
Definition: findQualityFiles.py:179
ParameterVector
reco::TrackBase::ParameterVector ParameterVector
Definition: TrackParameterAnalyzer.h:40
TrackParameterAnalyzer::edmSimVertexContainerToken_
edm::EDGetTokenT< edm::SimVertexContainer > edmSimVertexContainerToken_
Definition: TrackParameterAnalyzer.h:58
edm::Handle< edm::SimVertexContainer >
reco::TrackBase::i_dsz
Definition: TrackBase.h:83
TrackParameterAnalyzer::h1_par0_
TH1 * h1_par0_
Definition: TrackParameterAnalyzer.h:76
class-composition.Q
Q
Definition: class-composition.py:82
TrackParameterAnalyzer::h1_res3_
TH1 * h1_res3_
Definition: TrackParameterAnalyzer.h:72
funct::sin
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
MakerMacros.h
alignCSCRings.s
s
Definition: alignCSCRings.py:92
reco::TrackBase::i_qoverp
Definition: TrackBase.h:83
funct::cos
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
Track.h
TrackParameterAnalyzer::verbose_
bool verbose_
Definition: TrackParameterAnalyzer.h:82
mathSSE::sqrt
T sqrt(T t)
Definition: SSEVec.h:19
TrackParameterAnalyzer::h1_par3_
TH1 * h1_par3_
Definition: TrackParameterAnalyzer.h:79
reco::TrackBase::ParameterVector
math::Vector< dimension >::type ParameterVector
parameter vector
Definition: TrackBase.h:71
HLTMuonOfflineAnalyzer_cfi.z0
z0
Definition: HLTMuonOfflineAnalyzer_cfi.py:98
b
double b
Definition: hdecay.h:118
TrackParameterAnalyzer::h2_dvsphi_
TH2 * h2_dvsphi_
Definition: TrackParameterAnalyzer.h:75
AlCaHLTBitMon_QueryRunRegistry.string
string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
reco::TrackBase::i_lambda
Definition: TrackBase.h:83
TrackParameterAnalyzer::h1_res1_
TH1 * h1_res1_
Definition: TrackParameterAnalyzer.h:70
edm::ParameterSet
Definition: ParameterSet.h:47
GetReleaseVersion.h
TrackParameterAnalyzer::h1_Beff_
TH1 * h1_Beff_
Definition: TrackParameterAnalyzer.h:74
a
double a
Definition: hdecay.h:119
TrackParameterAnalyzer::match
bool match(const ParameterVector &a, const ParameterVector &b)
Definition: TrackParameterAnalyzer.cc:114
Event.h
iEvent
int iEvent
Definition: GenABIO.cc:224
M_PI
#define M_PI
Definition: BXVectorInputProducer.cc:49
submitPVResolutionJobs.q
q
Definition: submitPVResolutionJobs.py:84
edm::EventSetup
Definition: EventSetup.h:57
HltBtagPostValidation_cff.c
c
Definition: HltBtagPostValidation_cff.py:31
TrackParameterAnalyzer::rootFile_
TFile * rootFile_
Definition: TrackParameterAnalyzer.h:63
TrackParameterAnalyzer::endJob
void endJob() override
Definition: TrackParameterAnalyzer.cc:90
TrackParameterAnalyzer::recoTrackCollectionToken_
edm::EDGetTokenT< reco::TrackCollection > recoTrackCollectionToken_
Definition: TrackParameterAnalyzer.h:60
std
Definition: JetResolutionObject.h:76
TrackParameterAnalyzer::TrackParameterAnalyzer
TrackParameterAnalyzer(const edm::ParameterSet &)
Definition: TrackParameterAnalyzer.cc:36
edm::SimTrackContainer
std::vector< SimTrack > SimTrackContainer
Definition: SimTrackContainer.h:12
edm::getReleaseVersion
std::string getReleaseVersion()
Definition: GetReleaseVersion.cc:7
TrackParameterAnalyzer::h1_par1_
TH1 * h1_par1_
Definition: TrackParameterAnalyzer.h:77
TrackParameterAnalyzer::h1_res0_
TH1 * h1_res0_
Definition: TrackParameterAnalyzer.h:69
TrackParameterAnalyzer::outputFile_
std::string outputFile_
Definition: TrackParameterAnalyzer.h:62
funct::pow
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29
kappa
static const G4double kappa
Definition: UrbanMscModel93.cc:35
TrackParameterAnalyzer::h1_pull1_
TH1 * h1_pull1_
Definition: TrackParameterAnalyzer.h:65
reco::TrackBase::CovarianceMatrix
math::Error< dimension >::type CovarianceMatrix
5 parameter covariance matrix
Definition: TrackBase.h:74
edm::Event
Definition: Event.h:73
submitPVValidationJobs.t
string t
Definition: submitPVValidationJobs.py:644
hiRegitInitialStep_cff.recTracks
recTracks
Definition: hiRegitInitialStep_cff.py:14
TrackParameterAnalyzer::h1_pull3_
TH1 * h1_pull3_
Definition: TrackParameterAnalyzer.h:67
TrackParameterAnalyzer::simUnit_
double simUnit_
Definition: TrackParameterAnalyzer.h:81
edm::SimVertexContainer
std::vector< SimVertex > SimVertexContainer
Definition: SimVertexContainer.h:12
reco::TrackCollection
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:14
TrackParameterAnalyzer::h1_res4_
TH1 * h1_res4_
Definition: TrackParameterAnalyzer.h:73