1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
|
#include <memory>
#include <string>
#include <iostream>
#include <fstream>
#include <TMath.h>
#include "FWCore/PluginManager/interface/ModuleDef.h"
#include "FWCore/Framework/interface/MakerMacros.h"
#include "DataFormats/TrajectorySeed/interface/TrajectorySeedCollection.h"
#include "DataFormats/GeometryVector/interface/GlobalVector.h"
#include "DataFormats/GeometryVector/interface/LocalVector.h"
#include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
#include "Geometry/CommonDetUnit/interface/GeomDetType.h"
#include "Geometry/CommonDetUnit/interface/GeomDet.h"
#include "DataFormats/TrackerRecHit2D/interface/SiPixelRecHit.h"
#include "DataFormats/TrackReco/interface/TrackExtra.h"
#include "DataFormats/SiStripDetId/interface/StripSubdetector.h"
#include "DataFormats/TrackerRecHit2D/interface/SiStripMatchedRecHit2D.h"
#include "Geometry/CommonTopologies/interface/StripTopology.h"
#include "DataFormats/TrackerCommon/interface/TrackerTopology.h"
#include "Geometry/Records/interface/TrackerTopologyRcd.h"
#include "DataFormats/TrackReco/interface/TrackFwd.h"
#include "TrackingTools/Records/interface/TransientRecHitRecord.h"
#include "Geometry/CommonDetUnit/interface/GluedGeomDet.h"
#include "TrackingTools/TrackFitters/interface/TrajectoryStateCombiner.h"
#include "TrackingTools/PatternTools/interface/TrajTrackAssociation.h"
#include "TrackingTools/TransientTrack/interface/TransientTrack.h"
#include "SiPixelLorentzAngle.h"
using namespace std;
using namespace edm;
using namespace reco;
using namespace analyzer;
SiPixelLorentzAngle::SiPixelLorentzAngle(edm::ParameterSet const& conf)
: filename_(conf.getParameter<std::string>("fileName")),
filenameFit_(conf.getParameter<std::string>("fileNameFit")),
ptmin_(conf.getParameter<double>("ptMin")),
simData_(conf.getParameter<bool>("simData")),
normChi2Max_(conf.getParameter<double>("normChi2Max")),
clustSizeYMin_(conf.getParameter<int>("clustSizeYMin")),
residualMax_(conf.getParameter<double>("residualMax")),
clustChargeMax_(conf.getParameter<double>("clustChargeMax")),
hist_depth_(conf.getParameter<int>("binsDepth")),
hist_drift_(conf.getParameter<int>("binsDrift")),
trackerHitAssociatorConfig_(consumesCollector()) {
// anglefinder_=new TrackLocalAngle(conf);
hist_x_ = 50;
hist_y_ = 100;
min_x_ = -500.;
max_x_ = 500.;
min_y_ = -1500.;
max_y_ = 500.;
width_ = 0.0285;
min_depth_ = -100.;
max_depth_ = 400.;
min_drift_ = -1000.; //-200.;(conf.getParameter<double>("residualMax"))
max_drift_ = 1000.; //400.;
t_trajTrack = consumes<TrajTrackAssociationCollection>(conf.getParameter<edm::InputTag>("src"));
trackerTopoToken_ = esConsumes<TrackerTopology, TrackerTopologyRcd>();
trackerGeomToken_ = esConsumes<TrackerGeometry, TrackerDigiGeometryRecord>();
}
// Virtual destructor needed.
SiPixelLorentzAngle::~SiPixelLorentzAngle() = default;
void SiPixelLorentzAngle::beginJob() {
// cout << "started SiPixelLorentzAngle" << endl;
hFile_ = new TFile(filename_.c_str(), "RECREATE");
int bufsize = 64000;
// create tree structure
SiPixelLorentzAngleTree_ = new TTree("SiPixelLorentzAngleTree_", "SiPixel LorentzAngle tree", bufsize);
SiPixelLorentzAngleTree_->Branch("run", &run_, "run/I", bufsize);
SiPixelLorentzAngleTree_->Branch("event", &event_, "event/I", bufsize);
SiPixelLorentzAngleTree_->Branch("module", &module_, "module/I", bufsize);
SiPixelLorentzAngleTree_->Branch("ladder", &ladder_, "ladder/I", bufsize);
SiPixelLorentzAngleTree_->Branch("layer", &layer_, "layer/I", bufsize);
SiPixelLorentzAngleTree_->Branch("isflipped", &isflipped_, "isflipped/I", bufsize);
SiPixelLorentzAngleTree_->Branch("pt", &pt_, "pt/F", bufsize);
SiPixelLorentzAngleTree_->Branch("eta", &eta_, "eta/F", bufsize);
SiPixelLorentzAngleTree_->Branch("phi", &phi_, "phi/F", bufsize);
SiPixelLorentzAngleTree_->Branch("chi2", &chi2_, "chi2/D", bufsize);
SiPixelLorentzAngleTree_->Branch("ndof", &ndof_, "ndof/D", bufsize);
SiPixelLorentzAngleTree_->Branch("trackhit", &trackhit_, "x/F:y/F:alpha/D:beta/D:gamma_/D", bufsize);
SiPixelLorentzAngleTree_->Branch("simhit", &simhit_, "x/F:y/F:alpha/D:beta/D:gamma_/D", bufsize);
SiPixelLorentzAngleTree_->Branch("npix", &pixinfo_.npix, "npix/I", bufsize);
SiPixelLorentzAngleTree_->Branch("rowpix", pixinfo_.row, "row[npix]/F", bufsize);
SiPixelLorentzAngleTree_->Branch("colpix", pixinfo_.col, "col[npix]/F", bufsize);
SiPixelLorentzAngleTree_->Branch("adc", pixinfo_.adc, "adc[npix]/F", bufsize);
SiPixelLorentzAngleTree_->Branch("xpix", pixinfo_.x, "x[npix]/F", bufsize);
SiPixelLorentzAngleTree_->Branch("ypix", pixinfo_.y, "y[npix]/F", bufsize);
SiPixelLorentzAngleTree_->Branch(
"clust",
&clust_,
"x/F:y/F:charge/F:size_x/I:size_y/I:maxPixelCol/I:maxPixelRow:minPixelCol/I:minPixelRow/I",
bufsize);
SiPixelLorentzAngleTree_->Branch("rechit", &rechit_, "x/F:y/F", bufsize);
SiPixelLorentzAngleTreeForward_ =
new TTree("SiPixelLorentzAngleTreeForward_", "SiPixel LorentzAngle tree forward", bufsize);
SiPixelLorentzAngleTreeForward_->Branch("run", &run_, "run/I", bufsize);
SiPixelLorentzAngleTreeForward_->Branch("event", &event_, "event/I", bufsize);
SiPixelLorentzAngleTreeForward_->Branch("side", &sideF_, "side/I", bufsize);
SiPixelLorentzAngleTreeForward_->Branch("disk", &diskF_, "disk/I", bufsize);
SiPixelLorentzAngleTreeForward_->Branch("blade", &bladeF_, "blade/I", bufsize);
SiPixelLorentzAngleTreeForward_->Branch("panel", &panelF_, "panel/I", bufsize);
SiPixelLorentzAngleTreeForward_->Branch("module", &moduleF_, "module/I", bufsize);
SiPixelLorentzAngleTreeForward_->Branch("pt", &pt_, "pt/F", bufsize);
SiPixelLorentzAngleTreeForward_->Branch("eta", &eta_, "eta/F", bufsize);
SiPixelLorentzAngleTreeForward_->Branch("phi", &phi_, "phi/F", bufsize);
SiPixelLorentzAngleTreeForward_->Branch("chi2", &chi2_, "chi2/D", bufsize);
SiPixelLorentzAngleTreeForward_->Branch("ndof", &ndof_, "ndof/D", bufsize);
SiPixelLorentzAngleTreeForward_->Branch("trackhit", &trackhitF_, "x/F:y/F:alpha/D:beta/D:gamma_/D", bufsize);
SiPixelLorentzAngleTreeForward_->Branch("simhit", &simhitF_, "x/F:y/F:alpha/D:beta/D:gamma_/D", bufsize);
SiPixelLorentzAngleTreeForward_->Branch("npix", &pixinfoF_.npix, "npix/I", bufsize);
SiPixelLorentzAngleTreeForward_->Branch("rowpix", pixinfoF_.row, "row[npix]/F", bufsize);
SiPixelLorentzAngleTreeForward_->Branch("colpix", pixinfoF_.col, "col[npix]/F", bufsize);
SiPixelLorentzAngleTreeForward_->Branch("adc", pixinfoF_.adc, "adc[npix]/F", bufsize);
SiPixelLorentzAngleTreeForward_->Branch("xpix", pixinfoF_.x, "x[npix]/F", bufsize);
SiPixelLorentzAngleTreeForward_->Branch("ypix", pixinfoF_.y, "y[npix]/F", bufsize);
SiPixelLorentzAngleTreeForward_->Branch(
"clust",
&clustF_,
"x/F:y/F:charge/F:size_x/I:size_y/I:maxPixelCol/I:maxPixelRow:minPixelCol/I:minPixelRow/I",
bufsize);
SiPixelLorentzAngleTreeForward_->Branch("rechit", &rechitF_, "x/F:y/F", bufsize);
//book histograms
char name[128];
for (int i_module = 1; i_module <= 8; i_module++) {
for (int i_layer = 1; i_layer <= 3; i_layer++) {
sprintf(name, "h_drift_depth_adc_layer%i_module%i", i_layer, i_module);
_h_drift_depth_adc_[i_module + (i_layer - 1) * 8] =
new TH2F(name, name, hist_drift_, min_drift_, max_drift_, hist_depth_, min_depth_, max_depth_);
sprintf(name, "h_drift_depth_adc2_layer%i_module%i", i_layer, i_module);
_h_drift_depth_adc2_[i_module + (i_layer - 1) * 8] =
new TH2F(name, name, hist_drift_, min_drift_, max_drift_, hist_depth_, min_depth_, max_depth_);
sprintf(name, "h_drift_depth_noadc_layer%i_module%i", i_layer, i_module);
_h_drift_depth_noadc_[i_module + (i_layer - 1) * 8] =
new TH2F(name, name, hist_drift_, min_drift_, max_drift_, hist_depth_, min_depth_, max_depth_);
sprintf(name, "h_drift_depth_layer%i_module%i", i_layer, i_module);
_h_drift_depth_[i_module + (i_layer - 1) * 8] =
new TH2F(name, name, hist_drift_, min_drift_, max_drift_, hist_depth_, min_depth_, max_depth_);
sprintf(name, "h_mean_layer%i_module%i", i_layer, i_module);
_h_mean_[i_module + (i_layer - 1) * 8] = new TH1F(name, name, hist_depth_, min_depth_, max_depth_);
}
}
// just for some expaining plots
h_cluster_shape_adc_ = new TH2F(
"h_cluster_shape_adc", "cluster shape with adc weight", hist_x_, min_x_, max_x_, hist_y_, min_y_, max_y_);
h_cluster_shape_noadc_ = new TH2F(
"h_cluster_shape_noadc", "cluster shape without adc weight", hist_x_, min_x_, max_x_, hist_y_, min_y_, max_y_);
h_cluster_shape_ = new TH2F("h_cluster_shape", "cluster shape", hist_x_, min_x_, max_x_, hist_y_, min_y_, max_y_);
h_cluster_shape_adc_rot_ = new TH2F(
"h_cluster_shape_adc_rot", "cluster shape with adc weight", hist_x_, min_x_, max_x_, hist_y_, -max_y_, -min_y_);
h_cluster_shape_noadc_rot_ = new TH2F("h_cluster_shape_noadc_rot",
"cluster shape without adc weight",
hist_x_,
min_x_,
max_x_,
hist_y_,
-max_y_,
-min_y_);
h_cluster_shape_rot_ =
new TH2F("h_cluster_shape_rot", "cluster shape", hist_x_, min_x_, max_x_, hist_y_, -max_y_, -min_y_);
h_tracks_ = new TH1F("h_tracks", "h_tracks", 2, 0., 2.);
event_counter_ = 0;
trackEventsCounter_ = 0;
// trackcounter_ = 0;
hitCounter_ = 0;
usedHitCounter_ = 0;
pixelTracksCounter_ = 0;
}
// Functions that gets called by framework every event
void SiPixelLorentzAngle::analyze(const edm::Event& e, const edm::EventSetup& es) {
//Retrieve tracker topology from geometry
edm::ESHandle<TrackerTopology> tTopoHandle = es.getHandle(trackerTopoToken_);
const TrackerTopology* const tTopo = tTopoHandle.product();
event_counter_++;
// if(event_counter_ % 500 == 0) cout << "event number " << event_counter_ << endl;
cout << "event number " << event_counter_ << endl;
edm::ESHandle<TrackerGeometry> estracker = es.getHandle(trackerGeomToken_);
const TrackerGeometry* tracker = &(*estracker);
std::unique_ptr<TrackerHitAssociator> associate;
if (simData_)
associate = std::make_unique<TrackerHitAssociator>(e, trackerHitAssociatorConfig_);
// restet values
module_ = -1;
layer_ = -1;
ladder_ = -1;
isflipped_ = -1;
pt_ = -999;
eta_ = 999;
phi_ = 999;
pixinfo_.npix = 0;
run_ = e.id().run();
event_ = e.id().event();
// get the association map between tracks and trajectories
edm::Handle<TrajTrackAssociationCollection> trajTrackCollectionHandle;
e.getByToken(t_trajTrack, trajTrackCollectionHandle);
if (!trajTrackCollectionHandle->empty()) {
trackEventsCounter_++;
for (TrajTrackAssociationCollection::const_iterator it = trajTrackCollectionHandle->begin();
it != trajTrackCollectionHandle->end();
++it) {
const Track& track = *it->val;
const Trajectory& traj = *it->key;
// get the trajectory measurements
std::vector<TrajectoryMeasurement> tmColl = traj.measurements();
// TrajectoryStateOnSurface tsos = tsoscomb( itTraj->forwardPredictedState(), itTraj->backwardPredictedState() );
pt_ = track.pt();
eta_ = track.eta();
phi_ = track.phi();
chi2_ = traj.chiSquared();
ndof_ = traj.ndof();
if (pt_ < ptmin_)
continue;
// iterate over trajectory measurements
std::vector<PSimHit> matched;
h_tracks_->Fill(0);
bool pixeltrack = false;
for (std::vector<TrajectoryMeasurement>::const_iterator itTraj = tmColl.begin(); itTraj != tmColl.end();
itTraj++) {
if (!itTraj->updatedState().isValid())
continue;
TransientTrackingRecHit::ConstRecHitPointer recHit = itTraj->recHit();
if (!recHit->isValid() || recHit->geographicalId().det() != DetId::Tracker)
continue;
unsigned int subDetID = (recHit->geographicalId().subdetId());
if (subDetID == PixelSubdetector::PixelBarrel || subDetID == PixelSubdetector::PixelEndcap) {
if (!pixeltrack) {
h_tracks_->Fill(1);
pixelTracksCounter_++;
}
pixeltrack = true;
}
if (subDetID == PixelSubdetector::PixelBarrel) {
hitCounter_++;
DetId detIdObj = recHit->geographicalId();
const PixelGeomDetUnit* theGeomDet = dynamic_cast<const PixelGeomDetUnit*>(tracker->idToDet(detIdObj));
if (!theGeomDet)
continue;
const PixelTopology* topol = &(theGeomDet->specificTopology());
if (!topol)
continue;
layer_ = tTopo->pxbLayer(detIdObj);
ladder_ = tTopo->pxbLadder(detIdObj);
module_ = tTopo->pxbModule(detIdObj);
float tmp1 = theGeomDet->surface().toGlobal(Local3DPoint(0., 0., 0.)).perp();
float tmp2 = theGeomDet->surface().toGlobal(Local3DPoint(0., 0., 1.)).perp();
if (tmp2 < tmp1)
isflipped_ = 1;
else
isflipped_ = 0;
const SiPixelRecHit* recHitPix = dynamic_cast<const SiPixelRecHit*>((*recHit).hit());
if (!recHitPix)
continue;
rechit_.x = recHitPix->localPosition().x();
rechit_.y = recHitPix->localPosition().y();
SiPixelRecHit::ClusterRef const& cluster = recHitPix->cluster();
// fill entries in clust_
clust_.x = (cluster)->x();
clust_.y = (cluster)->y();
clust_.charge = (cluster->charge()) / 1000.;
clust_.size_x = cluster->sizeX();
clust_.size_y = cluster->sizeY();
clust_.maxPixelCol = cluster->maxPixelCol();
clust_.maxPixelRow = cluster->maxPixelRow();
clust_.minPixelCol = cluster->minPixelCol();
clust_.minPixelRow = cluster->minPixelRow();
// fill entries in pixinfo_:
fillPix(*cluster, topol, pixinfo_);
// fill the trackhit info
TrajectoryStateOnSurface tsos = itTraj->updatedState();
if (!tsos.isValid()) {
cout << "tsos not valid" << endl;
continue;
}
LocalVector trackdirection = tsos.localDirection();
LocalPoint trackposition = tsos.localPosition();
if (trackdirection.z() == 0)
continue;
// the local position and direction
trackhit_.alpha = atan2(trackdirection.z(), trackdirection.x());
trackhit_.beta = atan2(trackdirection.z(), trackdirection.y());
trackhit_.gamma = atan2(trackdirection.x(), trackdirection.y());
trackhit_.x = trackposition.x();
trackhit_.y = trackposition.y();
// fill entries in simhit_:
if (simData_) {
matched.clear();
matched = associate->associateHit((*recHitPix));
float dr_start = 9999.;
for (std::vector<PSimHit>::iterator isim = matched.begin(); isim != matched.end(); ++isim) {
DetId simdetIdObj((*isim).detUnitId());
if (simdetIdObj == detIdObj) {
float sim_x1 = (*isim).entryPoint().x(); // width (row index, in col direction)
float sim_y1 = (*isim).entryPoint().y(); // length (col index, in row direction)
float sim_x2 = (*isim).exitPoint().x();
float sim_y2 = (*isim).exitPoint().y();
float sim_xpos = 0.5 * (sim_x1 + sim_x2);
float sim_ypos = 0.5 * (sim_y1 + sim_y2);
float sim_px = (*isim).momentumAtEntry().x();
float sim_py = (*isim).momentumAtEntry().y();
float sim_pz = (*isim).momentumAtEntry().z();
float dr = (sim_xpos - (recHitPix->localPosition().x())) * (sim_xpos - recHitPix->localPosition().x()) +
(sim_ypos - recHitPix->localPosition().y()) * (sim_ypos - recHitPix->localPosition().y());
if (dr < dr_start) {
simhit_.x = sim_xpos;
simhit_.y = sim_ypos;
simhit_.alpha = atan2(sim_pz, sim_px);
simhit_.beta = atan2(sim_pz, sim_py);
simhit_.gamma = atan2(sim_px, sim_py);
dr_start = dr;
}
}
} // end of filling simhit_
}
// is one pixel in cluster a large pixel ? (hit will be excluded)
bool large_pix = false;
for (int j = 0; j < pixinfo_.npix; j++) {
int colpos = static_cast<int>(pixinfo_.col[j]);
if (pixinfo_.row[j] == 0 || pixinfo_.row[j] == 79 || pixinfo_.row[j] == 80 || pixinfo_.row[j] == 159 ||
colpos % 52 == 0 || colpos % 52 == 51) {
large_pix = true;
}
}
double residual = TMath::Sqrt((trackhit_.x - rechit_.x) * (trackhit_.x - rechit_.x) +
(trackhit_.y - rechit_.y) * (trackhit_.y - rechit_.y));
SiPixelLorentzAngleTree_->Fill();
if (!large_pix && (chi2_ / ndof_) < normChi2Max_ && cluster->sizeY() >= clustSizeYMin_ &&
residual < residualMax_ && (cluster->charge() < clustChargeMax_)) {
usedHitCounter_++;
// iterate over pixels in hit
for (int j = 0; j < pixinfo_.npix; j++) {
// use trackhits
float dx = (pixinfo_.x[j] - (trackhit_.x - width_ / 2. / TMath::Tan(trackhit_.alpha))) * 10000.;
float dy = (pixinfo_.y[j] - (trackhit_.y - width_ / 2. / TMath::Tan(trackhit_.beta))) * 10000.;
float depth = dy * tan(trackhit_.beta);
float drift = dx - dy * tan(trackhit_.gamma);
_h_drift_depth_adc_[module_ + (layer_ - 1) * 8]->Fill(drift, depth, pixinfo_.adc[j]);
_h_drift_depth_adc2_[module_ + (layer_ - 1) * 8]->Fill(drift, depth, pixinfo_.adc[j] * pixinfo_.adc[j]);
_h_drift_depth_noadc_[module_ + (layer_ - 1) * 8]->Fill(drift, depth);
if (layer_ == 3 && module_ == 1 && isflipped_) {
float dx_rot = dx * TMath::Cos(trackhit_.gamma) + dy * TMath::Sin(trackhit_.gamma);
float dy_rot = dy * TMath::Cos(trackhit_.gamma) - dx * TMath::Sin(trackhit_.gamma);
h_cluster_shape_adc_->Fill(dx, dy, pixinfo_.adc[j]);
h_cluster_shape_noadc_->Fill(dx, dy);
h_cluster_shape_adc_rot_->Fill(dx_rot, dy_rot, pixinfo_.adc[j]);
h_cluster_shape_noadc_rot_->Fill(dx_rot, dy_rot);
}
}
}
} else if (subDetID == PixelSubdetector::PixelEndcap) {
DetId detIdObj = recHit->geographicalId();
const PixelGeomDetUnit* theGeomDet = dynamic_cast<const PixelGeomDetUnit*>(tracker->idToDet(detIdObj));
if (!theGeomDet)
continue;
const PixelTopology* topol = &(theGeomDet->specificTopology());
if (!topol)
continue;
sideF_ = tTopo->pxfSide(detIdObj);
diskF_ = tTopo->pxfDisk(detIdObj);
bladeF_ = tTopo->pxfBlade(detIdObj);
panelF_ = tTopo->pxfPanel(detIdObj);
moduleF_ = tTopo->pxfModule(detIdObj);
//float tmp1 = theGeomDet->surface().toGlobal(Local3DPoint(0.,0.,0.)).perp();
//float tmp2 = theGeomDet->surface().toGlobal(Local3DPoint(0.,0.,1.)).perp();
//if ( tmp2<tmp1 ) isflipped_ = 1;
//else isflipped_ = 0;
const SiPixelRecHit* recHitPix = dynamic_cast<const SiPixelRecHit*>((*recHit).hit());
if (!recHitPix)
continue;
rechitF_.x = recHitPix->localPosition().x();
rechitF_.y = recHitPix->localPosition().y();
SiPixelRecHit::ClusterRef const& cluster = recHitPix->cluster();
// fill entries in clust_
clustF_.x = (cluster)->x();
clustF_.y = (cluster)->y();
clustF_.charge = (cluster->charge()) / 1000.;
clustF_.size_x = cluster->sizeX();
clustF_.size_y = cluster->sizeY();
clustF_.maxPixelCol = cluster->maxPixelCol();
clustF_.maxPixelRow = cluster->maxPixelRow();
clustF_.minPixelCol = cluster->minPixelCol();
clustF_.minPixelRow = cluster->minPixelRow();
// fill entries in pixinfo_:
fillPix(*cluster, topol, pixinfoF_);
// fill the trackhit info
TrajectoryStateOnSurface tsos = itTraj->updatedState();
if (!tsos.isValid()) {
cout << "tsos not valid" << endl;
continue;
}
LocalVector trackdirection = tsos.localDirection();
LocalPoint trackposition = tsos.localPosition();
if (trackdirection.z() == 0)
continue;
// the local position and direction
trackhitF_.alpha = atan2(trackdirection.z(), trackdirection.x());
trackhitF_.beta = atan2(trackdirection.z(), trackdirection.y());
trackhitF_.gamma = atan2(trackdirection.x(), trackdirection.y());
trackhitF_.x = trackposition.x();
trackhitF_.y = trackposition.y();
// fill entries in simhit_:
if (simData_) {
matched.clear();
matched = associate->associateHit((*recHitPix));
float dr_start = 9999.;
for (std::vector<PSimHit>::iterator isim = matched.begin(); isim != matched.end(); ++isim) {
DetId simdetIdObj((*isim).detUnitId());
if (simdetIdObj == detIdObj) {
float sim_x1 = (*isim).entryPoint().x(); // width (row index, in col direction)
float sim_y1 = (*isim).entryPoint().y(); // length (col index, in row direction)
float sim_x2 = (*isim).exitPoint().x();
float sim_y2 = (*isim).exitPoint().y();
float sim_xpos = 0.5 * (sim_x1 + sim_x2);
float sim_ypos = 0.5 * (sim_y1 + sim_y2);
float sim_px = (*isim).momentumAtEntry().x();
float sim_py = (*isim).momentumAtEntry().y();
float sim_pz = (*isim).momentumAtEntry().z();
float dr = (sim_xpos - (recHitPix->localPosition().x())) * (sim_xpos - recHitPix->localPosition().x()) +
(sim_ypos - recHitPix->localPosition().y()) * (sim_ypos - recHitPix->localPosition().y());
if (dr < dr_start) {
simhitF_.x = sim_xpos;
simhitF_.y = sim_ypos;
simhitF_.alpha = atan2(sim_pz, sim_px);
simhitF_.beta = atan2(sim_pz, sim_py);
simhitF_.gamma = atan2(sim_px, sim_py);
dr_start = dr;
}
}
} // end of filling simhit_
}
SiPixelLorentzAngleTreeForward_->Fill();
}
} //end iteration over trajectory measurements
} //end iteration over trajectories
}
}
void SiPixelLorentzAngle::endJob() {
// produce histograms with the average adc counts
for (int i_ring = 1; i_ring <= 24; i_ring++) {
_h_drift_depth_[i_ring]->Divide(_h_drift_depth_adc_[i_ring], _h_drift_depth_noadc_[i_ring]);
}
h_drift_depth_adc_slice_ =
new TH1F("h_drift_depth_adc_slice", "slice of adc histogram", hist_drift_, min_drift_, max_drift_);
TF1* f1 = new TF1("f1", "[0] + [1]*x", 50., 235.);
f1->SetParName(0, "p0");
f1->SetParName(1, "p1");
f1->SetParameter(0, 0);
f1->SetParameter(1, 0.4);
ofstream fLorentzFit(filenameFit_.c_str(), ios::trunc);
cout.precision(4);
fLorentzFit << "module"
<< "\t"
<< "layer"
<< "\t"
<< "offset"
<< "\t"
<< "error"
<< "\t"
<< "slope"
<< "\t"
<< "error"
<< "\t"
"rel.err"
<< "\t"
"pull"
<< "\t"
<< "chi2"
<< "\t"
<< "prob" << endl;
//loop over modlues and layers to fit the lorentz angle
for (int i_layer = 1; i_layer <= 3; i_layer++) {
for (int i_module = 1; i_module <= 8; i_module++) {
//loop over bins in depth (z-local-coordinate) (in order to fit slices)
for (int i = 1; i <= hist_depth_; i++) {
findMean(i, (i_module + (i_layer - 1) * 8));
} // end loop over bins in depth
_h_mean_[i_module + (i_layer - 1) * 8]->Fit(f1, "ERQ");
double p0 = f1->GetParameter(0);
double e0 = f1->GetParError(0);
double p1 = f1->GetParameter(1);
double e1 = f1->GetParError(1);
double chi2 = f1->GetChisquare();
double prob = f1->GetProb();
fLorentzFit << setprecision(4) << i_module << "\t" << i_layer << "\t" << p0 << "\t" << e0 << "\t" << p1
<< setprecision(3) << "\t" << e1 << "\t" << e1 / p1 * 100. << "\t" << (p1 - 0.424) / e1 << "\t"
<< chi2 << "\t" << prob << endl;
}
} // end loop over modules and layers
fLorentzFit.close();
hFile_->cd();
for (int i_module = 1; i_module <= 8; i_module++) {
for (int i_layer = 1; i_layer <= 3; i_layer++) {
_h_drift_depth_adc_[i_module + (i_layer - 1) * 8]->Write();
_h_drift_depth_adc2_[i_module + (i_layer - 1) * 8]->Write();
_h_drift_depth_noadc_[i_module + (i_layer - 1) * 8]->Write();
_h_drift_depth_[i_module + (i_layer - 1) * 8]->Write();
_h_mean_[i_module + (i_layer - 1) * 8]->Write();
}
}
h_cluster_shape_adc_->Write();
h_cluster_shape_noadc_->Write();
h_cluster_shape_adc_rot_->Write();
h_cluster_shape_noadc_rot_->Write();
h_tracks_->Write();
hFile_->Write();
hFile_->Close();
cout << "events: " << event_counter_ << endl;
cout << "events with tracks: " << trackEventsCounter_ << endl;
cout << "events with pixeltracks: " << pixelTracksCounter_ << endl;
cout << "hits in the pixel: " << hitCounter_ << endl;
cout << "number of used Hits: " << usedHitCounter_ << endl;
}
inline void SiPixelLorentzAngle::fillPix(const SiPixelCluster& LocPix, const PixelTopology* topol, Pixinfo& pixinfo)
{
const std::vector<SiPixelCluster::Pixel>& pixvector = LocPix.pixels();
pixinfo.npix = 0;
for (std::vector<SiPixelCluster::Pixel>::const_iterator itPix = pixvector.begin(); itPix != pixvector.end();
itPix++) {
// for(pixinfo.npix = 0; pixinfo.npix < static_cast<int>(pixvector.size()); ++pixinfo.npix) {
pixinfo.row[pixinfo.npix] = itPix->x;
pixinfo.col[pixinfo.npix] = itPix->y;
pixinfo.adc[pixinfo.npix] = itPix->adc;
LocalPoint lp = topol->localPosition(MeasurementPoint(itPix->x + 0.5, itPix->y + 0.5));
pixinfo.x[pixinfo.npix] = lp.x();
pixinfo.y[pixinfo.npix] = lp.y();
pixinfo.npix++;
}
}
void SiPixelLorentzAngle::findMean(int i, int i_ring) {
double nentries = 0;
h_drift_depth_adc_slice_->Reset("ICE");
// determine sigma and sigma^2 of the adc counts and average adc counts
//loop over bins in drift width
for (int j = 1; j <= hist_drift_; j++) {
if (_h_drift_depth_noadc_[i_ring]->GetBinContent(j, i) >= 1) {
double adc_error2 =
(_h_drift_depth_adc2_[i_ring]->GetBinContent(j, i) - _h_drift_depth_adc_[i_ring]->GetBinContent(j, i) *
_h_drift_depth_adc_[i_ring]->GetBinContent(j, i) /
_h_drift_depth_noadc_[i_ring]->GetBinContent(j, i)) /
_h_drift_depth_noadc_[i_ring]->GetBinContent(j, i);
_h_drift_depth_adc_[i_ring]->SetBinError(j, i, sqrt(adc_error2));
double error2 = adc_error2 / (_h_drift_depth_noadc_[i_ring]->GetBinContent(j, i) - 1.);
_h_drift_depth_[i_ring]->SetBinError(j, i, sqrt(error2));
} else {
_h_drift_depth_[i_ring]->SetBinError(j, i, 0);
_h_drift_depth_adc_[i_ring]->SetBinError(j, i, 0);
}
h_drift_depth_adc_slice_->SetBinContent(j, _h_drift_depth_adc_[i_ring]->GetBinContent(j, i));
h_drift_depth_adc_slice_->SetBinError(j, _h_drift_depth_adc_[i_ring]->GetBinError(j, i));
nentries += _h_drift_depth_noadc_[i_ring]->GetBinContent(j, i);
} // end loop over bins in drift width
double mean = h_drift_depth_adc_slice_->GetMean(1);
double error = 0;
if (nentries != 0) {
error = h_drift_depth_adc_slice_->GetRMS(1) / sqrt(nentries);
}
_h_mean_[i_ring]->SetBinContent(i, mean);
_h_mean_[i_ring]->SetBinError(i, error);
}
DEFINE_FWK_MODULE(SiPixelLorentzAngle);
|