CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
SiPixelTrackResidualSource.cc
Go to the documentation of this file.
1 // Package: SiPixelMonitorTrack
2 // Class: SiPixelTrackResidualSource
3 //
4 // class SiPixelTrackResidualSource SiPixelTrackResidualSource.cc
5 // DQM/SiPixelMonitorTrack/src/SiPixelTrackResidualSource.cc
6 //
7 // Description: SiPixel hit-to-track residual data quality monitoring modules
8 // Implementation: prototype -> improved -> never final - end of the 1st step
9 //
10 // Original Author: Shan-Huei Chuang
11 // Created: Fri Mar 23 18:41:42 CET 2007
12 // Updated by Lukas Wehrli (plots for clusters on/off track added)
13 
14 
15 #include <iostream>
16 #include <map>
17 #include <string>
18 #include <vector>
19 #include <utility>
20 
23 
26 
30 //#include "DataFormats/SiPixelDetId/interface/PixelBarrelNameWrapper.h"
34 
37 
39 
45 
49 
50 //Claudia new libraries
60 
61 using namespace std;
62 using namespace edm;
63 
64 
66  pSet_(pSet),
67  modOn( pSet.getUntrackedParameter<bool>("modOn",true) ),
68  reducedSet( pSet.getUntrackedParameter<bool>("reducedSet",true) ),
69  ladOn( pSet.getUntrackedParameter<bool>("ladOn",false) ),
70  layOn( pSet.getUntrackedParameter<bool>("layOn",false) ),
71  phiOn( pSet.getUntrackedParameter<bool>("phiOn",false) ),
72  ringOn( pSet.getUntrackedParameter<bool>("ringOn",false) ),
73  bladeOn( pSet.getUntrackedParameter<bool>("bladeOn",false) ),
74  diskOn( pSet.getUntrackedParameter<bool>("diskOn",false) ),
75  isUpgrade( pSet.getUntrackedParameter<bool>("isUpgrade",false) ),
76  noOfLayers(0),
77  noOfDisks(0)
78 {
79  pSet_ = pSet;
80  debug_ = pSet_.getUntrackedParameter<bool>("debug", false);
83  tracksrc_ = pSet_.getParameter<edm::InputTag>("trajectoryInput");
84  ttrhbuilder_ = pSet_.getParameter<std::string>("TTRHBuilder");
85  ptminres_= pSet.getUntrackedParameter<double>("PtMinRes",4.0) ;
86  beamSpotToken_ = consumes<reco::BeamSpot>(std::string("offlineBeamSpot"));
87  vtxsrc_=pSet_.getUntrackedParameter<std::string>("vtxsrc", "offlinePrimaryVertices");
88  offlinePrimaryVerticesToken_ = consumes<reco::VertexCollection>(vtxsrc_);// consumes<reco::VertexCollection>(std::string("hiSelectedVertex")); //"offlinePrimaryVertices"));
89  generalTracksToken_ = consumes<reco::TrackCollection>(pSet_.getParameter<edm::InputTag>("tracksrc"));
90  tracksrcToken_ = consumes<std::vector<Trajectory> >(pSet_.getParameter<edm::InputTag>("trajectoryInput"));
91  trackToken_ = consumes<std::vector<reco::Track> >(pSet_.getParameter<edm::InputTag>("trajectoryInput"));
92  trackAssociationToken_ = consumes<TrajTrackAssociationCollection>(pSet_.getParameter<edm::InputTag>("trajectoryInput"));
93  clustersrcToken_ = consumes<edmNew::DetSetVector<SiPixelCluster> >(pSet_.getParameter<edm::InputTag>("clustersrc"));
94 
95  LogInfo("PixelDQM") << "SiPixelTrackResidualSource constructor" << endl;
96  LogInfo ("PixelDQM") << "Mod/Lad/Lay/Phi " << modOn << "/" << ladOn << "/"
97  << layOn << "/" << phiOn << std::endl;
98  LogInfo ("PixelDQM") << "Blade/Disk/Ring" << bladeOn << "/" << diskOn << "/"
99  << ringOn << std::endl;
100 
101  topFolderName_ = pSet_.getParameter<std::string>("TopFolderName");
102 
103  firstRun = true;
104  NTotal=0;
105  NLowProb=0;
106 }
107 
108 
110  LogInfo("PixelDQM") << "SiPixelTrackResidualSource destructor" << endl;
111 
112  std::map<uint32_t,SiPixelTrackResidualModule*>::iterator struct_iter;
113  for (struct_iter = theSiPixelStructure.begin() ; struct_iter != theSiPixelStructure.end() ; struct_iter++){
114  delete struct_iter->second;
115  struct_iter->second = 0;
116  }
117 }
118 
119 
121  LogInfo("PixelDQM") << "SiPixelTrackResidualSource beginRun()" << endl;
122  // retrieve TrackerGeometry for pixel dets
124  iSetup.get<TrackerDigiGeometryRecord>().get(TG);
125  if (debug_) LogVerbatim("PixelDQM") << "TrackerGeometry "<< &(*TG) <<" size is "<< TG->dets().size() << endl;
126  edm::ESHandle<TrackerTopology> tTopoHandle;
127  iSetup.get<TrackerTopologyRcd>().get(tTopoHandle);
128  const TrackerTopology *pTT = tTopoHandle.product();
129 
130  // build theSiPixelStructure with the pixel barrel and endcap dets from TrackerGeometry
131  for (TrackerGeometry::DetContainer::const_iterator pxb = TG->detsPXB().begin();
132  pxb!=TG->detsPXB().end(); pxb++) {
133  if (dynamic_cast<PixelGeomDetUnit const *>((*pxb))!=0) {
134  SiPixelTrackResidualModule* module = new SiPixelTrackResidualModule((*pxb)->geographicalId().rawId());
135  theSiPixelStructure.insert(pair<uint32_t, SiPixelTrackResidualModule*>((*pxb)->geographicalId().rawId(), module));
136  //int DBlayer = PixelBarrelNameWrapper(pSet_, DetId((*pxb)->geographicalId())).layerName();
137  int DBlayer = PixelBarrelName(DetId((*pxb)->geographicalId()),pTT,isUpgrade).layerName();
138  if (noOfLayers < DBlayer) noOfLayers = DBlayer;
139  }
140  }
141  for (TrackerGeometry::DetContainer::const_iterator pxf = TG->detsPXF().begin();
142  pxf!=TG->detsPXF().end(); pxf++) {
143  if (dynamic_cast<PixelGeomDetUnit const *>((*pxf))!=0) {
144  SiPixelTrackResidualModule* module = new SiPixelTrackResidualModule((*pxf)->geographicalId().rawId());
145  theSiPixelStructure.insert(pair<uint32_t, SiPixelTrackResidualModule*>((*pxf)->geographicalId().rawId(), module));
146  int DBdisk;
147  DBdisk = PixelEndcapName(DetId((*pxf)->geographicalId()),pTT,isUpgrade).diskName();
148  if (noOfDisks < DBdisk) noOfDisks = DBdisk;
149  }
150  }
151  LogInfo("PixelDQM") << "SiPixelStructure size is " << theSiPixelStructure.size() << endl;
152 }
153 
155 
156  // book residual histograms in theSiPixelFolder - one (x,y) pair of histograms per det
157  SiPixelFolderOrganizer theSiPixelFolder(false);
158  std::stringstream nameX, titleX, nameY, titleY;
159 
160  if(ladOn){
161  iBooker.setCurrentFolder(topFolderName_+"/Barrel");
162  for (int i = 1; i <= noOfLayers; i++){
163  nameX.str(std::string()); nameX <<"siPixelTrackResidualsX_SummedLayer_" << i;
164  titleX.str(std::string()); titleX <<"Layer"<< i << "Hit-to-Track Residual in r-phi";
165  meResidualXSummedLay.push_back(iBooker.book1D(nameX.str(),titleX.str(),100,-150,150));
166  meResidualXSummedLay.at(i-1)->setAxisTitle("hit-to-track residual in r-phi (um)",1);
167  nameY.str(std::string()); nameY <<"siPixelTrackResidualsY_SummedLayer_" << i;
168  titleY.str(std::string()); titleY <<"Layer"<< i << "Hit-to-Track Residual in Z";
169  meResidualYSummedLay.push_back(iBooker.book1D(nameY.str(),titleY.str(),100,-300,300));
170  meResidualYSummedLay.at(i-1)->setAxisTitle("hit-to-track residual in z (um)",1);
171  }
172  }
173 
174  for (std::map<uint32_t, SiPixelTrackResidualModule*>::iterator pxd = theSiPixelStructure.begin();
175  pxd!=theSiPixelStructure.end(); pxd++){
176 
177  if(modOn){
178  if (theSiPixelFolder.setModuleFolder(iBooker,(*pxd).first,0,isUpgrade)) (*pxd).second->book(pSet_,iSetup,iBooker,reducedSet,0,isUpgrade);
179  else throw cms::Exception("LogicError") << "SiPixelTrackResidualSource Folder Creation Failed! ";
180  }
181  if(ladOn){
182  if (theSiPixelFolder.setModuleFolder(iBooker,(*pxd).first,1,isUpgrade)) {
183 
184  (*pxd).second->book(pSet_,iSetup,iBooker,reducedSet,1,isUpgrade);
185  }
186  else throw cms::Exception("LogicError") << "SiPixelTrackResidualSource ladder Folder Creation Failed! ";
187  }
188  if(layOn){
189  if (theSiPixelFolder.setModuleFolder(iBooker,(*pxd).first,2,isUpgrade)) (*pxd).second->book(pSet_,iSetup,iBooker,reducedSet,2,isUpgrade);
190  else throw cms::Exception("LogicError") << "SiPixelTrackResidualSource layer Folder Creation Failed! ";
191  }
192  if(phiOn){
193  if (theSiPixelFolder.setModuleFolder(iBooker,(*pxd).first,3,isUpgrade)) (*pxd).second->book(pSet_,iSetup,iBooker,reducedSet,3,isUpgrade);
194  else throw cms::Exception("LogicError") << "SiPixelTrackResidualSource phi Folder Creation Failed! ";
195  }
196  if(bladeOn){
197  if (theSiPixelFolder.setModuleFolder(iBooker,(*pxd).first,4,isUpgrade)) (*pxd).second->book(pSet_,iSetup,iBooker,reducedSet,4,isUpgrade);
198  else throw cms::Exception("LogicError") << "SiPixelTrackResidualSource Blade Folder Creation Failed! ";
199  }
200  if(diskOn){
201  if (theSiPixelFolder.setModuleFolder(iBooker,(*pxd).first,5,isUpgrade)) (*pxd).second->book(pSet_,iSetup,iBooker,reducedSet,5,isUpgrade);
202  else throw cms::Exception("LogicError") << "SiPixelTrackResidualSource Disk Folder Creation Failed! ";
203  }
204  if(ringOn){
205  if (theSiPixelFolder.setModuleFolder(iBooker,(*pxd).first,6,isUpgrade)) (*pxd).second->book(pSet_,iSetup,iBooker,reducedSet,6,isUpgrade);
206  else throw cms::Exception("LogicError") << "SiPixelTrackResidualSource Ring Folder Creation Failed! ";
207  }
208  }
209 
210 // edm::InputTag tracksrc = pSet_.getParameter<edm::InputTag>("trajectoryInput");
211 // edm::InputTag clustersrc = pSet_.getParameter<edm::InputTag>("clustersrc");
212 
213  //number of tracks
214  iBooker.setCurrentFolder(topFolderName_+"/Tracks");
215  meNofTracks_ = iBooker.book1D("ntracks_" + tracksrc_.label(),"Number of Tracks",4,0,4);
216  meNofTracks_->setAxisTitle("Number of Tracks",1);
217  meNofTracks_->setBinLabel(1,"All");
218  meNofTracks_->setBinLabel(2,"Pixel");
219  meNofTracks_->setBinLabel(3,"BPix");
220  meNofTracks_->setBinLabel(4,"FPix");
221 
222  //number of tracks in pixel fiducial volume
223  iBooker.setCurrentFolder(topFolderName_+"/Tracks");
224  meNofTracksInPixVol_ = iBooker.book1D("ntracksInPixVol_" + tracksrc_.label(),"Number of Tracks crossing Pixel fiducial Volume",2,0,2);
225  meNofTracksInPixVol_->setAxisTitle("Number of Tracks",1);
226  meNofTracksInPixVol_->setBinLabel(1,"With Hits");
227  meNofTracksInPixVol_->setBinLabel(2,"Without Hits");
228 
229  //number of clusters (associated to track / not associated)
230  iBooker.setCurrentFolder(topFolderName_+"/Clusters/OnTrack");
231  meNofClustersOnTrack_ = iBooker.book1D("nclusters_" + clustersrc_.label() + "_tot","Number of Clusters (on track)",3,0,3);
232  meNofClustersOnTrack_->setAxisTitle("Number of Clusters on Track",1);
236  iBooker.setCurrentFolder(topFolderName_+"/Clusters/OffTrack");
237  meNofClustersNotOnTrack_ = iBooker.book1D("nclusters_" + clustersrc_.label() + "_tot","Number of Clusters (off track)",3,0,3);
238  meNofClustersNotOnTrack_->setAxisTitle("Number of Clusters off Track",1);
242 
243  //cluster charge and size
244  //charge
245  //on track
246  iBooker.setCurrentFolder(topFolderName_+"/Clusters/OnTrack");
247  std::stringstream ss1, ss2;
248  meClChargeOnTrack_all = iBooker.book1D("charge_" + clustersrc_.label(),"Charge (on track)",500,0.,500.);
249  meClChargeOnTrack_all->setAxisTitle("Charge size (in ke)",1);
250  meClChargeOnTrack_bpix = iBooker.book1D("charge_" + clustersrc_.label() + "_Barrel","Charge (on track, barrel)",500,0.,500.);
251  meClChargeOnTrack_bpix->setAxisTitle("Charge size (in ke)",1);
252  meClChargeOnTrack_fpix = iBooker.book1D("charge_" + clustersrc_.label() + "_Endcap","Charge (on track, endcap)",500,0.,500.);
253  meClChargeOnTrack_fpix->setAxisTitle("Charge size (in ke)",1);
254  for (int i = 1; i <= noOfLayers; i++)
255  {
256  ss1.str(std::string()); ss1 << "charge_" + clustersrc_.label() + "_Layer_" << i;
257  ss2.str(std::string()); ss2 << "Charge (on track, layer" << i << ")";
258  meClChargeOnTrack_layers.push_back(iBooker.book1D(ss1.str(),ss2.str(),500,0.,500.));
259  meClChargeOnTrack_layers.at(i-1)->setAxisTitle("Charge size (in ke)",1);
260  }
261  for (int i = 1; i <= noOfDisks; i++)
262  {
263  ss1.str(std::string()); ss1 << "charge_" + clustersrc_.label() + "_Disk_p" << i;
264  ss2.str(std::string()); ss2 << "Charge (on track, diskp" << i << ")";
265  meClChargeOnTrack_diskps.push_back(iBooker.book1D(ss1.str(),ss2.str(),500,0.,500.));
266  meClChargeOnTrack_diskps.at(i-1)->setAxisTitle("Charge size (in ke)",1);
267  }
268  for (int i = 1; i <= noOfDisks; i++)
269  {
270  ss1.str(std::string()); ss1 << "charge_" + clustersrc_.label() + "_Disk_m" << i;
271  ss2.str(std::string()); ss2 << "Charge (on track, diskm" << i << ")";
272  meClChargeOnTrack_diskms.push_back(iBooker.book1D(ss1.str(),ss2.str(),500,0.,500.));
273  meClChargeOnTrack_diskms.at(i-1)->setAxisTitle("Charge size (in ke)",1);
274  }
275  //off track
276  iBooker.setCurrentFolder(topFolderName_+"/Clusters/OffTrack");
277  meClChargeNotOnTrack_all = iBooker.book1D("charge_" + clustersrc_.label(),"Charge (off track)",500,0.,500.);
278  meClChargeNotOnTrack_all->setAxisTitle("Charge size (in ke)",1);
279  meClChargeNotOnTrack_bpix = iBooker.book1D("charge_" + clustersrc_.label() + "_Barrel","Charge (off track, barrel)",500,0.,500.);
280  meClChargeNotOnTrack_bpix->setAxisTitle("Charge size (in ke)",1);
281  meClChargeNotOnTrack_fpix = iBooker.book1D("charge_" + clustersrc_.label() + "_Endcap","Charge (off track, endcap)",500,0.,500.);
282  meClChargeNotOnTrack_fpix->setAxisTitle("Charge size (in ke)",1);
283  for (int i = 1; i <= noOfLayers; i++)
284  {
285  ss1.str(std::string()); ss1 << "charge_" + clustersrc_.label() + "_Layer_" << i;
286  ss2.str(std::string()); ss2 << "Charge (off track, layer" << i << ")";
287  meClChargeNotOnTrack_layers.push_back(iBooker.book1D(ss1.str(),ss2.str(),500,0.,500.));
288  meClChargeNotOnTrack_layers.at(i-1)->setAxisTitle("Charge size (in ke)",1);
289  }
290  for (int i = 1; i <= noOfDisks; i++)
291  {
292  ss1.str(std::string()); ss1 << "charge_" + clustersrc_.label() + "_Disk_p" << i;
293  ss2.str(std::string()); ss2 << "Charge (off track, diskp" << i << ")";
294  meClChargeNotOnTrack_diskps.push_back(iBooker.book1D(ss1.str(),ss2.str(),500,0.,500.));
295  meClChargeNotOnTrack_diskps.at(i-1)->setAxisTitle("Charge size (in ke)",1);
296  }
297  for (int i = 1; i <= noOfDisks; i++)
298  {
299  ss1.str(std::string()); ss1 << "charge_" + clustersrc_.label() + "_Disk_m" << i;
300  ss2.str(std::string()); ss2 << "Charge (off track, diskm" << i << ")";
301  meClChargeNotOnTrack_diskms.push_back(iBooker.book1D(ss1.str(),ss2.str(),500,0.,500.));
302  meClChargeNotOnTrack_diskms.at(i-1)->setAxisTitle("Charge size (in ke)",1);
303  }
304 
305  //size
306  //on track
307  iBooker.setCurrentFolder(topFolderName_+"/Clusters/OnTrack");
308  meClSizeOnTrack_all = iBooker.book1D("size_" + clustersrc_.label(),"Size (on track)",100,0.,100.);
309  meClSizeOnTrack_all->setAxisTitle("Cluster size (in pixels)",1);
310  meClSizeOnTrack_bpix = iBooker.book1D("size_" + clustersrc_.label() + "_Barrel","Size (on track, barrel)",100,0.,100.);
311  meClSizeOnTrack_bpix->setAxisTitle("Cluster size (in pixels)",1);
312  meClSizeOnTrack_fpix = iBooker.book1D("size_" + clustersrc_.label() + "_Endcap","Size (on track, endcap)",100,0.,100.);
313  meClSizeOnTrack_fpix->setAxisTitle("Cluster size (in pixels)",1);
314  for (int i = 1; i <= noOfLayers; i++)
315  {
316  ss1.str(std::string()); ss1 << "size_" + clustersrc_.label() + "_Layer_" << i;
317  ss2.str(std::string()); ss2 << "Size (on track, layer" << i << ")";
318  meClSizeOnTrack_layers.push_back(iBooker.book1D(ss1.str(),ss2.str(),500,0.,500.));
319  meClSizeOnTrack_layers.at(i-1)->setAxisTitle("Cluster size (in pixels)",1);
320  }
321  for (int i = 1; i <= noOfDisks; i++)
322  {
323  ss1.str(std::string()); ss1 << "size_" + clustersrc_.label() + "_Disk_p" << i;
324  ss2.str(std::string()); ss2 << "Size (on track, diskp" << i << ")";
325  meClSizeOnTrack_diskps.push_back(iBooker.book1D(ss1.str(),ss2.str(),500,0.,500.));
326  meClSizeOnTrack_diskps.at(i-1)->setAxisTitle("Cluster size (in pixels)",1);
327  }
328  for (int i = 1; i <= noOfDisks; i++)
329  {
330  ss1.str(std::string()); ss1 << "size_" + clustersrc_.label() + "_Disk_m1" << i;
331  ss2.str(std::string()); ss2 << "Size (on track, diskm" << i << ")";
332  meClSizeOnTrack_diskms.push_back(iBooker.book1D(ss1.str(),ss2.str(),500,0.,500.));
333  meClSizeOnTrack_diskms.at(i-1)->setAxisTitle("Cluster size (in pixels)",1);
334  }
335  meClSizeXOnTrack_all = iBooker.book1D("sizeX_" + clustersrc_.label(),"SizeX (on track)",100,0.,100.);
336  meClSizeXOnTrack_all->setAxisTitle("Cluster sizeX (in pixels)",1);
337  meClSizeXOnTrack_bpix = iBooker.book1D("sizeX_" + clustersrc_.label() + "_Barrel","SizeX (on track, barrel)",100,0.,100.);
338  meClSizeXOnTrack_bpix->setAxisTitle("Cluster sizeX (in pixels)",1);
339  meClSizeXOnTrack_fpix = iBooker.book1D("sizeX_" + clustersrc_.label() + "_Endcap","SizeX (on track, endcap)",100,0.,100.);
340  meClSizeXOnTrack_fpix->setAxisTitle("Cluster sizeX (in pixels)",1);
341  for (int i = 1; i <= noOfLayers; i++)
342  {
343  ss1.str(std::string()); ss1 << "sizeX_" + clustersrc_.label() + "_Layer_" << i;
344  ss2.str(std::string()); ss2 << "SizeX (on track, layer" << i << ")";
345  meClSizeXOnTrack_layers.push_back(iBooker.book1D(ss1.str(),ss2.str(),500,0.,500.));
346  meClSizeXOnTrack_layers.at(i-1)->setAxisTitle("Cluster size (in pixels)",1);
347  }
348  for (int i = 1; i <= noOfDisks; i++)
349  {
350  ss1.str(std::string()); ss1 << "sizeX_" + clustersrc_.label() + "_Disk_p" << i;
351  ss2.str(std::string()); ss2 << "SizeX (on track, diskp" << i << ")";
352  meClSizeXOnTrack_diskps.push_back(iBooker.book1D(ss1.str(),ss2.str(),500,0.,500.));
353  meClSizeXOnTrack_diskps.at(i-1)->setAxisTitle("Cluster size (in pixels)",1);
354  }
355  for (int i = 1; i <= noOfDisks; i++)
356  {
357  ss1.str(std::string()); ss1 << "sizeX_" + clustersrc_.label() + "_Disk_m" << i;
358  ss2.str(std::string()); ss2 << "SizeX (on track, diskm" << i << ")";
359  meClSizeXOnTrack_diskms.push_back(iBooker.book1D(ss1.str(),ss2.str(),500,0.,500.));
360  meClSizeXOnTrack_diskms.at(i-1)->setAxisTitle("Cluster size (in pixels)",1);
361  }
362  meClSizeYOnTrack_all = iBooker.book1D("sizeY_" + clustersrc_.label(),"SizeY (on track)",100,0.,100.);
363  meClSizeYOnTrack_all->setAxisTitle("Cluster sizeY (in pixels)",1);
364  meClSizeYOnTrack_bpix = iBooker.book1D("sizeY_" + clustersrc_.label() + "_Barrel","SizeY (on track, barrel)",100,0.,100.);
365  meClSizeYOnTrack_bpix->setAxisTitle("Cluster sizeY (in pixels)",1);
366  meClSizeYOnTrack_fpix = iBooker.book1D("sizeY_" + clustersrc_.label() + "_Endcap","SizeY (on track, endcap)",100,0.,100.);
367  meClSizeYOnTrack_fpix->setAxisTitle("Cluster sizeY (in pixels)",1);
368  for (int i = 1; i <= noOfLayers; i++)
369  {
370  ss1.str(std::string()); ss1 << "sizeY_" + clustersrc_.label() + "_Layer_" << i;
371  ss2.str(std::string()); ss2 << "SizeY (on track, layer" << i << ")";
372  meClSizeYOnTrack_layers.push_back(iBooker.book1D(ss1.str(),ss2.str(),500,0.,500.));
373  meClSizeYOnTrack_layers.at(i-1)->setAxisTitle("Cluster size (in pixels)",1);
374  }
375  for (int i = 1; i <= noOfDisks; i++)
376  {
377  ss1.str(std::string()); ss1 << "sizeY_" + clustersrc_.label() + "_Disk_p" << i;
378  ss2.str(std::string()); ss2 << "SizeY (on track, diskp" << i << ")";
379  meClSizeYOnTrack_diskps.push_back(iBooker.book1D(ss1.str(),ss2.str(),500,0.,500.));
380  meClSizeYOnTrack_diskps.at(i-1)->setAxisTitle("Cluster size (in pixels)",1);
381  }
382  for (int i = 1; i <= noOfDisks; i++)
383  {
384  ss1.str(std::string()); ss1 << "sizeY_" + clustersrc_.label() + "_Disk_m" << i;
385  ss2.str(std::string()); ss2 << "SizeY (on track, diskm" << i << ")";
386  meClSizeYOnTrack_diskms.push_back(iBooker.book1D(ss1.str(),ss2.str(),500,0.,500.));
387  meClSizeYOnTrack_diskms.at(i-1)->setAxisTitle("Cluster size (in pixels)",1);
388  }
389  //off track
390  iBooker.setCurrentFolder(topFolderName_+"/Clusters/OffTrack");
391  meClSizeNotOnTrack_all = iBooker.book1D("size_" + clustersrc_.label(),"Size (off track)",100,0.,100.);
392  meClSizeNotOnTrack_all->setAxisTitle("Cluster size (in pixels)",1);
393  meClSizeNotOnTrack_bpix = iBooker.book1D("size_" + clustersrc_.label() + "_Barrel","Size (off track, barrel)",100,0.,100.);
394  meClSizeNotOnTrack_bpix->setAxisTitle("Cluster size (in pixels)",1);
395  meClSizeNotOnTrack_fpix = iBooker.book1D("size_" + clustersrc_.label() + "_Endcap","Size (off track, endcap)",100,0.,100.);
396  meClSizeNotOnTrack_fpix->setAxisTitle("Cluster size (in pixels)",1);
397  for (int i = 1; i <= noOfLayers; i++)
398  {
399  ss1.str(std::string()); ss1 << "size_" + clustersrc_.label() + "_Layer_" << i;
400  ss2.str(std::string()); ss2 << "Size (off track, layer" << i << ")";
401  meClSizeNotOnTrack_layers.push_back(iBooker.book1D(ss1.str(),ss2.str(),500,0.,500.));
402  meClSizeNotOnTrack_layers.at(i-1)->setAxisTitle("Cluster size (in pixels)",1);
403  }
404  for (int i = 1; i <= noOfDisks; i++)
405  {
406  ss1.str(std::string()); ss1 << "size_" + clustersrc_.label() + "_Disk_p" << i;
407  ss2.str(std::string()); ss2 << "Size (off track, diskp" << i << ")";
408  meClSizeNotOnTrack_diskps.push_back(iBooker.book1D(ss1.str(),ss2.str(),500,0.,500.));
409  meClSizeNotOnTrack_diskps.at(i-1)->setAxisTitle("Cluster size (in pixels)",1);
410  }
411  for (int i = 1; i <= noOfDisks; i++)
412  {
413  ss1.str(std::string()); ss1 << "size_" + clustersrc_.label() + "_Disk_m" << i;
414  ss2.str(std::string()); ss2 << "Size (off track, diskm" << i << ")";
415  meClSizeNotOnTrack_diskms.push_back(iBooker.book1D(ss1.str(),ss2.str(),500,0.,500.));
416  meClSizeNotOnTrack_diskms.at(i-1)->setAxisTitle("Cluster size (in pixels)",1);
417  }
418  meClSizeXNotOnTrack_all = iBooker.book1D("sizeX_" + clustersrc_.label(),"SizeX (off track)",100,0.,100.);
419  meClSizeXNotOnTrack_all->setAxisTitle("Cluster sizeX (in pixels)",1);
420  meClSizeXNotOnTrack_bpix = iBooker.book1D("sizeX_" + clustersrc_.label() + "_Barrel","SizeX (off track, barrel)",100,0.,100.);
421  meClSizeXNotOnTrack_bpix->setAxisTitle("Cluster sizeX (in pixels)",1);
422  meClSizeXNotOnTrack_fpix = iBooker.book1D("sizeX_" + clustersrc_.label() + "_Endcap","SizeX (off track, endcap)",100,0.,100.);
423  meClSizeXNotOnTrack_fpix->setAxisTitle("Cluster sizeX (in pixels)",1);
424  for (int i = 1; i <= noOfLayers; i++)
425  {
426  ss1.str(std::string()); ss1 << "sizeX_" + clustersrc_.label() + "_Layer_" << i;
427  ss2.str(std::string()); ss2 << "SizeX (off track, layer" << i << ")";
428  meClSizeXNotOnTrack_layers.push_back(iBooker.book1D(ss1.str(),ss2.str(),500,0.,500.));
429  meClSizeXNotOnTrack_layers.at(i-1)->setAxisTitle("Cluster size (in pixels)",1);
430  }
431  for (int i = 1; i <= noOfDisks; i++)
432  {
433  ss1.str(std::string()); ss1 << "sizeX_" + clustersrc_.label() + "_Disk_p" << i;
434  ss2.str(std::string()); ss2 << "SizeX (off track, diskp" << i << ")";
435  meClSizeXNotOnTrack_diskps.push_back(iBooker.book1D(ss1.str(),ss2.str(),500,0.,500.));
436  meClSizeXNotOnTrack_diskps.at(i-1)->setAxisTitle("Cluster size (in pixels)",1);
437  }
438  for (int i = 1; i <= noOfDisks; i++)
439  {
440  ss1.str(std::string()); ss1 << "sizeX_" + clustersrc_.label() + "_Disk_m" << i;
441  ss2.str(std::string()); ss2 << "SizeX (off track, diskm" << i << ")";
442  meClSizeXNotOnTrack_diskms.push_back(iBooker.book1D(ss1.str(),ss2.str(),500,0.,500.));
443  meClSizeXNotOnTrack_diskms.at(i-1)->setAxisTitle("Cluster size (in pixels)",1);
444  }
445  meClSizeYNotOnTrack_all = iBooker.book1D("sizeY_" + clustersrc_.label(),"SizeY (off track)",100,0.,100.);
446  meClSizeYNotOnTrack_all->setAxisTitle("Cluster sizeY (in pixels)",1);
447  meClSizeYNotOnTrack_bpix = iBooker.book1D("sizeY_" + clustersrc_.label() + "_Barrel","SizeY (off track, barrel)",100,0.,100.);
448  meClSizeYNotOnTrack_bpix->setAxisTitle("Cluster sizeY (in pixels)",1);
449  meClSizeYNotOnTrack_fpix = iBooker.book1D("sizeY_" + clustersrc_.label() + "_Endcap","SizeY (off track, endcap)",100,0.,100.);
450  meClSizeYNotOnTrack_fpix->setAxisTitle("Cluster sizeY (in pixels)",1);
451  for (int i = 1; i <= noOfLayers; i++)
452  {
453  ss1.str(std::string()); ss1 << "sizeY_" + clustersrc_.label() + "_Layer_" << i;
454  ss2.str(std::string()); ss2 << "SizeY (off track, layer" << i << ")";
455  meClSizeYNotOnTrack_layers.push_back(iBooker.book1D(ss1.str(),ss2.str(),500,0.,500.));
456  meClSizeYNotOnTrack_layers.at(i-1)->setAxisTitle("Cluster size (in pixels)",1);
457  }
458  for (int i = 1; i <= noOfDisks; i++)
459  {
460  ss1.str(std::string()); ss1 << "sizeY_" + clustersrc_.label() + "_Disk_p" << i;
461  ss2.str(std::string()); ss2 << "SizeY (off track, diskp" << i << ")";
462  meClSizeYNotOnTrack_diskps.push_back(iBooker.book1D(ss1.str(),ss2.str(),500,0.,500.));
463  meClSizeYNotOnTrack_diskps.at(i-1)->setAxisTitle("Cluster size (in pixels)",1);
464  }
465  for (int i = 1; i <= noOfDisks; i++)
466  {
467  ss1.str(std::string()); ss1 << "sizeY_" + clustersrc_.label() + "_Disk_m" << i;
468  ss2.str(std::string()); ss2 << "SizeY (off track, diskm" << i << ")";
469  meClSizeYNotOnTrack_diskms.push_back(iBooker.book1D(ss1.str(),ss2.str(),500,0.,500.));
470  meClSizeYNotOnTrack_diskms.at(i-1)->setAxisTitle("Cluster size (in pixels)",1);
471  }
472 
473  //cluster global position
474  //on track
475  iBooker.setCurrentFolder(topFolderName_+"/Clusters/OnTrack");
476  //bpix
477  for (int i = 1; i <= noOfLayers; i++)
478  {
479  ss1.str(std::string()); ss1 << "position_" + clustersrc_.label() + "_Layer_" << i;
480  ss2.str(std::string()); ss2 << "Clusters Layer" << i << " (on track)";
481  meClPosLayersOnTrack.push_back(iBooker.book2D(ss1.str(),ss2.str(),200,-30.,30.,128,-3.2,3.2));
482  meClPosLayersOnTrack.at(i-1)->setAxisTitle("Global Z (cm)",1);
483  meClPosLayersOnTrack.at(i-1)->setAxisTitle("Global #phi",2);
484 
485  int ybins = -1; float ymin = 0.; float ymax = 0.;
486  if (i==1) { ybins = 23; ymin = -11.5; ymax = 11.5; }
487  if (i==2) { ybins = 33; ymin = -17.5; ymax = 17.5; }
488  if (i==3) { ybins = 45; ymin = -24.5; ymax = 24.5; }
489  ss1.str(std::string()); ss1 << "position_" + clustersrc_.label() + "_LadvsMod_Layer_" << i;
490  ss2.str(std::string()); ss2 << "Clusters Layer" << i << "_LadvsMod (on track)";
491  meClPosLayersLadVsModOnTrack.push_back(iBooker.book2D(ss1.str(),ss2.str(),11,-5.5,5.5,ybins,ymin,ymax));
492  meClPosLayersLadVsModOnTrack.at(i-1)->setAxisTitle("z-module",1);
493  meClPosLayersLadVsModOnTrack.at(i-1)->setAxisTitle("Ladder",2);
494 
495  }
496  //fpix
497  for (int i = 1; i <= noOfDisks; i++)
498  {
499  ss1.str(std::string()); ss1 << "position_" + clustersrc_.label() + "_pz_Disk_" << i;
500  ss2.str(std::string()); ss2 << "Clusters +Z Disk" << i << " (on track)";
501  meClPosDiskspzOnTrack.push_back(iBooker.book2D(ss1.str(),ss2.str(),80,-20.,20.,80,-20.,20.));
502  meClPosDiskspzOnTrack.at(i-1)->setAxisTitle("Global X (cm)",1);
503  meClPosDiskspzOnTrack.at(i-1)->setAxisTitle("Global Y (cm)",2);
504  }
505  for (int i = 1; i <= noOfDisks; i++)
506  {
507  ss1.str(std::string()); ss1 << "position_" + clustersrc_.label() + "_mz_Disk_" << i;
508  ss2.str(std::string()); ss2 << "Clusters -Z Disk" << i << " (on track)";
509  meClPosDisksmzOnTrack.push_back(iBooker.book2D(ss1.str(),ss2.str(),80,-20.,20.,80,-20.,20.));
510  meClPosDisksmzOnTrack.at(i-1)->setAxisTitle("Global X (cm)",1);
511  meClPosDisksmzOnTrack.at(i-1)->setAxisTitle("Global Y (cm)",2);
512  }
513  meNClustersOnTrack_all = iBooker.book1D("nclusters_" + clustersrc_.label(),"Number of Clusters (on Track)",50,0.,50.);
514  meNClustersOnTrack_all->setAxisTitle("Number of Clusters",1);
515  meNClustersOnTrack_bpix = iBooker.book1D("nclusters_" + clustersrc_.label() + "_Barrel","Number of Clusters (on track, barrel)",50,0.,50.);
516  meNClustersOnTrack_bpix->setAxisTitle("Number of Clusters",1);
517  meNClustersOnTrack_fpix = iBooker.book1D("nclusters_" + clustersrc_.label() + "_Endcap","Number of Clusters (on track, endcap)",50,0.,50.);
518  meNClustersOnTrack_fpix->setAxisTitle("Number of Clusters",1);
519  for (int i = 1; i <= noOfLayers; i++)
520  {
521  ss1.str(std::string()); ss1 << "nclusters_" + clustersrc_.label() + "_Layer_" << i;
522  ss2.str(std::string()); ss2 << "Number of Clusters (on track, layer" << i << ")";
523  meNClustersOnTrack_layers.push_back(iBooker.book1D(ss1.str(),ss2.str(),500,0.,500.));
524  meNClustersOnTrack_layers.at(i-1)->setAxisTitle("Number of Clusters",1);
525  }
526  for (int i = 1; i <= noOfDisks; i++)
527  {
528  ss1.str(std::string()); ss1 << "nclusters_" + clustersrc_.label() + "_Disk_p" << i;
529  ss2.str(std::string()); ss2 << "Number of Clusters (on track, diskp" << i << ")";
530  meNClustersOnTrack_diskps.push_back(iBooker.book1D(ss1.str(),ss2.str(),50,0.,50.));
531  meNClustersOnTrack_diskps.at(i-1)->setAxisTitle("Number of Clusters",1);
532  }
533  for (int i = 1; i <= noOfDisks; i++)
534  {
535  ss1.str(std::string()); ss1 << "nclusters_" + clustersrc_.label() + "_Disk_m" << i;
536  ss2.str(std::string()); ss2 << "Number of Clusters (on track, diskm" << i << ")";
537  meNClustersOnTrack_diskms.push_back(iBooker.book1D(ss1.str(),ss2.str(),500,0.,500.));
538  meNClustersOnTrack_diskms.at(i-1)->setAxisTitle("Number of Clusters",1);
539  }
540 
541  //not on track
542  iBooker.setCurrentFolder(topFolderName_+"/Clusters/OffTrack");
543  //bpix
544  for (int i = 1; i <= noOfLayers; i++)
545  {
546  ss1.str(std::string()); ss1 << "position_" + clustersrc_.label() + "_Layer_" << i;
547  ss2.str(std::string()); ss2 << "Clusters Layer" << i << " (off track)";
548  meClPosLayersNotOnTrack.push_back(iBooker.book2D(ss1.str(),ss2.str(),200,-30.,30.,128,-3.2,3.2));
549  meClPosLayersNotOnTrack.at(i-1)->setAxisTitle("Global Z (cm)",1);
550  meClPosLayersNotOnTrack.at(i-1)->setAxisTitle("Global #phi",2);
551  }
552  //fpix
553  for (int i = 1; i <= noOfDisks; i++)
554  {
555  ss1.str(std::string()); ss1 << "position_" + clustersrc_.label() + "_pz_Disk_" << i;
556  ss2.str(std::string()); ss2 << "Clusters +Z Disk" << i << " (off track)";
557  meClPosDiskspzNotOnTrack.push_back(iBooker.book2D(ss1.str(),ss2.str(),80,-20.,20.,80,-20.,20.));
558  meClPosDiskspzNotOnTrack.at(i-1)->setAxisTitle("Global X (cm)",1);
559  meClPosDiskspzNotOnTrack.at(i-1)->setAxisTitle("Global Y (cm)",2);
560  }
561  for (int i = 1; i <= noOfDisks; i++)
562  {
563  ss1.str(std::string()); ss1 << "position_" + clustersrc_.label() + "_mz_Disk_" << i;
564  ss2.str(std::string()); ss2 << "Clusters -Z Disk" << i << " (off track)";
565  meClPosDisksmzNotOnTrack.push_back(iBooker.book2D(ss1.str(),ss2.str(),80,-20.,20.,80,-20.,20.));
566  meClPosDisksmzNotOnTrack.at(i-1)->setAxisTitle("Global X (cm)",1);
567  meClPosDisksmzNotOnTrack.at(i-1)->setAxisTitle("Global Y (cm)",2);
568  }
569  meNClustersNotOnTrack_all = iBooker.book1D("nclusters_" + clustersrc_.label(),"Number of Clusters (off Track)",50,0.,50.);
570  meNClustersNotOnTrack_all->setAxisTitle("Number of Clusters",1);
571  meNClustersNotOnTrack_bpix = iBooker.book1D("nclusters_" + clustersrc_.label() + "_Barrel","Number of Clusters (off track, barrel)",50,0.,50.);
572  meNClustersNotOnTrack_bpix->setAxisTitle("Number of Clusters",1);
573  meNClustersNotOnTrack_fpix = iBooker.book1D("nclusters_" + clustersrc_.label() + "_Endcap","Number of Clusters (off track, endcap)",50,0.,50.);
574  meNClustersNotOnTrack_fpix->setAxisTitle("Number of Clusters",1);
575  for (int i = 1; i <= noOfLayers; i++)
576  {
577  ss1.str(std::string()); ss1 << "nclusters_" + clustersrc_.label() + "_Layer_" << i;
578  ss2.str(std::string()); ss2 << "Number of Clusters (off track, layer" << i << ")";
579  meNClustersNotOnTrack_layers.push_back(iBooker.book1D(ss1.str(),ss2.str(),500,0.,500.));
580  meNClustersNotOnTrack_layers.at(i-1)->setAxisTitle("Number of Clusters",1);
581  }
582  for (int i = 1; i <= noOfDisks; i++)
583  {
584  ss1.str(std::string()); ss1 << "nclusters_" + clustersrc_.label() + "_Disk_p" << i;
585  ss2.str(std::string()); ss2 << "Number of Clusters (off track, diskp" << i << ")";
586  meNClustersNotOnTrack_diskps.push_back(iBooker.book1D(ss1.str(),ss2.str(),500,0.,500.));
587  meNClustersNotOnTrack_diskps.at(i-1)->setAxisTitle("Number of Clusters",1);
588  }
589  for (int i = 1; i <= noOfDisks; i++)
590  {
591  ss1.str(std::string()); ss1 << "nclusters_" + clustersrc_.label() + "_Disk_m" << i;
592  ss2.str(std::string()); ss2 << "Number of Clusters (off track, diskm" << i << ")";
593  meNClustersNotOnTrack_diskms.push_back(iBooker.book1D(ss1.str(),ss2.str(),500,0.,500.));
594  meNClustersNotOnTrack_diskms.at(i-1)->setAxisTitle("Number of Clusters",1);
595 
596  }
597  //HitProbability
598  //on track
599  iBooker.setCurrentFolder(topFolderName_+"/Clusters/OnTrack");
600  meHitProbability = iBooker.book1D("FractionLowProb","Fraction of hits with low probability;FractionLowProb;#HitsOnTrack",100,0.,1.);
601 
602  if (debug_) {
603  // book summary residual histograms in a debugging folder - one (x,y) pair of histograms per subdetector
604  iBooker.setCurrentFolder("debugging");
605  char hisID[80];
606  for (int s=0; s<3; s++) {
607  sprintf(hisID,"residual_x_subdet_%i",s);
608  meSubdetResidualX[s] = iBooker.book1D(hisID,"Pixel Hit-to-Track Residual in X",500,-5.,5.);
609 
610  sprintf(hisID,"residual_y_subdet_%i",s);
611  meSubdetResidualY[s] = iBooker.book1D(hisID,"Pixel Hit-to-Track Residual in Y",500,-5.,5.);
612  }
613  }
614 
615  firstRun = false;
616 }
617 
618 
620  //Retrieve tracker topology from geometry
621  edm::ESHandle<TrackerTopology> tTopoHandle;
622  iSetup.get<TrackerTopologyRcd>().get(tTopoHandle);
623  const TrackerTopology* const tTopo = tTopoHandle.product();
624 
625  // retrieve TrackerGeometry again and MagneticField for use in transforming
626  // a TrackCandidate's P(ersistent)TrajectoryStateoOnDet (PTSoD) to a TrajectoryStateOnSurface (TSoS)
628  iSetup.get<TrackerDigiGeometryRecord>().get(TG);
629  const TrackerGeometry* theTrackerGeometry = TG.product();
630 
631  //analytic triplet method to calculate the track residuals in the pixe barrel detector
632 
633  //--------------------------------------------------------------------
634  // beam spot:
635  //
637  //iEvent.getByLabel( "offlineBeamSpot", rbs );
638  iEvent.getByToken( beamSpotToken_, rbs );
639  math::XYZPoint bsP = math::XYZPoint(0,0,0);
640  if( !rbs.failedToGet() && rbs.isValid() )
641  {
642  bsP = math::XYZPoint( rbs->x0(), rbs->y0(), rbs->z0() );
643  }
644 
645  //--------------------------------------------------------------------
646  // primary vertices:
647  //
649  //iEvent.getByLabel("offlinePrimaryVertices", vertices );
650  iEvent.getByToken( offlinePrimaryVerticesToken_, vertices );
651 
652  if( vertices.failedToGet() ) return;
653  if( !vertices.isValid() ) return;
654 
655  math::XYZPoint vtxN = math::XYZPoint(0,0,0);
656  math::XYZPoint vtxP = math::XYZPoint(0,0,0);
657 
658  double bestNdof = 0.0;
659  double maxSumPt = 0.0;
660  reco::Vertex bestPvx;
661  for(reco::VertexCollection::const_iterator iVertex = vertices->begin();
662  iVertex != vertices->end(); ++iVertex ) {
663  if( iVertex->ndof() > bestNdof ) {
664  bestNdof = iVertex->ndof();
665  vtxN = math::XYZPoint( iVertex->x(), iVertex->y(), iVertex->z() );
666  }//ndof
667  if( iVertex->p4().pt() > maxSumPt ) {
668  maxSumPt = iVertex->p4().pt();
669  vtxP = math::XYZPoint( iVertex->x(), iVertex->y(), iVertex->z() );
670  bestPvx = *iVertex;
671  }//sumpt
672 
673  }//vertex
674 
675  if( maxSumPt < 1.0 ) vtxP = vtxN;
676 
677  //---------------------------------------------
678  //get Tracks
679  //
681  //iEvent.getByLabel( "generalTracks", TracksForRes );
682  iEvent.getByToken( generalTracksToken_, TracksForRes );
683 
684  //
685  // transient track builder, needs B-field from data base (global tag in .py)
686  //
688  iSetup.get<TransientTrackRecord>().get( "TransientTrackBuilder", theB );
689 
690  //get the TransienTrackingRecHitBuilder needed for extracting the global position of the hits in the pixel
691  edm::ESHandle<TransientTrackingRecHitBuilder> theTrackerRecHitBuilder;
692  iSetup.get<TransientRecHitRecord>().get(ttrhbuilder_,theTrackerRecHitBuilder);
693 
694  //check that tracks are valid
695  if( TracksForRes.failedToGet() ) return;
696  if( !TracksForRes.isValid() ) return;
697 
698  //get tracker geometry
700  iSetup.get<TrackerDigiGeometryRecord>().get(pDD);
701 
702  if( !pDD.isValid() ) {
703  cout << "Unable to find TrackerDigiGeometry. Return\n";
704  return;
705  }
706 
707  int kk = -1;
708  //----------------------------------------------------------------------------
709  // Residuals:
710  //
711  for( reco::TrackCollection::const_iterator iTrack = TracksForRes->begin();
712  iTrack != TracksForRes->end(); ++iTrack ) {
713  //count
714  kk++;
715  //Calculate minimal track pt before curling
716  // cpt = cqRB = 0.3*R[m]*B[T] = 1.14*R[m] for B=3.8T
717  // D = 2R = 2*pt/1.14
718  // calo: D = 1.3 m => pt = 0.74 GeV/c
719  double pt = iTrack->pt();
720  if( pt < 0.75 ) continue;// curls up
721  if( abs( iTrack->dxy(vtxP) ) > 5*iTrack->dxyError() ) continue; // not prompt
722 
723  double charge = iTrack->charge();
724 
725  reco::TransientTrack tTrack = theB->build(*iTrack);
726  //get curvature of the track, needed for the residuals
727  double kap = tTrack.initialFreeState().transverseCurvature();
728  //needed for the TransienTrackingRecHitBuilder
730  if( iTrack->extra().isNonnull() &&iTrack->extra().isAvailable() ){
731 
732  double x1 = 0;
733  double y1 = 0;
734  double z1 = 0;
735  double x2 = 0;
736  double y2 = 0;
737  double z2 = 0;
738  double x3 = 0;
739  double y3 = 0;
740  double z3 = 0;
741  int n1 = 0;
742  int n2 = 0;
743  int n3 = 0;
744 
745  //for saving the pixel barrel hits
746  vector<TransientTrackingRecHit::RecHitPointer> GoodPixBarrelHits;
747  //looping through the RecHits of the track
748  for( trackingRecHit_iterator irecHit = iTrack->recHitsBegin();
749  irecHit != iTrack->recHitsEnd(); ++irecHit){
750 
751  if( (*irecHit)->isValid() ){
752  DetId detId = (*irecHit)->geographicalId();
753  // enum Detector { Tracker=1, Muon=2, Ecal=3, Hcal=4, Calo=5 };
754  if( detId.det() != 1 ){
755  if(debug_){
756  cout << "rec hit ID = " << detId.det() << " not in tracker!?!?\n";
757  }
758  continue;
759  }
760  uint32_t subDet = detId.subdetId();
761 
762  // enum SubDetector{ PixelBarrel=1, PixelEndcap=2 };
763  // enum SubDetector{ TIB=3, TID=4, TOB=5, TEC=6 };
764 
765  TransientTrackingRecHit::RecHitPointer trecHit = theTrackerRecHitBuilder->build( &*(*irecHit), initialTSOS);
766 
767 
768  double gX = trecHit->globalPosition().x();
769  double gY = trecHit->globalPosition().y();
770  double gZ = trecHit->globalPosition().z();
771 
772 
773  if( subDet == PixelSubdetector::PixelBarrel ) {
774 
775  int ilay = tTopo->pxbLayer(detId);
776 
777  if( ilay == 1 ){
778  n1++;
779  x1 = gX;
780  y1 = gY;
781  z1 = gZ;
782 
783  GoodPixBarrelHits.push_back((trecHit));
784  }//PXB1
785  if( ilay == 2 ){
786 
787  n2++;
788  x2 = gX;
789  y2 = gY;
790  z2 = gZ;
791 
792  GoodPixBarrelHits.push_back((trecHit));
793 
794  }//PXB2
795  if( ilay == 3 ){
796 
797  n3++;
798  x3 = gX;
799  y3 = gY;
800  z3 = gZ;
801  GoodPixBarrelHits.push_back((trecHit));
802  }
803  }//PXB
804 
805 
806  }//valid
807  }//loop rechits
808 
809  //CS extra plots
810 
811 
812  if( n1+n2+n3 == 3 && n1*n2*n3 > 0) {
813  for( unsigned int i = 0; i < GoodPixBarrelHits.size(); i++){
814 
815  if( GoodPixBarrelHits[i]->isValid() ){
816  DetId detId = GoodPixBarrelHits[i]->geographicalId().rawId();
817  int ilay = tTopo->pxbLayer(detId);
818  if(pt > ptminres_){
819 
820  double dca2 = 0.0, dz2=0.0;
821  double ptsig = pt;
822  if(charge<0.) ptsig = -pt;
823  //Filling the histograms in modules
824 
825  MeasurementPoint Test;
826  MeasurementPoint Test2;
827  Test=MeasurementPoint(0,0);
828  Test2=MeasurementPoint(0,0);
829  Measurement2DVector residual;
830 
831  if( ilay == 1 ){
832 
833  triplets(x2,y2,z2,x1,y1,z1,x3,y3,z3,ptsig,dca2,dz2, kap);
834 
835  Test=MeasurementPoint(dca2*1E4,dz2*1E4);
836  residual=Test-Test2;
837  }
838 
839  if( ilay == 2 ){
840 
841  triplets(x1,y1,z1,x2,y2,z2,x3,y3,z3,ptsig,dca2,dz2, kap);
842 
843  Test=MeasurementPoint(dca2*1E4,dz2*1E4);
844  residual=Test-Test2;
845 
846  }
847 
848  if( ilay == 3 ){
849 
850  triplets(x1,y1,z1,x3,y3,z3,x2,y2,z2,ptsig,dca2,dz2, kap);
851 
852  Test=MeasurementPoint(dca2*1E4,dz2*1E4);
853  residual=Test-Test2;
854  }
855  // fill the residual histograms
856 
857  std::map<uint32_t, SiPixelTrackResidualModule*>::iterator pxd = theSiPixelStructure.find(detId);
858  if (pxd!=theSiPixelStructure.end()) (*pxd).second->fill(residual, reducedSet, modOn, ladOn, layOn, phiOn, bladeOn, diskOn, ringOn);
859 
860  if(ladOn&&pxd!=theSiPixelStructure.end()){
861  meResidualXSummedLay.at(PixelBarrelNameUpgrade((*pxd).first).layerName()-1)->Fill(residual.x());
862  meResidualYSummedLay.at(PixelBarrelNameUpgrade((*pxd).first).layerName()-1)->Fill(residual.y());
863  }
864 
865  }//three hits
866  }//is valid
867  }//rechits loop
868  }//pt 4
869  }
870 
871 
872  }//-----Tracks
874  //get trajectories
875  edm::Handle<std::vector<Trajectory> > trajCollectionHandle;
876  //iEvent.getByLabel(tracksrc_,trajCollectionHandle);
877  iEvent.getByToken ( tracksrcToken_, trajCollectionHandle );
878  auto const & trajColl = *(trajCollectionHandle.product());
879 
880  //get tracks
881  edm::Handle<std::vector<reco::Track> > trackCollectionHandle;
882  //iEvent.getByLabel(tracksrc_,trackCollectionHandle);
883  iEvent.getByToken( trackToken_, trackCollectionHandle );
884  auto const & trackColl = *(trackCollectionHandle.product());
885 
886  //get the map
888  //iEvent.getByLabel(tracksrc_,match);
889  iEvent.getByToken( trackAssociationToken_, match);
890  auto const & ttac = *(match.product());
891 
892  // get clusters
894  //iEvent.getByLabel( clustersrc_, clusterColl );
895  iEvent.getByToken( clustersrcToken_, clusterColl );
896  auto const & clustColl = *(clusterColl.product());
897 
898 
899  if(debug_){
900  std::cout << "Trajectories\t : " << trajColl.size() << std::endl;
901  std::cout << "recoTracks \t : " << trackColl.size() << std::endl;
902  std::cout << "Map entries \t : " << ttac.size() << std::endl;
903  }
904 
905  std::set<SiPixelCluster> clusterSet;
906  TrajectoryStateCombiner tsoscomb;
907  int tracks=0, pixeltracks=0, bpixtracks=0, fpixtracks=0;
908  int trackclusters=0, barreltrackclusters=0, endcaptrackclusters=0;
909  int otherclusters=0, barrelotherclusters=0, endcapotherclusters=0;
910 
911  //Loop over map entries
912  for(TrajTrackAssociationCollection::const_iterator it = ttac.begin();it != ttac.end(); ++it){
913  const edm::Ref<std::vector<Trajectory> > traj_iterator = it->key;
914  // Trajectory Map, extract Trajectory for this track
915  reco::TrackRef trackref = it->val;
916  tracks++;
917 
918  bool isBpixtrack = false, isFpixtrack = false, crossesPixVol=false;
919 
920  //find out whether track crosses pixel fiducial volume (for cosmic tracks)
921 
922  double d0 = (*trackref).d0(), dz = (*trackref).dz();
923 
924  if(abs(d0)<15 && abs(dz)<50) crossesPixVol = true;
925 
926  const std::vector<TrajectoryMeasurement>& tmeasColl =traj_iterator->measurements();
927  std::vector<TrajectoryMeasurement>::const_iterator tmeasIt;
928  //loop on measurements to find out whether there are bpix and/or fpix hits
929  for(tmeasIt = tmeasColl.begin();tmeasIt!=tmeasColl.end();tmeasIt++){
930  if(! tmeasIt->updatedState().isValid()) continue;
931  TransientTrackingRecHit::ConstRecHitPointer testhit = tmeasIt->recHit();
932  if(! testhit->isValid() || testhit->geographicalId().det() != DetId::Tracker) continue;
933  uint testSubDetID = (testhit->geographicalId().subdetId());
934  if(testSubDetID==PixelSubdetector::PixelBarrel) isBpixtrack = true;
935  if(testSubDetID==PixelSubdetector::PixelEndcap) isFpixtrack = true;
936  }//end loop on measurements
937  if(isBpixtrack) {
938  bpixtracks++;
939  if(debug_) std::cout << "bpixtrack\n";
940  }
941  if(isFpixtrack) {
942  fpixtracks++;
943  if(debug_) std::cout << "fpixtrack\n";
944  }
945  if(isBpixtrack || isFpixtrack){
946  pixeltracks++;
947 
948  if(crossesPixVol) meNofTracksInPixVol_->Fill(0,1);
949 
950  const std::vector<TrajectoryMeasurement>& tmeasColl = traj_iterator->measurements();
951  for(std::vector<TrajectoryMeasurement>::const_iterator tmeasIt = tmeasColl.begin(); tmeasIt!=tmeasColl.end(); tmeasIt++){
952  if(! tmeasIt->updatedState().isValid()) continue;
953 
954  TrajectoryStateOnSurface tsos = tsoscomb( tmeasIt->forwardPredictedState(), tmeasIt->backwardPredictedState() );
955  if (!tsos.isValid()) continue; // Happens rarely, due to singular matrix or similar
956 
958  if(! hit->isValid() || hit->geographicalId().det() != DetId::Tracker ) {
959  continue;
960  } else {
961 
962 // //residual
963  const DetId & hit_detId = hit->geographicalId();
964  //uint IntRawDetID = (hit_detId.rawId());
965  uint IntSubDetID = (hit_detId.subdetId());
966 
967  if(IntSubDetID == 0 ) continue; // don't look at SiStrip hits!
968 
969  // get the enclosed persistent hit
970  const TrackingRecHit *persistentHit = hit->hit();
971  // check if it's not null, and if it's a valid pixel hit
972  if ((persistentHit != 0) && (typeid(*persistentHit) == typeid(SiPixelRecHit))) {
973  // tell the C++ compiler that the hit is a pixel hit
974  const SiPixelRecHit* pixhit = static_cast<const SiPixelRecHit*>( hit->hit() );
975  //Hit probability:
976  float hit_prob = -1.;
977  if(pixhit->hasFilledProb()){
978  hit_prob = pixhit->clusterProbability(0);
979  //std::cout<<"HITPROB= "<<hit_prob<<std::endl;
980  if(hit_prob<pow(10.,-15.)) NLowProb++;
981  NTotal++;
982  if(NTotal>0) meHitProbability->Fill(float(NLowProb/NTotal));
983  }
984 
985  // get the edm::Ref to the cluster
986  edm::Ref<edmNew::DetSetVector<SiPixelCluster>, SiPixelCluster> const& clust = (*pixhit).cluster();
987 
988  // check if the ref is not null
989  if (clust.isNonnull()) {
990 
991  //define tracker and pixel geometry and topology
992  const TrackerGeometry& theTracker(*theTrackerGeometry);
993  const PixelGeomDetUnit* theGeomDet = static_cast<const PixelGeomDetUnit*> (theTracker.idToDet(hit_detId) );
994 
995  //test if PixelGeomDetUnit exists
996  if(theGeomDet == 0) {
997  if(debug_) std::cout << "NO THEGEOMDET\n";
998  continue;
999  }
1000 
1001  const PixelTopology * topol = &(theGeomDet->specificTopology());
1002  //fill histograms for clusters on tracks
1003  //correct SiPixelTrackResidualModule
1004  std::map<uint32_t, SiPixelTrackResidualModule*>::iterator pxd = theSiPixelStructure.find((*hit).geographicalId().rawId());
1005 
1006  //CHARGE CORRECTION (for track impact angle)
1007  // calculate alpha and beta from cluster position
1009  LocalVector localDir = ltp.momentum()/ltp.momentum().mag();
1010 
1011  float clust_alpha = atan2(localDir.z(), localDir.x());
1012  float clust_beta = atan2(localDir.z(), localDir.y());
1013  double corrCharge = clust->charge() * sqrt( 1.0 / ( 1.0/pow( tan(clust_alpha), 2 ) +
1014  1.0/pow( tan(clust_beta ), 2 ) +
1015  1.0 )
1016  )/1000.;
1017 
1018  if (pxd!=theSiPixelStructure.end()) (*pxd).second->fill((*clust), true, corrCharge, reducedSet, modOn, ladOn, layOn, phiOn, bladeOn, diskOn, ringOn);
1019 
1020 
1021  trackclusters++;
1022  //CORR CHARGE
1023  meClChargeOnTrack_all->Fill(corrCharge);
1024  meClSizeOnTrack_all->Fill((*clust).size());
1025  meClSizeXOnTrack_all->Fill((*clust).sizeX());
1026  meClSizeYOnTrack_all->Fill((*clust).sizeY());
1027  clusterSet.insert(*clust);
1028 
1029  //find cluster global position (rphi, z)
1030  // get cluster center of gravity (of charge)
1031  float xcenter = clust->x();
1032  float ycenter = clust->y();
1033  // get the cluster position in local coordinates (cm)
1034  LocalPoint clustlp = topol->localPosition( MeasurementPoint(xcenter, ycenter) );
1035  // get the cluster position in global coordinates (cm)
1036  GlobalPoint clustgp = theGeomDet->surface().toGlobal( clustlp );
1037 
1038  //find location of hit (barrel or endcap, same for cluster)
1039  bool barrel = DetId((*hit).geographicalId()).subdetId() == static_cast<int>(PixelSubdetector::PixelBarrel);
1040  bool endcap = DetId((*hit).geographicalId()).subdetId() == static_cast<int>(PixelSubdetector::PixelEndcap);
1041  if(barrel) {
1042  barreltrackclusters++;
1043  //CORR CHARGE
1044  meClChargeOnTrack_bpix->Fill(corrCharge);
1045  meClSizeOnTrack_bpix->Fill((*clust).size());
1046  meClSizeXOnTrack_bpix->Fill((*clust).sizeX());
1047  meClSizeYOnTrack_bpix->Fill((*clust).sizeY());
1048  int DBlayer;
1049  DBlayer = PixelBarrelName(DetId((*hit).geographicalId()), tTopo, isUpgrade).layerName();
1050  float phi = clustgp.phi();
1051  float z = clustgp.z();
1052 
1053  PixelBarrelName pbn(DetId((*hit).geographicalId()), tTopo, isUpgrade);
1054  int ladder = pbn.ladderName();
1055  int module = pbn.moduleName();
1056 
1057  PixelBarrelName::Shell sh = pbn.shell(); //enum
1058  int ladderSigned=ladder;
1059  int moduleSigned=module;
1060  // Shell { mO = 1, mI = 2 , pO =3 , pI =4 };
1061  int shell = int(sh);
1062  // change the module sign for z<0
1063  if(shell==1 || shell==2) { moduleSigned = -module; }
1064  // change ladeer sign for Outer )x<0)
1065  if(shell==1 || shell==3) { ladderSigned = -ladder; }
1066 
1067  for (int i = 0; i < noOfLayers; i++)
1068  {
1069  if (DBlayer == i + 1) {
1070  meClPosLayersOnTrack.at(i)->Fill(z,phi);
1071  meClPosLayersLadVsModOnTrack.at(i)->Fill(moduleSigned,ladderSigned);
1072  meClChargeOnTrack_layers.at(i)->Fill(corrCharge);
1073  meClSizeOnTrack_layers.at(i)->Fill((*clust).size());
1074  meClSizeXOnTrack_layers.at(i)->Fill((*clust).sizeX());
1075  meClSizeYOnTrack_layers.at(i)->Fill((*clust).sizeY());
1076  }
1077  }
1078  }
1079  if(endcap) {
1080  endcaptrackclusters++;
1081  //CORR CHARGE
1082  meClChargeOnTrack_fpix->Fill(corrCharge);
1083  meClSizeOnTrack_fpix->Fill((*clust).size());
1084  meClSizeXOnTrack_fpix->Fill((*clust).sizeX());
1085  meClSizeYOnTrack_fpix->Fill((*clust).sizeY());
1086  int DBdisk = 0;
1087  DBdisk = PixelEndcapName(DetId((*hit).geographicalId()), tTopo, isUpgrade).diskName();
1088  float x = clustgp.x();
1089  float y = clustgp.y();
1090  float z = clustgp.z();
1091  if(z>0){
1092  for (int i = 0; i < noOfDisks; i++)
1093  {
1094  if (DBdisk == i + 1) {
1095  meClPosDiskspzOnTrack.at(i)->Fill(x,y);
1096  meClChargeOnTrack_diskps.at(i)->Fill(corrCharge);
1097  meClSizeOnTrack_diskps.at(i)->Fill((*clust).size());
1098  meClSizeXOnTrack_diskps.at(i)->Fill((*clust).sizeX());
1099  meClSizeYOnTrack_diskps.at(i)->Fill((*clust).sizeY());
1100  }
1101  }
1102  }
1103  else{
1104  for (int i = 0; i < noOfDisks; i++)
1105  {
1106  if (DBdisk == i + 1) {
1107  meClPosDisksmzOnTrack.at(i)->Fill(x,y);
1108  meClChargeOnTrack_diskms.at(i)->Fill(corrCharge);
1109  meClSizeOnTrack_diskms.at(i)->Fill((*clust).size());
1110  meClSizeXOnTrack_diskms.at(i)->Fill((*clust).sizeX());
1111  meClSizeYOnTrack_diskms.at(i)->Fill((*clust).sizeY());
1112  }
1113  }
1114  }
1115  }
1116 
1117  }//end if (cluster exists)
1118 
1119  }//end if (persistent hit exists and is pixel hit)
1120 
1121  }//end of else
1122 
1123 
1124  }//end for (all traj measurements of pixeltrack)
1125  }//end if (is pixeltrack)
1126  else {
1127  if(debug_) std::cout << "no pixeltrack:\n";
1128  if(crossesPixVol) meNofTracksInPixVol_->Fill(1,1);
1129  }
1130 
1131  }//end loop on map entries
1132 
1133  //find clusters that are NOT on track
1134  //edmNew::DetSet<SiPixelCluster>::const_iterator di;
1135  if(debug_) std::cout << "clusters not on track: (size " << clustColl.size() << ") ";
1136 
1137  for(TrackerGeometry::DetContainer::const_iterator it = TG->dets().begin(); it != TG->dets().end(); it++){
1138  //if(dynamic_cast<PixelGeomDetUnit const *>((*it))!=0){
1139  DetId detId = (*it)->geographicalId();
1140  if(detId>=302055684 && detId<=352477708){ // make sure it's a Pixel module WITHOUT using dynamic_cast!
1141  int nofclOnTrack = 0, nofclOffTrack=0;
1142  float z=0.;
1143  edmNew::DetSetVector<SiPixelCluster>::const_iterator isearch = clustColl.find(detId);
1144  if( isearch != clustColl.end() ) { // Not an empty iterator
1146  for(di=isearch->begin(); di!=isearch->end(); di++){
1147  unsigned int temp = clusterSet.size();
1148  clusterSet.insert(*di);
1149  //check if cluster is off track
1150  if(clusterSet.size()>temp) {
1151  otherclusters++;
1152  nofclOffTrack++;
1153  //fill histograms for clusters off tracks
1154  //correct SiPixelTrackResidualModule
1155  std::map<uint32_t, SiPixelTrackResidualModule*>::iterator pxd = theSiPixelStructure.find((*it)->geographicalId().rawId());
1156 
1157  if (pxd!=theSiPixelStructure.end()) (*pxd).second->fill((*di), false, -1., reducedSet, modOn, ladOn, layOn, phiOn, bladeOn, diskOn, ringOn);
1158 
1159 
1160 
1161  meClSizeNotOnTrack_all->Fill((*di).size());
1162  meClSizeXNotOnTrack_all->Fill((*di).sizeX());
1163  meClSizeYNotOnTrack_all->Fill((*di).sizeY());
1164  meClChargeNotOnTrack_all->Fill((*di).charge()/1000);
1165 
1167  //find cluster global position (rphi, z) get cluster
1168  //define tracker and pixel geometry and topology
1169  const TrackerGeometry& theTracker(*theTrackerGeometry);
1170  const PixelGeomDetUnit* theGeomDet = static_cast<const PixelGeomDetUnit*> (theTracker.idToDet(detId) );
1171  //test if PixelGeomDetUnit exists
1172  if(theGeomDet == 0) {
1173  if(debug_) std::cout << "NO THEGEOMDET\n";
1174  continue;
1175  }
1176  const PixelTopology * topol = &(theGeomDet->specificTopology());
1177 
1178  //center of gravity (of charge)
1179  float xcenter = di->x();
1180  float ycenter = di->y();
1181  // get the cluster position in local coordinates (cm)
1182  LocalPoint clustlp = topol->localPosition( MeasurementPoint(xcenter, ycenter) );
1183  // get the cluster position in global coordinates (cm)
1184  GlobalPoint clustgp = theGeomDet->surface().toGlobal( clustlp );
1185 
1187 
1188  //barrel
1189  if(DetId(detId).subdetId() == 1) {
1190  meClSizeNotOnTrack_bpix->Fill((*di).size());
1191  meClSizeXNotOnTrack_bpix->Fill((*di).sizeX());
1192  meClSizeYNotOnTrack_bpix->Fill((*di).sizeY());
1193  meClChargeNotOnTrack_bpix->Fill((*di).charge()/1000);
1194  barrelotherclusters++;
1195  int DBlayer = PixelBarrelName(DetId(detId), tTopo, isUpgrade).layerName();
1196  float phi = clustgp.phi();
1197  //float r = clustgp.perp();
1198  z = clustgp.z();
1199  for (int i = 0; i < noOfLayers; i++)
1200  {
1201  if (DBlayer == i + 1) {
1202  meClPosLayersNotOnTrack.at(i)->Fill(z,phi);
1203  meClSizeNotOnTrack_layers.at(i)->Fill((*di).size());
1204  meClSizeXNotOnTrack_layers.at(i)->Fill((*di).sizeX());
1205  meClSizeYNotOnTrack_layers.at(i)->Fill((*di).sizeY());
1206  meClChargeNotOnTrack_layers.at(i)->Fill((*di).charge()/1000);
1207  }
1208  }
1209  }
1210  //endcap
1211  if(DetId(detId).subdetId() == 2) {
1212  meClSizeNotOnTrack_fpix->Fill((*di).size());
1213  meClSizeXNotOnTrack_fpix->Fill((*di).sizeX());
1214  meClSizeYNotOnTrack_fpix->Fill((*di).sizeY());
1215  meClChargeNotOnTrack_fpix->Fill((*di).charge()/1000);
1216  endcapotherclusters++;
1217  int DBdisk = PixelEndcapName(DetId(detId), tTopo, isUpgrade).diskName();
1218  float x = clustgp.x();
1219  float y = clustgp.y();
1220  z = clustgp.z();
1221  if(z>0){
1222  for (int i = 0; i < noOfDisks; i++)
1223  {
1224  if (DBdisk == i + 1) {
1225  meClPosDiskspzNotOnTrack.at(i)->Fill(x,y);
1226  meClSizeNotOnTrack_diskps.at(i)->Fill((*di).size());
1227  meClSizeXNotOnTrack_diskps.at(i)->Fill((*di).sizeX());
1228  meClSizeYNotOnTrack_diskps.at(i)->Fill((*di).sizeY());
1229  meClChargeNotOnTrack_diskps.at(i)->Fill((*di).charge()/1000);
1230  }
1231  }
1232  }
1233  else{
1234  for (int i = 0; i < noOfDisks; i++)
1235  {
1236  if (DBdisk == i + 1) {
1237  meClPosDisksmzNotOnTrack.at(i)->Fill(x,y);
1238  meClSizeNotOnTrack_diskms.at(i)->Fill((*di).size());
1239  meClSizeXNotOnTrack_diskms.at(i)->Fill((*di).sizeX());
1240  meClSizeYNotOnTrack_diskms.at(i)->Fill((*di).sizeY());
1241  meClChargeNotOnTrack_diskms.at(i)->Fill((*di).charge()/1000);
1242  }
1243  }
1244  }
1245 
1246  }
1247  }// end "if cluster off track"
1248  else {
1249  nofclOnTrack++;
1250  if(z == 0){
1251  //find cluster global position (rphi, z) get cluster
1252  //define tracker and pixel geometry and topology
1253  const TrackerGeometry& theTracker(*theTrackerGeometry);
1254  const PixelGeomDetUnit* theGeomDet = static_cast<const PixelGeomDetUnit*> (theTracker.idToDet(detId) );
1255  //test if PixelGeomDetUnit exists
1256  if(theGeomDet == 0) {
1257  if(debug_) std::cout << "NO THEGEOMDET\n";
1258  continue;
1259  }
1260  const PixelTopology * topol = &(theGeomDet->specificTopology());
1261  //center of gravity (of charge)
1262  float xcenter = di->x();
1263  float ycenter = di->y();
1264  // get the cluster position in local coordinates (cm)
1265  LocalPoint clustlp = topol->localPosition( MeasurementPoint(xcenter, ycenter) );
1266  // get the cluster position in global coordinates (cm)
1267  GlobalPoint clustgp = theGeomDet->surface().toGlobal( clustlp );
1268  z = clustgp.z();
1269  }
1270  }
1271  }
1272  }
1273  //++ fill the number of clusters on a module
1274  std::map<uint32_t, SiPixelTrackResidualModule*>::iterator pxd = theSiPixelStructure.find((*it)->geographicalId().rawId());
1275  if (pxd!=theSiPixelStructure.end()) (*pxd).second->nfill(nofclOnTrack, nofclOffTrack, reducedSet, modOn, ladOn, layOn, phiOn, bladeOn, diskOn, ringOn);
1276  if(nofclOnTrack!=0) meNClustersOnTrack_all->Fill(nofclOnTrack);
1277  if(nofclOffTrack!=0) meNClustersNotOnTrack_all->Fill(nofclOffTrack);
1278  //barrel
1279  if(DetId(detId).subdetId() == 1){
1280  if(nofclOnTrack!=0) meNClustersOnTrack_bpix->Fill(nofclOnTrack);
1281  if(nofclOffTrack!=0) meNClustersNotOnTrack_bpix->Fill(nofclOffTrack);
1282  int DBlayer = PixelBarrelName(DetId(detId), tTopo, isUpgrade).layerName();
1283  for (int i = 0; i < noOfLayers; i++)
1284  {
1285  if (DBlayer == i + 1) {
1286  if(nofclOnTrack!=0) meNClustersOnTrack_layers.at(i)->Fill(nofclOnTrack);
1287  if(nofclOffTrack!=0) meNClustersNotOnTrack_layers.at(i)->Fill(nofclOffTrack);
1288  }
1289  }
1290  }//end barrel
1291  //endcap
1292  if(DetId(detId).subdetId() == 2) {
1293  int DBdisk = PixelEndcapName(DetId(detId )).diskName();
1294  //z = clustgp.z();
1295  if(nofclOnTrack!=0) meNClustersOnTrack_fpix->Fill(nofclOnTrack);
1296  if(nofclOffTrack!=0) meNClustersNotOnTrack_fpix->Fill(nofclOffTrack);
1297  if(z>0){
1298  for (int i = 0; i < noOfDisks; i++)
1299  {
1300  if (DBdisk == i + 1) {
1301  if(nofclOnTrack!=0) meNClustersOnTrack_diskps.at(i)->Fill(nofclOnTrack);
1302  if(nofclOffTrack!=0) meNClustersNotOnTrack_diskps.at(i)->Fill(nofclOffTrack);
1303  }
1304  }
1305  }
1306  if(z<0){
1307  for (int i = 0; i < noOfDisks; i++)
1308  {
1309  if (DBdisk == i + 1) {
1310  if(nofclOnTrack!=0) meNClustersOnTrack_diskms.at(i)->Fill(nofclOnTrack);
1311  if(nofclOffTrack!=0) meNClustersNotOnTrack_diskms.at(i)->Fill(nofclOffTrack);
1312  }
1313  }
1314  }
1315  }
1316 
1317  }//end if it's a Pixel module
1318  }//end for loop over tracker detector geometry modules
1319 
1320  if(trackclusters>0) (meNofClustersOnTrack_)->Fill(0,trackclusters);
1321  if(barreltrackclusters>0)(meNofClustersOnTrack_)->Fill(1,barreltrackclusters);
1322  if(endcaptrackclusters>0)(meNofClustersOnTrack_)->Fill(2,endcaptrackclusters);
1323  if(otherclusters>0)(meNofClustersNotOnTrack_)->Fill(0,otherclusters);
1324  if(barrelotherclusters>0)(meNofClustersNotOnTrack_)->Fill(1,barrelotherclusters);
1325  if(endcapotherclusters>0)(meNofClustersNotOnTrack_)->Fill(2,endcapotherclusters);
1326  if(tracks>0)(meNofTracks_)->Fill(0,tracks);
1327  if(pixeltracks>0)(meNofTracks_)->Fill(1,pixeltracks);
1328  if(bpixtracks>0)(meNofTracks_)->Fill(2,bpixtracks);
1329  if(fpixtracks>0)(meNofTracks_)->Fill(3,fpixtracks);
1330 }
1331 void SiPixelTrackResidualSource::triplets(double x1,double y1,double z1,double x2,double y2,double z2,double x3,double y3,double z3,
1332  double ptsig, double & dca2,double & dz2, double kap) {
1333 
1334  //Define some constants
1335  using namespace std;
1336 
1337  //Curvature kap from global Track
1338 
1339  //inverse of the curvature is the radius in the transverse plane
1340  double rho = 1/kap;
1341  //Check that the hits are in the correct layers
1342  double r1 = sqrt( x1*x1 + y1*y1 );
1343  double r3 = sqrt( x3*x3 + y3*y3 );
1344 
1345  if( r3-r1 < 2.0 ) cout << "warn r1 = " << r1 << ", r3 = " << r3 << endl;
1346 
1347  // Calculate the centre of the helix in xy-projection with radius rho from the track.
1348  //start with a line (sekante) connecting the two points (x1,y1) and (x3,y3) vec_L = vec_x3-vec_x1
1349  //with L being the length of that vector.
1350  double L=sqrt((x3-x1)*(x3-x1)+(y3-y1)*(y3-y1));
1351  //lam is the line from the middel point of vec_q towards the center of the circle X0,Y0
1352  // we already have kap and rho = 1/kap
1353  double lam = sqrt(rho*rho - L*L/4);
1354 
1355  // There are two solutions, the sign of kap gives the information
1356  // which of them is correct.
1357  //
1358  if( kap > 0 ) lam = -lam;
1359 
1360  //
1361  // ( X0, Y0 ) is the centre of the circle that describes the helix in xy-projection.
1362  //
1363  double x0 = 0.5*( x1 + x3 ) + lam/L * ( -y1 + y3 );
1364  double y0 = 0.5*( y1 + y3 ) + lam/L * ( x1 - x3 );
1365 
1366  // Calculate the dipangle in z direction (needed later for z residual) :
1367  //Starting from the heliz equation whihc has to hold for both points z1,z3
1368  double num = ( y3 - y0 ) * ( x1 - x0 ) - ( x3 - x0 ) * ( y1 - y0 );
1369  double den = ( x1 - x0 ) * ( x3 - x0 ) + ( y1 - y0 ) * ( y3 - y0 );
1370  double tandip = kap * ( z3 - z1 ) / atan( num / den );
1371 
1372 
1373  // angle from first hit to dca point:
1374  //
1375  double dphi = atan( ( ( x1 - x0 ) * y0 - ( y1 - y0 ) * x0 )
1376  / ( ( x1 - x0 ) * x0 + ( y1 - y0 ) * y0 ) );
1377  //z position of the track based on the middle of the circle
1378  //track equation for the z component
1379  double uz0 = z1 + tandip * dphi * rho;
1380 
1382  //RESIDUAL IN R-PHI
1384  //Calculate distance dca2 from point (x2,y2) to the circle which is given by
1385  //the distance of the point to the middlepoint dcM = sqrt((x0-x2)^2+(y0-y2)) and rho
1386  //dca = rho +- dcM
1387  if(kap>0) dca2=rho-sqrt((x0-x2)*(x0-x2)+(y0-y2)*(y0-y2));
1388  else dca2=rho+sqrt((-x0+x2)*(-x0+x2)+(-y0+y2)*(-y0+y2));
1389 
1391  //RESIDUAL IN Z
1393  double xx =0 ;
1394  double yy =0 ;
1395  //sign of kappa determines the calculation
1396  //xx and yy are the new coordinates starting from x2, y2 that are on the track itself
1397  //vec_X2+-dca2*vec(X0-X2)/|(X0-X2)|
1398  if(kap<0){
1399  xx = x2+(dca2*((x0-x2))/sqrt((x0-x2)*(x0-x2)+(y0-y2)*(y0-y2)));
1400  yy = y2+(dca2*((y0-y2))/sqrt((x0-x2)*(x0-x2)+(y0-y2)*(y0-y2)));
1401  }
1402  else if(kap>=0){
1403  xx = x2-(dca2*((x0-x2))/sqrt((x0-x2)*(x0-x2)+(y0-y2)*(y0-y2)));
1404  yy = y2-(dca2*((y0-y2))/sqrt((x0-x2)*(x0-x2)+(y0-y2)*(y0-y2)));
1405  }
1406 
1407  //to get residual in z start with calculating the new uz2 position if one has moved to xx, yy
1408  //on the track. First calculate the change in phi2 with respect to the center X0, Y0
1409  double dphi2 = atan( ( ( xx - x0 ) * y0 - ( yy - y0 ) * x0 )
1410  / ( ( xx - x0 ) * x0 + ( yy - y0 ) * y0 ) );
1411  //Solve track equation for this new z depending on the dip angle of the track (see above
1412  //calculated based on X1, X3 and X0, use uz0 as reference point again.
1413  double uz2= uz0 - dphi2*tandip*rho;
1414 
1415  //subtract new z position from the old one
1416  dz2=z2-uz2;
1417 
1418  //if we are interested in the arclength this is unsigned though
1419  // double cosphi2 = (x2*xx+y2*yy)/(sqrt(x2*x2+y2*y2)*sqrt(xx*xx+yy*yy));
1420  //double arcdca2=sqrt(x2*x2+y2*y2)*acos(cosphi2);
1421 
1422 }
1423 
1424 
1425 DEFINE_FWK_MODULE(SiPixelTrackResidualSource); // define this as a plug-in
GlobalPoint toGlobal(const Point2DBase< Scalar, LocalTag > lp) const
Definition: Surface.h:106
T getParameter(std::string const &) const
edm::EDGetTokenT< reco::TrackCollection > generalTracksToken_
T getUntrackedParameter(std::string const &, T const &) const
int i
Definition: DBlmapReader.cc:9
bool hasFilledProb() const
boost::transform_iterator< IterHelp, const_IdIter > const_iterator
const_iterator end(bool update=false) const
std::vector< MonitorElement * > meClChargeOnTrack_layers
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:252
float clusterProbability(unsigned int flags=0) const
T y() const
Definition: PV2DBase.h:46
std::vector< MonitorElement * > meClSizeYNotOnTrack_diskps
edm::EDGetTokenT< std::vector< Trajectory > > tracksrcToken_
virtual LocalPoint localPosition(const MeasurementPoint &) const =0
const LocalTrajectoryParameters & localParameters() const
int moduleName() const
module id (index in z)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:462
std::vector< MonitorElement * > meClSizeXNotOnTrack_diskms
virtual void dqmBeginRun(const edm::Run &r, edm::EventSetup const &iSetup)
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
virtual void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
std::vector< MonitorElement * > meClSizeXOnTrack_layers
std::vector< MonitorElement * > meNClustersNotOnTrack_diskps
Geom::Phi< T > phi() const
Definition: PV3DBase.h:69
T y() const
Definition: PV3DBase.h:63
std::vector< MonitorElement * > meClSizeYOnTrack_diskms
void setBinLabel(int bin, const std::string &label, int axis=1)
set bin label for x, y or z axis (axis=1, 2, 3 respectively)
std::vector< MonitorElement * > meClSizeNotOnTrack_diskms
std::vector< MonitorElement * > meClPosLayersOnTrack
std::vector< MonitorElement * > meClPosDisksmzOnTrack
std::vector< MonitorElement * > meClSizeOnTrack_layers
std::vector< MonitorElement * > meClPosDiskspzOnTrack
data_type const * const_iterator
Definition: DetSetNew.h:30
std::vector< MonitorElement * > meClChargeNotOnTrack_diskps
key_type key() const
Accessor for product key.
Definition: Ref.h:264
std::vector< MonitorElement * > meClChargeOnTrack_diskms
std::vector< MonitorElement * > meClPosDiskspzNotOnTrack
const Plane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:40
SiPixelTrackResidualSource(const edm::ParameterSet &)
std::vector< MonitorElement * > meClPosLayersNotOnTrack
std::vector< MonitorElement * > meResidualXSummedLay
tuple gX
Definition: corrVsCorr.py:109
TrajectoryStateOnSurface innermostMeasurementState() const
void Fill(long long x)
uint32_t rawId() const
get the raw id
Definition: DetId.h:43
std::vector< MonitorElement * > meClSizeXNotOnTrack_layers
Measurement2DPoint MeasurementPoint
Measurement points are two-dimensional by default.
int iEvent
Definition: GenABIO.cc:230
T mag() const
Definition: PV3DBase.h:67
std::vector< MonitorElement * > meClChargeOnTrack_diskps
std::shared_ptr< TrackingRecHit const > ConstRecHitPointer
void Fill(HcalDetId &id, double val, std::vector< TH2F > &depth)
T sqrt(T t)
Definition: SSEVec.h:18
std::vector< MonitorElement * > meResidualYSummedLay
bool setModuleFolder(const uint32_t &rawdetid=0, int type=0, bool isUpgrade=false)
Set folder name for a module or plaquette.
std::vector< MonitorElement * > meClSizeYNotOnTrack_diskms
T z() const
Definition: PV3DBase.h:64
void triplets(double x1, double y1, double z1, double x2, double y2, double z2, double x3, double y3, double z3, double ptsig, double &dc, double &dz, double kap)
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:115
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
LocalVector momentum() const
Momentum vector in the local frame.
std::vector< MonitorElement * > meClPosDisksmzNotOnTrack
edm::EDGetTokenT< std::vector< reco::Track > > trackToken_
std::vector< MonitorElement * > meClSizeYOnTrack_layers
std::vector< MonitorElement * > meNClustersOnTrack_diskps
bool isValid() const
Definition: HandleBase.h:75
int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:37
edm::EDGetTokenT< edmNew::DetSetVector< SiPixelCluster > > clustersrcToken_
std::shared_ptr< TrackingRecHit const > RecHitPointer
double ndof() const
Definition: Vertex.h:95
std::vector< MonitorElement * > meClSizeNotOnTrack_diskps
std::vector< MonitorElement * > meClSizeXOnTrack_diskps
unsigned int pxbLayer(const DetId &id) const
std::vector< MonitorElement * > meClSizeYNotOnTrack_layers
std::vector< MonitorElement * > meClSizeYOnTrack_diskps
FreeTrajectoryState initialFreeState() const
bool failedToGet() const
Definition: HandleBase.h:79
Definition: DetId.h:18
virtual TrackingRecHit const * hit() const
int ladderName() const
ladder id (index in phi)
std::vector< MonitorElement * > meClPosLayersLadVsModOnTrack
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:276
std::map< uint32_t, SiPixelTrackResidualModule * > theSiPixelStructure
T const * product() const
Definition: Handle.h:81
MonitorElement * book2D(Args &&...args)
Definition: DQMStore.h:133
tuple tracks
Definition: testEve_cfg.py:39
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
std::vector< MonitorElement * > meClChargeNotOnTrack_diskms
std::vector< MonitorElement * > meClSizeOnTrack_diskps
const T & get() const
Definition: EventSetup.h:56
virtual const PixelTopology & specificTopology() const
Returns a reference to the pixel proxy topology.
int layerName() const
layer id
T const * product() const
Definition: ESHandle.h:86
edm::EDGetTokenT< reco::VertexCollection > offlinePrimaryVerticesToken_
Shell shell() const
tuple gY
Definition: corrVsCorr.py:110
double transverseCurvature() const
std::string const & label() const
Definition: InputTag.h:36
std::vector< MonitorElement * > meNClustersOnTrack_diskms
Pixel cluster – collection of neighboring pixels above threshold.
std::vector< MonitorElement * > meNClustersNotOnTrack_diskms
tuple cout
Definition: gather_cfg.py:145
std::vector< MonitorElement * > meClSizeXNotOnTrack_diskps
int diskName() const
disk id
std::pair< typename Association::data_type::first_type, double > match(Reference key, Association association, bool bestMatchByMaxValue)
Generic matching function.
Definition: Utils.h:10
std::vector< MonitorElement * > meClChargeNotOnTrack_layers
virtual void analyze(const edm::Event &, const edm::EventSetup &)
std::vector< MonitorElement * > meNClustersOnTrack_layers
volatile std::atomic< bool > shutdown_flag false
std::vector< MonitorElement * > meClSizeNotOnTrack_layers
bool isValid() const
Definition: ESHandle.h:47
std::vector< MonitorElement * > meNClustersNotOnTrack_layers
Detector det() const
get the detector field from this detid
Definition: DetId.h:35
std::vector< MonitorElement * > meClSizeOnTrack_diskms
T x() const
Definition: PV2DBase.h:45
edm::EDGetTokenT< reco::BeamSpot > beamSpotToken_
T x() const
Definition: PV3DBase.h:62
void setAxisTitle(const std::string &title, int axis=1)
set x-, y- or z-axis title (axis=1, 2, 3 respectively)
Definition: vlib.h:208
edm::EDGetTokenT< TrajTrackAssociationCollection > trackAssociationToken_
int layerName() const
layer id
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
std::vector< MonitorElement * > meClSizeXOnTrack_diskms
const_iterator begin(bool update=false) const
TrackingRecHitCollection::base::const_iterator trackingRecHit_iterator
iterator over a vector of reference to TrackingRecHit in the same collection
Definition: Run.h:43
virtual const TrackerGeomDet * idToDet(DetId) const
Our base class.
Definition: SiPixelRecHit.h:23