CMS 3D CMS Logo

List of all members | Public Member Functions | Protected Member Functions | Private Member Functions | Private Attributes | Static Private Attributes
CTPPSPixelDQMSource Class Reference
Inheritance diagram for CTPPSPixelDQMSource:

Public Member Functions

 CTPPSPixelDQMSource (const edm::ParameterSet &ps)
 
 ~CTPPSPixelDQMSource () override
 

Protected Member Functions

void analyze (edm::Event const &e, edm::EventSetup const &eSetup) override
 
void bookHistograms (DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
 
void dqmBeginRun (edm::Run const &, edm::EventSetup const &) override
 
void endRun (edm::Run const &run, edm::EventSetup const &eSetup) override
 

Private Member Functions

int getDet (int id)
 
int getPixPlane (int id)
 
int getPlaneIndex (int arm, int station, int rp, int plane)
 
int getRPglobalBin (int arm, int stn)
 
int getRPindex (int arm, int station, int rp)
 
int getRPInStationBin (int rp)
 
int prIndex (int rp, int plane)
 

Private Attributes

int cluSizeMax
 
int ClusMultPlane [RPotsTotalNumber][NplaneMAX]
 
MonitorElementh2AllPlanesActive
 
MonitorElementh2CluSize [NArms][NStationMAX]
 
MonitorElementh2HitsMultipl [NArms][NStationMAX]
 
MonitorElementh2HitsMultROC [RPotsTotalNumber]
 
MonitorElementh2trackXY0 [RPotsTotalNumber]
 
MonitorElementh2xyHits [RPotsTotalNumber][NplaneMAX]
 
MonitorElementh2xyROCHits [RPotsTotalNumber *NplaneMAX][NROCsMAX]
 
MonitorElementhBX
 
MonitorElementhBXshort
 
MonitorElementhHitsMult [RPotsTotalNumber][NplaneMAX]
 
int HitsMultPlane [RPotsTotalNumber][NplaneMAX]
 
int HitsMultROC [RPotsTotalNumber *NplaneMAX][NROCsMAX]
 
MonitorElementhp2HitsMultROC_LS [RPotsTotalNumber]
 
MonitorElementhp2xyADC [RPotsTotalNumber][NplaneMAX]
 
MonitorElementhpixLTrack
 
MonitorElementhpRPactive
 
MonitorElementhROCadc [RPotsTotalNumber *NplaneMAX][NROCsMAX]
 
MonitorElementhRPotActivBX [RPotsTotalNumber]
 
MonitorElementhRPotActivBXall [RPotsTotalNumber]
 
MonitorElementhRPotActivBXroc [RPotsTotalNumber]
 
MonitorElementhRPotActivPlanes [RPotsTotalNumber]
 
MonitorElementhtrackHits [RPotsTotalNumber]
 
MonitorElementhtrackMult [RPotsTotalNumber]
 
const int IndexNotValid = 0
 
int multHitsMax
 
long int nEvents = 0
 
int pixColMAX = 156
 
int pixRowMAX = 160
 
int ROCSizeInX = pixRowMAX/2
 
int ROCSizeInY = pixColMAX/3
 
int RPindexValid [RPotsTotalNumber]
 
int RPstatus [StationIDMAX][RPotsIDMAX]
 
unsigned int rpStatusWord = 0x8008
 
int StationStatus [StationIDMAX]
 
CTPPSPixelIndices thePixIndices
 
edm::EDGetTokenT< edm::DetSetVector< CTPPSPixelCluster > > tokenCluster
 
edm::EDGetTokenT< edm::DetSetVector< CTPPSPixelDigi > > tokenDigi
 
edm::EDGetTokenT< edm::DetSetVector< CTPPSPixelLocalTrack > > tokenTrack
 
int TrackFitDimension = 4
 
unsigned int verbosity
 
float x0_MAX
 
float x0_MIN
 
float y0_MAX
 
float y0_MIN
 

Static Private Attributes

static int ADCMax = 256
 
static int ClusMultMAX = 10
 
static int ClusterSizeMax = 9
 
static int hitMultMAX = 50
 
static int NArms =2
 
static int NLocalTracksMAX =20
 
static int NPlaneBins = NplaneMAX*NRPotBinsInStation
 
static int NplaneMAX =6
 
static int NROCsMAX = 6
 
static int NRPglobalBins = 4
 
static int NRPotBinsInStation = RPn_last-RPn_first
 
static int NRPotsMAX =6
 
static int NStationMAX =3
 
static int RPn_first = 3
 
static int RPn_last = 4
 
static int RPotsIDMAX =8
 
static int RPotsTotalNumber =NArms*NStationMAX*NRPotsMAX
 
static int StationIDMAX =4
 

Detailed Description

Definition at line 33 of file CTPPSPixelDQMSource.cc.

Constructor & Destructor Documentation

CTPPSPixelDQMSource::CTPPSPixelDQMSource ( const edm::ParameterSet ps)

Definition at line 152 of file CTPPSPixelDQMSource.cc.

References edm::ParameterSet::getParameter(), tokenCluster, tokenDigi, and tokenTrack.

152  :
153  verbosity(ps.getUntrackedParameter<unsigned int>("verbosity", 0)),
154  rpStatusWord(ps.getUntrackedParameter<unsigned int>("RPStatusWord",0x8008))
155 {
156  tokenDigi = consumes<DetSetVector<CTPPSPixelDigi> >(ps.getParameter<edm::InputTag>("tagRPixDigi"));
157  tokenCluster=consumes<DetSetVector<CTPPSPixelCluster>>(ps.getParameter<edm::InputTag>("tagRPixCluster"));
158 tokenTrack=consumes<DetSetVector<CTPPSPixelLocalTrack>>(ps.getParameter<edm::InputTag>("tagRPixLTrack"));
159 }
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
edm::EDGetTokenT< edm::DetSetVector< CTPPSPixelDigi > > tokenDigi
edm::EDGetTokenT< edm::DetSetVector< CTPPSPixelCluster > > tokenCluster
edm::EDGetTokenT< edm::DetSetVector< CTPPSPixelLocalTrack > > tokenTrack
CTPPSPixelDQMSource::~CTPPSPixelDQMSource ( )
override

Definition at line 163 of file CTPPSPixelDQMSource.cc.

164 {
165 }

Member Function Documentation

void CTPPSPixelDQMSource::analyze ( edm::Event const &  e,
edm::EventSetup const &  eSetup 
)
overrideprotected

Definition at line 378 of file CTPPSPixelDQMSource.cc.

References ecalMGPA::adc(), CTPPSDetId::arm(), edm::DetSet< T >::begin(), edm::EventBase::bunchCrossing(), cluSizeMax, ClusMultPlane, ClusterSizeMax, cuy::col, MonitorElement::Fill(), objects.autophobj::float, getDet(), getPixPlane(), getPlaneIndex(), getRPglobalBin(), getRPindex(), h2AllPlanesActive, h2CluSize, h2HitsMultipl, h2HitsMultROC, h2trackXY0, h2xyHits, hBX, hBXshort, hHitsMult, hitMultMAX, HitsMultPlane, HitsMultROC, hp2HitsMultROC_LS, hp2xyADC, hpixLTrack, hpRPactive, hRPotActivBX, hRPotActivBXall, hRPotActivBXroc, hRPotActivPlanes, htrackHits, htrackMult, mps_fire::i, edm::HandleBase::isValid(), SiStripPI::max, multHitsMax, nEvents, NLocalTracksMAX, np, NplaneMAX, NROCsMAX, NRPotsMAX, NStationMAX, AlCaHLTBitMon_ParallelJobs::p, prIndex(), alignCSCRings::r, CTPPSDetId::rp(), RPindexValid, RPotsTotalNumber, RPstatus, CTPPSDetId::station(), relativeConstraints::station, StationStatus, thePixIndices, tokenCluster, tokenDigi, tokenTrack, TrackFitDimension, CTPPSPixelIndices::transformToROC(), verbosity, DetId::VeryForward, x0_MAX, x0_MIN, y0_MAX, and y0_MIN.

379 {
380  ++nEvents;
381  int lumiId = event.getLuminosityBlock().id().luminosityBlock();
382  if(lumiId<0) lumiId=0;
383 
384  int RPactivity[RPotsTotalNumber], RPdigiSize[RPotsTotalNumber];
385  int pixRPTracks[RPotsTotalNumber];
386 
387  for(int rp=0; rp<RPotsTotalNumber; rp++)
388  { RPactivity[rp] = RPdigiSize[rp] = pixRPTracks[rp] = 0;}
389 
390  for(int ind=0; ind<RPotsTotalNumber; ind++) {
391  for(int p=0; p<NplaneMAX; p++) {
392  HitsMultPlane[ind][p] = 0;
393  ClusMultPlane[ind][p] = 0;
394  }
395  }
396  for(int ind=0; ind<RPotsTotalNumber*NplaneMAX; ind++) {
397  for(int roc=0; roc<NROCsMAX; roc++) {
398  HitsMultROC[ind][roc] = 0;
399  }
400  }
402  event.getByToken(tokenDigi, pixDigi);
403 
405  event.getByToken(tokenCluster, pixClus);
406 
408  event.getByToken(tokenTrack, pixTrack);
409 
410  hBX->Fill(event.bunchCrossing());
411  hBXshort->Fill(event.bunchCrossing());
412 
413  if(pixTrack.isValid()) {
414  for(const auto &ds_tr : *pixTrack)
415  {
416  int idet = getDet(ds_tr.id);
417  if(idet != DetId::VeryForward) {
418  if(verbosity>1) LogPrint("CTPPSPixelDQMSource") <<"not CTPPS: ds_tr.id"<<ds_tr.id;
419  continue;
420  }
421  CTPPSDetId theId(ds_tr.id);
422  int arm = theId.arm()&0x1;
423  int station = theId.station()&0x3;
424  int rpot = theId.rp()&0x7;
425  int rpInd = getRPindex(arm,station,rpot);
426 
428  dit != ds_tr.end(); ++dit) {
429  ++pixRPTracks[rpInd];
430  int nh_tr = (dit->getNDF() + TrackFitDimension)/2;
431  for(int i=0; i<=NplaneMAX; i++) {
432  if(i==nh_tr) htrackHits[rpInd]->Fill(nh_tr,1.);
433  else htrackHits[rpInd]->Fill(i,0.);
434  }
435  float x0 = dit->getX0();
436  float y0 = dit->getY0();
437  h2trackXY0[rpInd]->Fill(x0,y0);
438  if(x0_MAX < x0) x0_MAX = x0;
439  if(y0_MAX < y0) y0_MAX = y0;
440  if(x0_MIN > x0) x0_MIN = x0;
441  if(y0_MIN > y0) y0_MIN = y0;
442  }
443  }
444  } // end if(pixTrack.isValid())
445 
446 
447  bool valid = false;
448  valid |= pixDigi.isValid();
449 // valid |= Clus.isValid();
450 
451  if(!valid && verbosity) LogPrint("CTPPSPixelDQMSource") <<"No valid data in Event "<<nEvents;
452 
453  if(pixDigi.isValid()) {
454  for(const auto &ds_digi : *pixDigi)
455  {
456  int idet = getDet(ds_digi.id);
457  if(idet != DetId::VeryForward) {
458  if(verbosity>1) LogPrint("CTPPSPixelDQMSource") <<"not CTPPS: ds_digi.id"<<ds_digi.id;
459  continue;
460  }
461  // int subdet = getSubdet(ds_digi.id);
462 
463  int plane = getPixPlane(ds_digi.id);
464 
465  CTPPSDetId theId(ds_digi.id);
466  int arm = theId.arm()&0x1;
467  int station = theId.station()&0x3;
468  int rpot = theId.rp()&0x7;
469  int rpInd = getRPindex(arm,station,rpot);
470  RPactivity[rpInd] = 1;
471  ++RPdigiSize[rpInd];
472 
473  if(StationStatus[station] && RPstatus[station][rpot]) {
474 
475  h2HitsMultipl[arm][station]->Fill(prIndex(rpot,plane),ds_digi.data.size());
476  h2AllPlanesActive->Fill(plane,getRPglobalBin(arm,station));
477 
478  int index = getRPindex(arm,station,rpot);
479  HitsMultPlane[index][plane] += ds_digi.data.size();
480  if(RPindexValid[index]) {
481  int nh = ds_digi.data.size();
482  if(nh > hitMultMAX) nh = hitMultMAX;
483  hHitsMult[index][plane]->Fill(nh);
484  }
485  int rocHistIndex = getPlaneIndex(arm,station,rpot,plane);
486 
487  for(DetSet<CTPPSPixelDigi>::const_iterator dit = ds_digi.begin();
488  dit != ds_digi.end(); ++dit) {
489  int row = dit->row();
490  int col = dit->column();
491  int adc = dit->adc();
492 
493  if(RPindexValid[index]) {
494  h2xyHits[index][plane]->Fill(col,row);
495  hp2xyADC[index][plane]->Fill(col,row,adc);
496  int colROC, rowROC;
497  int trocId;
498  if(!thePixIndices.transformToROC(col,row, trocId, colROC, rowROC)) {
499  if(trocId>=0 && trocId<NROCsMAX) {
500  ++HitsMultROC[rocHistIndex][trocId];
501  }
502  }
503  } //end if(RPindexValid[index]) {
504  }
505  if(int(ds_digi.data.size()) > multHitsMax) multHitsMax = ds_digi.data.size();
506 
507  } // end if(StationStatus[station]) {
508  } // end for(const auto &ds_digi : *pixDigi)
509  } //if(pixDigi.isValid()) {
510 
511  if(pixClus.isValid())
512  for(const auto &ds : *pixClus)
513  {
514  int idet = getDet(ds.id);
515  if(idet != DetId::VeryForward) {
516  if(verbosity>1) LogPrint("CTPPSPixelDQMSource") <<"not CTPPS: cluster.id"<<ds.id;
517  continue;
518  }
519  // int subdet = getSubdet(ds.id);
520 
521  int plane = getPixPlane(ds.id);
522 
523  CTPPSDetId theId(ds.id);
524  int arm = theId.arm()&0x1;
525  int station = theId.station()&0x3;
526  int rpot = theId.rp()&0x7;
527 
528  if((StationStatus[station]==0) || (RPstatus[station][rpot]==0)) continue;
529 
530  int index = getRPindex(arm,station,rpot);
531  ++ClusMultPlane[index][plane];
532 
533  for (const auto &p : ds) {
534  int clusize = p.size();
535  h2CluSize[arm][station]->Fill(prIndex(rpot,plane),clusize);
536  if(cluSizeMax < clusize) cluSizeMax = clusize;
537  if(clusize > ClusterSizeMax) clusize = ClusterSizeMax;
538  }
539  } // end if(pixClus.isValid()) for(const auto &ds : *pixClus)
540 
541  bool allRPactivity = false;
542  for(int rp=0; rp<RPotsTotalNumber; rp++) if(RPactivity[rp]>0) allRPactivity=true;
543  for(int arm=0; arm<2; arm++) {
544  for(int stn=0; stn<NStationMAX; stn++) {
545  for(int rp=0; rp<NRPotsMAX; rp++) {
546  int index = getRPindex(arm,stn,rp);
547  if(RPindexValid[index]==0) continue;
548  hpRPactive->Fill(getRPglobalBin(arm,stn),RPactivity[index]);
549 // if(RPactivity[index]==0) continue;
550  if(!allRPactivity) continue;
551  hpixLTrack->Fill(getRPglobalBin(arm,stn),pixRPTracks[index]);
552  int ntr = pixRPTracks[index];
553  if(ntr > NLocalTracksMAX) ntr = NLocalTracksMAX;
554  for(int i=0; i<=NLocalTracksMAX; i++) {
555  if(i==ntr) htrackMult[index]->Fill(ntr,1.);
556  else htrackMult[index]->Fill(i,0.);
557  }
558 
559  int np = 0;
560  for(int p=0; p<NplaneMAX; p++) if(HitsMultPlane[index][p]>0) np++;
561  for(int p=0; p<=NplaneMAX; p++) {
562  if(p == np) hRPotActivPlanes[index]->Fill(p,1.);
563  else hRPotActivPlanes[index]->Fill(p,0.);
564  }
565  if(np>=5) hRPotActivBX[index]->Fill(event.bunchCrossing());
566  hRPotActivBXall[index]->Fill(event.bunchCrossing(),float(RPdigiSize[index]));
567 
568  int rocf[NplaneMAX];
569  for(int r=0; r<NROCsMAX; r++) rocf[r]=0;
570  for(int p=0; p<NplaneMAX; p++) {
571  int indp = getPlaneIndex(arm,stn,rp,p);
572  for(int r=0; r<NROCsMAX; r++) if(HitsMultROC[indp][r] > 0) ++rocf[r];
573  for(int r=0; r<NROCsMAX; r++) {
574  h2HitsMultROC[index]->Fill(p,r,HitsMultROC[indp][r]);
575  hp2HitsMultROC_LS[index]->Fill(lumiId,p*NROCsMAX+r,HitsMultROC[indp][r]);
576  }
577  }
578  int max = 0;
579  for(int r=0; r<NROCsMAX; r++)
580  if(max < rocf[r]) { max = rocf[r]; }
581  if(max >= 4) hRPotActivBXroc[index]->Fill(event.bunchCrossing());
582  } //end for(int rp=0; rp<NRPotsMAX; rp++) {
583  }
584  } //end for(int arm=0; arm<2; arm++) {
585 
586  if((nEvents % 100)) return;
587  if(verbosity) LogPrint("CTPPSPixelDQMSource")<<"analyze event "<<nEvents;
588 }
int adc(sample_type sample)
get the ADC sample (12 bits)
int StationStatus[StationIDMAX]
MonitorElement * hHitsMult[RPotsTotalNumber][NplaneMAX]
int HitsMultPlane[RPotsTotalNumber][NplaneMAX]
MonitorElement * hp2xyADC[RPotsTotalNumber][NplaneMAX]
MonitorElement * hpRPactive
MonitorElement * hp2HitsMultROC_LS[RPotsTotalNumber]
MonitorElement * h2CluSize[NArms][NStationMAX]
edm::EDGetTokenT< edm::DetSetVector< CTPPSPixelDigi > > tokenDigi
MonitorElement * htrackMult[RPotsTotalNumber]
MonitorElement * h2HitsMultROC[RPotsTotalNumber]
MonitorElement * h2trackXY0[RPotsTotalNumber]
MonitorElement * hpixLTrack
CTPPSPixelIndices thePixIndices
edm::EDGetTokenT< edm::DetSetVector< CTPPSPixelCluster > > tokenCluster
void Fill(long long x)
edm::EDGetTokenT< edm::DetSetVector< CTPPSPixelLocalTrack > > tokenTrack
MonitorElement * htrackHits[RPotsTotalNumber]
int prIndex(int rp, int plane)
int RPindexValid[RPotsTotalNumber]
int np
Definition: AMPTWrapper.h:33
MonitorElement * hRPotActivBX[RPotsTotalNumber]
bool isValid() const
Definition: HandleBase.h:74
int getRPindex(int arm, int station, int rp)
MonitorElement * hRPotActivBXroc[RPotsTotalNumber]
iterator begin()
Definition: DetSet.h:59
int ClusMultPlane[RPotsTotalNumber][NplaneMAX]
MonitorElement * hRPotActivBXall[RPotsTotalNumber]
Base class for CTPPS detector IDs.
Definition: CTPPSDetId.h:32
int transformToROC(const int col, const int row, int &rocId, int &colROC, int &rowROC) const
col
Definition: cuy.py:1009
MonitorElement * h2AllPlanesActive
int getRPglobalBin(int arm, int stn)
int HitsMultROC[RPotsTotalNumber *NplaneMAX][NROCsMAX]
int getPlaneIndex(int arm, int station, int rp, int plane)
MonitorElement * h2xyHits[RPotsTotalNumber][NplaneMAX]
MonitorElement * hRPotActivPlanes[RPotsTotalNumber]
int RPstatus[StationIDMAX][RPotsIDMAX]
MonitorElement * hBXshort
MonitorElement * h2HitsMultipl[NArms][NStationMAX]
Definition: event.py:1
void CTPPSPixelDQMSource::bookHistograms ( DQMStore::IBooker ibooker,
edm::Run const &  ,
edm::EventSetup const &   
)
overrideprotected

Definition at line 206 of file CTPPSPixelDQMSource.cc.

References CTPPSDetId::armName(), DQMStore::IBooker::book1D(), DQMStore::IBooker::book1DD(), DQMStore::IBooker::book2D(), DQMStore::IBooker::book2DD(), DQMStore::IBooker::bookProfile(), DQMStore::IBooker::bookProfile2D(), DQMStore::IBooker::cd(), ClusterSizeMax, MillePedeFileConverter_cfg::e, getRPglobalBin(), CTPPSDetId::getRPId(), getRPindex(), CTPPSDetId::getStationId(), MonitorElement::getTH2D(), MonitorElement::getTH2F(), MonitorElement::getTProfile(), MonitorElement::getTProfile2D(), h2AllPlanesActive, h2CluSize, h2HitsMultipl, h2HitsMultROC, h2trackXY0, h2xyHits, hBX, hBXshort, hHitsMult, hitMultMAX, hp2HitsMultROC_LS, hp2xyADC, hpixLTrack, hpRPactive, hRPotActivBX, hRPotActivBXall, hRPotActivBXroc, hRPotActivPlanes, htrackHits, htrackMult, pileupCalc::nbins, CTPPSDetId::nFull, NLocalTracksMAX, CTPPSDetId::nPath, NPlaneBins, NplaneMAX, NROCsMAX, NRPglobalBins, CTPPSDetId::nShort, NStationMAX, AlCaHLTBitMon_ParallelJobs::p, pixRowMAX, alignCSCRings::r, ROCSizeInX, ROCSizeInY, RPindexValid, RPn_first, RPn_last, RPstatus, alignCSCRings::s, sd, CTPPSDetId::sdTrackingPixel, DQMStore::IBooker::setCurrentFolder(), CTPPSDetId::setRP(), CTPPSDetId::setStation(), StationStatus, and AlCaHLTBitMon_QueryRunRegistry::string.

208 {
209  ibooker.cd();
210  ibooker.setCurrentFolder("CTPPS/TrackingPixel");
211  char s[50];
212  string armTitleShort, stnTitleShort;
213  hBX = ibooker.book1D("events per BX", "ctpps_pixel;Event.BX", 4002, -1.5, 4000. +0.5);
214  hBXshort = ibooker.book1D("events per BX(short)", "ctpps_pixel;Event.BX", 102, -1.5, 100. + 0.5);
215 
216  string str1st = "Pixel planes activity";
217  h2AllPlanesActive = ibooker.book2DD(str1st,str1st+"(digi task);Plane #",
219  TH2D *h1st = h2AllPlanesActive->getTH2D();
220  h1st->SetOption("colz");
221  TAxis *yah1st = h1st->GetYaxis();
222 
223  string str2 = "Pixel RP active";
224  hpRPactive = ibooker.bookProfile(str2,str2+" per event(digi task)",
225  NRPglobalBins, 0.5,NRPglobalBins+0.5, -0.1, 1.1,"" );
226  TAxis *xaRPact = hpRPactive->getTProfile()->GetXaxis();
227  hpRPactive->getTProfile()->SetOption("hist");
228  hpRPactive->getTProfile()->SetMinimum(0.);
229  hpRPactive->getTProfile()->SetMaximum(1.1);
230 
231  str2 = "Pixel Local Tracks";
232  hpixLTrack = ibooker.bookProfile(str2,str2+" per event",
233  NRPglobalBins, 0.5, NRPglobalBins+0.5,-0.1,NLocalTracksMAX,"");
234 
235  TAxis *xah1trk = hpixLTrack->getTProfile()->GetXaxis();
236  hpixLTrack->getTProfile()->GetYaxis()->SetTitle("average number of tracks per event");
237  hpixLTrack->getTProfile()->SetOption("hist");
238 
239  for(int arm=0; arm<2; arm++) {
241  string sd, armTitle;
242  ID.armName(sd, CTPPSDetId::nPath);
243  ID.armName(armTitle, CTPPSDetId::nFull);
244  ID.armName(armTitleShort, CTPPSDetId::nShort);
245 
246  ibooker.setCurrentFolder(sd);
247 
248  for(int stn=0; stn<NStationMAX; stn++) {
249  if(StationStatus[stn]==0) continue;
250  ID.setStation(stn);
251  string stnd, stnTitle;
252 
253  CTPPSDetId(ID.getStationId()).stationName(stnd, CTPPSDetId::nPath);
254  CTPPSDetId(ID.getStationId()).stationName(stnTitle, CTPPSDetId::nFull);
255  CTPPSDetId(ID.getStationId()).stationName(stnTitleShort, CTPPSDetId::nShort);
256 
257  ibooker.setCurrentFolder(stnd);
258 //--------- RPots ---
259  int pixBinW = 4;
260  for(int rp=RPn_first; rp<RPn_last; rp++) { // only installed pixel pots
261  ID.setRP(rp);
262  string rpd, rpTitle;
263  CTPPSDetId(ID.getRPId()).rpName(rpTitle, CTPPSDetId::nShort);
264  string rpBinName = armTitleShort + "_" + stnTitleShort+"_"+rpTitle;
265  yah1st->SetBinLabel(getRPglobalBin(arm,stn), rpBinName.c_str());
266  xah1trk->SetBinLabel(getRPglobalBin(arm,stn), rpBinName.c_str());
267  xaRPact->SetBinLabel(getRPglobalBin(arm,stn), rpBinName.c_str());
268 
269  if(RPstatus[stn][rp]==0) continue;
270  int indexP = getRPindex(arm,stn,rp);
271  RPindexValid[indexP] = 1;
272 
273  CTPPSDetId(ID.getRPId()).rpName(rpTitle, CTPPSDetId::nFull);
274  CTPPSDetId(ID.getRPId()).rpName(rpd, CTPPSDetId::nPath);
275 
276  ibooker.setCurrentFolder(rpd);
277 
278  string st2 = ": " + stnTitle;
279 
280  string st = "hit multiplicity in planes";
281  string st3 = ";PlaneIndex(=pixelPot*PlaneMAX + plane)";
282  h2HitsMultipl[arm][stn]= ibooker.book2DD(st,st+st2+st3+";multiplicity",
284  h2HitsMultipl[arm][stn]->getTH2D()->SetOption("colz");
285 
286  st = "cluster size in planes";
287  h2CluSize[arm][stn] = ibooker.book2D(st,st+st2+st3+";Cluster size",
289  h2CluSize[arm][stn]->getTH2F()->SetOption("colz");
290 
291  const float x0Maximum = 70.;
292  const float y0Maximum = 15.;
293  st = "track intercept point";
294  h2trackXY0[indexP] = ibooker.book2D(st,st+st2+";x0;y0",
295  int(x0Maximum)*2, 0.,x0Maximum,int(y0Maximum)*4,-y0Maximum,y0Maximum);
296  h2trackXY0[indexP]->getTH2F()->SetOption("colz");
297 
298  st = "number of tracks per event";
299  htrackMult[indexP] = ibooker.bookProfile(st,rpTitle+";number of tracks",
300  NLocalTracksMAX+1, -0.5,NLocalTracksMAX+0.5, -0.5,NLocalTracksMAX+0.5,"");
301  htrackMult[indexP]->getTProfile()->SetOption("hist");
302 
303  st = "number of hits per track";
304  htrackHits[indexP] = ibooker.bookProfile(st,rpTitle+";number of hits",
305  5,1.5,6.5, -0.1,1.1, "");
306  htrackHits[indexP]->getTProfile()->SetOption("hist");
307 
308  hRPotActivPlanes[indexP] =
309  ibooker.bookProfile("number of fired planes per event", rpTitle+";nPlanes;Probability",
310  NplaneMAX+1, -0.5, NplaneMAX+0.5, -0.5,NplaneMAX+0.5,"");
311  hRPotActivPlanes[indexP]->getTProfile()->SetOption("hist");
312 
313  h2HitsMultROC[indexP] = ibooker.bookProfile2D("ROCs hits multiplicity per event",
314  rpTitle+";plane # ;ROC #", NplaneMAX,-0.5,NplaneMAX-0.5, NROCsMAX,-0.5,NROCsMAX-0.5,
315  0.,ROCSizeInX*ROCSizeInY,"");
316  h2HitsMultROC[indexP]->getTProfile2D()->SetOption("colztext");
317  h2HitsMultROC[indexP]->getTProfile2D()->SetMinimum(1.e-10);
318 
319 
320  hp2HitsMultROC_LS[indexP]=ibooker.bookProfile2D("ROCs_hits_multiplicity_per_event vs LS",
321  rpTitle+";LumiSection;Plane#___ROC#", 1000,0.,1000.,
322  NplaneMAX*NROCsMAX,0.,double(NplaneMAX*NROCsMAX),0.,ROCSizeInX*ROCSizeInY,"");
323  hp2HitsMultROC_LS[indexP]->getTProfile2D()->SetOption("colz");
324  hp2HitsMultROC_LS[indexP]->getTProfile2D()->SetMinimum(1.0e-10);
325  hp2HitsMultROC_LS[indexP]->getTProfile2D()->SetCanExtend(TProfile2D::kXaxis);
326  TAxis *yahp2 = hp2HitsMultROC_LS[indexP]->getTProfile2D()->GetYaxis();
327  for(int p=0; p<NplaneMAX; p++) {
328  sprintf(s,"plane%d_0",p);
329  yahp2->SetBinLabel(p*NplaneMAX+1,s);
330  for(int r=1; r<NROCsMAX; r++) {
331  sprintf(s," %d_%d",p,r);
332  yahp2->SetBinLabel(p*NplaneMAX+r+1,s);
333  }
334  }
335 
336  ibooker.setCurrentFolder(rpd+"/latency");
337  hRPotActivBX[indexP] =
338  ibooker.book1D("5 fired planes per BX", rpTitle+";Event.BX",4002,-1.5, 4000.+0.5);
339 
340  hRPotActivBXroc[indexP] =
341  ibooker.book1D("4 fired ROCs per BX", rpTitle+";Event.BX", 4002, -1.5, 4000.+0.5);
342 
343  hRPotActivBXall[indexP] =
344  ibooker.book1D("hits per BX", rpTitle+";Event.BX", 4002, -1.5, 4000.+0.5);
345 
346  int nbins = pixRowMAX/pixBinW;
347 
348  for(int p=0; p<NplaneMAX; p++) {
349  sprintf(s,"plane_%d",p);
350  string pd = rpd+"/"+string(s);
351  ibooker.setCurrentFolder(pd);
352  string st1 = ": "+rpTitle+"_"+string(s);
353 
354  st = "hits position";
355  h2xyHits[indexP][p]=ibooker.book2DD(st,st1+";pix col;pix row",
357  h2xyHits[indexP][p]->getTH2D()->SetOption("colz");
358 
359  st = "adc average value";
360  hp2xyADC[indexP][p]=ibooker.bookProfile2D(st,st1+";pix col;pix row",
361  nbins,0,pixRowMAX,nbins,0,pixRowMAX, 0.,512.,"");
362  hp2xyADC[indexP][p]->getTProfile2D()->SetOption("colz");
363 
364  st = "hits multiplicity";
365  hHitsMult[indexP][p]=ibooker.book1DD(st,st1+";number of hits;N / 1 hit",
366  hitMultMAX+1,-0.5,hitMultMAX+0.5);
367  } // end of for(int p=0; p<NplaneMAX;..
368 
369  } // end for(int rp=0; rp<NRPotsMAX;...
370  } // end of for(int stn=0; stn<
371  } //end of for(int arm=0; arm<2;...
372 
373  return;
374 }
TProfile * getTProfile() const
int StationStatus[StationIDMAX]
MonitorElement * hHitsMult[RPotsTotalNumber][NplaneMAX]
MonitorElement * hp2xyADC[RPotsTotalNumber][NplaneMAX]
MonitorElement * hpRPactive
MonitorElement * hp2HitsMultROC_LS[RPotsTotalNumber]
MonitorElement * h2CluSize[NArms][NStationMAX]
TProfile2D * getTProfile2D() const
MonitorElement * bookProfile(Args &&...args)
Definition: DQMStore.h:160
MonitorElement * htrackMult[RPotsTotalNumber]
MonitorElement * h2HitsMultROC[RPotsTotalNumber]
uint32_t ID
Definition: Definitions.h:26
MonitorElement * h2trackXY0[RPotsTotalNumber]
MonitorElement * hpixLTrack
MonitorElement * htrackHits[RPotsTotalNumber]
int RPindexValid[RPotsTotalNumber]
MonitorElement * bookProfile2D(Args &&...args)
Definition: DQMStore.h:166
MonitorElement * hRPotActivBX[RPotsTotalNumber]
TH2D * getTH2D() const
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:118
int getRPindex(int arm, int station, int rp)
TH2F * getTH2F() const
MonitorElement * hRPotActivBXroc[RPotsTotalNumber]
MonitorElement * hRPotActivBXall[RPotsTotalNumber]
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:274
MonitorElement * book2D(Args &&...args)
Definition: DQMStore.h:136
double sd
Base class for CTPPS detector IDs.
Definition: CTPPSDetId.h:32
MonitorElement * h2AllPlanesActive
MonitorElement * book2DD(Args &&...args)
Definition: DQMStore.h:148
int getRPglobalBin(int arm, int stn)
MonitorElement * h2xyHits[RPotsTotalNumber][NplaneMAX]
MonitorElement * hRPotActivPlanes[RPotsTotalNumber]
MonitorElement * book1DD(Args &&...args)
Definition: DQMStore.h:130
int RPstatus[StationIDMAX][RPotsIDMAX]
MonitorElement * hBXshort
MonitorElement * h2HitsMultipl[NArms][NStationMAX]
void CTPPSPixelDQMSource::dqmBeginRun ( edm::Run const &  run,
edm::EventSetup const &   
)
overrideprotected

Definition at line 169 of file CTPPSPixelDQMSource.cc.

References cluSizeMax, CTPPSPixelLocalTrack::dimension, CTPPSPixelIndices::getDefaultColDetSize(), CTPPSPixelIndices::getDefaultRowDetSize(), multHitsMax, nEvents, NRPotsMAX, NStationMAX, pixColMAX, pixRowMAX, ROCSizeInX, ROCSizeInY, RPindexValid, RPotsIDMAX, RPstatus, rpStatusWord, StationIDMAX, StationStatus, thePixIndices, TrackFitDimension, verbosity, x0_MAX, x0_MIN, y0_MAX, and y0_MIN.

170 {
171  if(verbosity) LogPrint("CTPPSPixelDQMSource") <<"RPstatusWord= "<<rpStatusWord;
172  nEvents = 0;
173 
174  CTPPSPixelLocalTrack thePixelLocalTrack;
175  TrackFitDimension = thePixelLocalTrack.dimension;
176 
179  ROCSizeInX = pixRowMAX/2; // ROC row size in pixels = 80
180  ROCSizeInY = pixColMAX/3;
181 
182  for(int stn=0; stn<StationIDMAX; stn++) {
183  StationStatus[stn]=0;
184  for(int rp=0; rp<RPotsIDMAX; rp++) RPstatus[stn][rp]=0;
185  }
186 
187  unsigned int rpSts = rpStatusWord<<1;
188  for(int stn=0; stn<NStationMAX; stn++) {
189  int stns = 0;
190  for(int rp=0; rp<NRPotsMAX; rp++) {
191  rpSts = (rpSts >> 1); RPstatus[stn][rp] = rpSts&1;
192  if(RPstatus[stn][rp] > 0) stns = 1;
193  }
194  StationStatus[stn]=stns;
195  }
196 
197  for(int ind=0; ind<2*3*NRPotsMAX; ind++) RPindexValid[ind] = 0;
198 
199  multHitsMax = cluSizeMax = -1;
200  x0_MIN = y0_MIN = 1.0e06;
201  x0_MAX = y0_MAX = -1.0e06;
202 }
int StationStatus[StationIDMAX]
int getDefaultRowDetSize() const
CTPPSPixelIndices thePixIndices
int RPindexValid[RPotsTotalNumber]
int getDefaultColDetSize() const
int RPstatus[StationIDMAX][RPotsIDMAX]
void CTPPSPixelDQMSource::endRun ( edm::Run const &  run,
edm::EventSetup const &  eSetup 
)
overrideprotected

Definition at line 591 of file CTPPSPixelDQMSource.cc.

References cluSizeMax, DEFINE_FWK_MODULE, multHitsMax, nEvents, edm::RunBase::run(), verbosity, x0_MAX, x0_MIN, y0_MAX, and y0_MIN.

592 {
593  if(!verbosity) return;
594  LogPrint("CTPPSPixelDQMSource")
595  <<"end of Run "<<run.run()<<": "<<nEvents<<" events\n"
596  <<"multHitsMax= "<<multHitsMax <<" cluSizeMax= "<<cluSizeMax
597  <<"\nx0: "<< x0_MIN<<"/" << x0_MAX<<"y0: "<<y0_MIN<<"/"<<y0_MAX<<"\n";
598 }
int CTPPSPixelDQMSource::getDet ( int  id)
inlineprivate

Definition at line 135 of file CTPPSPixelDQMSource.cc.

References DetId::kDetOffset.

Referenced by analyze().

136  { return (id>>DetId::kDetOffset)&0xF; }
static const int kDetOffset
Definition: DetId.h:22
int CTPPSPixelDQMSource::getPixPlane ( int  id)
inlineprivate

Definition at line 137 of file CTPPSPixelDQMSource.cc.

Referenced by analyze().

138  { return ((id>>16)&0x7); }
int CTPPSPixelDQMSource::getPlaneIndex ( int  arm,
int  station,
int  rp,
int  plane 
)
inlineprivate

Definition at line 117 of file CTPPSPixelDQMSource.cc.

References getRPindex().

Referenced by analyze().

117  {
118  if(plane<0 || plane>=NplaneMAX) return(IndexNotValid);
119  int rc = getRPindex(arm, station, rp);
120  if(rc == IndexNotValid) return(IndexNotValid);
121  return(rc*NplaneMAX + plane);
122  }
int getRPindex(int arm, int station, int rp)
int CTPPSPixelDQMSource::getRPglobalBin ( int  arm,
int  stn 
)
inlineprivate

Definition at line 128 of file CTPPSPixelDQMSource.cc.

References constexpr, and NStationMAX.

Referenced by analyze(), and bookHistograms().

128  {
129  static constexpr int stationBinOrder[NStationMAX] = {0, 4, 1};
130  return( arm*2 + stationBinOrder[stn] +1 );
131  }
#define constexpr
int CTPPSPixelDQMSource::getRPindex ( int  arm,
int  station,
int  rp 
)
inlineprivate

Definition at line 110 of file CTPPSPixelDQMSource.cc.

References relativeConstraints::station.

Referenced by analyze(), bookHistograms(), and getPlaneIndex().

110  {
111  if(arm<0 || station<0 || rp<0) return(IndexNotValid);
112  if(arm>1 || station>=NStationMAX || rp>=NRPotsMAX) return(IndexNotValid);
113  int rc = (arm*NStationMAX+station)*NRPotsMAX + rp;
114  return(rc);
115  }
int CTPPSPixelDQMSource::getRPInStationBin ( int  rp)
inlineprivate

Definition at line 124 of file CTPPSPixelDQMSource.cc.

124 { return(rp - RPn_first +1); }
int CTPPSPixelDQMSource::prIndex ( int  rp,
int  plane 
)
inlineprivate

Definition at line 133 of file CTPPSPixelDQMSource.cc.

Referenced by analyze().

134  {return((rp - RPn_first)*NplaneMAX + plane);}

Member Data Documentation

int CTPPSPixelDQMSource::ADCMax = 256
staticprivate

Definition at line 58 of file CTPPSPixelDQMSource.cc.

int CTPPSPixelDQMSource::cluSizeMax
private

Definition at line 141 of file CTPPSPixelDQMSource.cc.

Referenced by analyze(), dqmBeginRun(), and endRun().

int CTPPSPixelDQMSource::ClusMultMAX = 10
staticprivate

Definition at line 63 of file CTPPSPixelDQMSource.cc.

int CTPPSPixelDQMSource::ClusMultPlane[RPotsTotalNumber][NplaneMAX]
private

Definition at line 102 of file CTPPSPixelDQMSource.cc.

Referenced by analyze().

int CTPPSPixelDQMSource::ClusterSizeMax = 9
staticprivate

Definition at line 64 of file CTPPSPixelDQMSource.cc.

Referenced by analyze(), and bookHistograms().

MonitorElement * CTPPSPixelDQMSource::h2AllPlanesActive
private

Definition at line 77 of file CTPPSPixelDQMSource.cc.

Referenced by analyze(), and bookHistograms().

MonitorElement* CTPPSPixelDQMSource::h2CluSize[NArms][NStationMAX]
private

Definition at line 81 of file CTPPSPixelDQMSource.cc.

Referenced by analyze(), and bookHistograms().

MonitorElement* CTPPSPixelDQMSource::h2HitsMultipl[NArms][NStationMAX]
private

Definition at line 80 of file CTPPSPixelDQMSource.cc.

Referenced by analyze(), and bookHistograms().

MonitorElement* CTPPSPixelDQMSource::h2HitsMultROC[RPotsTotalNumber]
private

Definition at line 92 of file CTPPSPixelDQMSource.cc.

Referenced by analyze(), and bookHistograms().

MonitorElement* CTPPSPixelDQMSource::h2trackXY0[RPotsTotalNumber]
private

Definition at line 86 of file CTPPSPixelDQMSource.cc.

Referenced by analyze(), and bookHistograms().

MonitorElement* CTPPSPixelDQMSource::h2xyHits[RPotsTotalNumber][NplaneMAX]
private

Definition at line 95 of file CTPPSPixelDQMSource.cc.

Referenced by analyze(), and bookHistograms().

MonitorElement* CTPPSPixelDQMSource::h2xyROCHits[RPotsTotalNumber *NplaneMAX][NROCsMAX]
private

Definition at line 97 of file CTPPSPixelDQMSource.cc.

MonitorElement* CTPPSPixelDQMSource::hBX
private

Definition at line 77 of file CTPPSPixelDQMSource.cc.

Referenced by analyze(), and bookHistograms().

MonitorElement * CTPPSPixelDQMSource::hBXshort
private

Definition at line 77 of file CTPPSPixelDQMSource.cc.

Referenced by analyze(), and bookHistograms().

MonitorElement* CTPPSPixelDQMSource::hHitsMult[RPotsTotalNumber][NplaneMAX]
private

Definition at line 94 of file CTPPSPixelDQMSource.cc.

Referenced by analyze(), and bookHistograms().

int CTPPSPixelDQMSource::hitMultMAX = 50
staticprivate

Definition at line 62 of file CTPPSPixelDQMSource.cc.

Referenced by analyze(), and bookHistograms().

int CTPPSPixelDQMSource::HitsMultPlane[RPotsTotalNumber][NplaneMAX]
private

Definition at line 101 of file CTPPSPixelDQMSource.cc.

Referenced by analyze().

int CTPPSPixelDQMSource::HitsMultROC[RPotsTotalNumber *NplaneMAX][NROCsMAX]
private

Definition at line 100 of file CTPPSPixelDQMSource.cc.

Referenced by analyze().

MonitorElement* CTPPSPixelDQMSource::hp2HitsMultROC_LS[RPotsTotalNumber]
private

Definition at line 93 of file CTPPSPixelDQMSource.cc.

Referenced by analyze(), and bookHistograms().

MonitorElement* CTPPSPixelDQMSource::hp2xyADC[RPotsTotalNumber][NplaneMAX]
private

Definition at line 96 of file CTPPSPixelDQMSource.cc.

Referenced by analyze(), and bookHistograms().

MonitorElement * CTPPSPixelDQMSource::hpixLTrack
private

Definition at line 77 of file CTPPSPixelDQMSource.cc.

Referenced by analyze(), and bookHistograms().

MonitorElement* CTPPSPixelDQMSource::hpRPactive
private

Definition at line 78 of file CTPPSPixelDQMSource.cc.

Referenced by analyze(), and bookHistograms().

MonitorElement* CTPPSPixelDQMSource::hROCadc[RPotsTotalNumber *NplaneMAX][NROCsMAX]
private

Definition at line 98 of file CTPPSPixelDQMSource.cc.

MonitorElement* CTPPSPixelDQMSource::hRPotActivBX[RPotsTotalNumber]
private

Definition at line 90 of file CTPPSPixelDQMSource.cc.

Referenced by analyze(), and bookHistograms().

MonitorElement* CTPPSPixelDQMSource::hRPotActivBXall[RPotsTotalNumber]
private

Definition at line 99 of file CTPPSPixelDQMSource.cc.

Referenced by analyze(), and bookHistograms().

MonitorElement* CTPPSPixelDQMSource::hRPotActivBXroc[RPotsTotalNumber]
private

Definition at line 91 of file CTPPSPixelDQMSource.cc.

Referenced by analyze(), and bookHistograms().

MonitorElement* CTPPSPixelDQMSource::hRPotActivPlanes[RPotsTotalNumber]
private

Definition at line 89 of file CTPPSPixelDQMSource.cc.

Referenced by analyze(), and bookHistograms().

MonitorElement* CTPPSPixelDQMSource::htrackHits[RPotsTotalNumber]
private

Definition at line 88 of file CTPPSPixelDQMSource.cc.

Referenced by analyze(), and bookHistograms().

MonitorElement* CTPPSPixelDQMSource::htrackMult[RPotsTotalNumber]
private

Definition at line 87 of file CTPPSPixelDQMSource.cc.

Referenced by analyze(), and bookHistograms().

const int CTPPSPixelDQMSource::IndexNotValid = 0
private

Definition at line 108 of file CTPPSPixelDQMSource.cc.

int CTPPSPixelDQMSource::multHitsMax
private

Definition at line 141 of file CTPPSPixelDQMSource.cc.

Referenced by analyze(), dqmBeginRun(), and endRun().

int CTPPSPixelDQMSource::NArms =2
staticprivate

Definition at line 52 of file CTPPSPixelDQMSource.cc.

long int CTPPSPixelDQMSource::nEvents = 0
private

Definition at line 47 of file CTPPSPixelDQMSource.cc.

Referenced by analyze(), dqmBeginRun(), endRun(), and looper.Looper::loop().

int CTPPSPixelDQMSource::NLocalTracksMAX =20
staticprivate

Definition at line 61 of file CTPPSPixelDQMSource.cc.

Referenced by analyze(), and bookHistograms().

int CTPPSPixelDQMSource::NPlaneBins = NplaneMAX*NRPotBinsInStation
staticprivate

Definition at line 75 of file CTPPSPixelDQMSource.cc.

Referenced by bookHistograms().

int CTPPSPixelDQMSource::NplaneMAX =6
staticprivate

Definition at line 55 of file CTPPSPixelDQMSource.cc.

Referenced by analyze(), and bookHistograms().

int CTPPSPixelDQMSource::NROCsMAX = 6
staticprivate

Definition at line 56 of file CTPPSPixelDQMSource.cc.

Referenced by analyze(), and bookHistograms().

int CTPPSPixelDQMSource::NRPglobalBins = 4
staticprivate

Definition at line 127 of file CTPPSPixelDQMSource.cc.

Referenced by bookHistograms().

int CTPPSPixelDQMSource::NRPotBinsInStation = RPn_last-RPn_first
staticprivate

Definition at line 74 of file CTPPSPixelDQMSource.cc.

int CTPPSPixelDQMSource::NRPotsMAX =6
staticprivate

Definition at line 54 of file CTPPSPixelDQMSource.cc.

Referenced by analyze(), and dqmBeginRun().

int CTPPSPixelDQMSource::NStationMAX =3
staticprivate

Definition at line 53 of file CTPPSPixelDQMSource.cc.

Referenced by analyze(), bookHistograms(), dqmBeginRun(), and getRPglobalBin().

int CTPPSPixelDQMSource::pixColMAX = 156
private

Definition at line 68 of file CTPPSPixelDQMSource.cc.

Referenced by dqmBeginRun().

int CTPPSPixelDQMSource::pixRowMAX = 160
private

Definition at line 67 of file CTPPSPixelDQMSource.cc.

Referenced by bookHistograms(), and dqmBeginRun().

int CTPPSPixelDQMSource::ROCSizeInX = pixRowMAX/2
private

Definition at line 69 of file CTPPSPixelDQMSource.cc.

Referenced by bookHistograms(), and dqmBeginRun().

int CTPPSPixelDQMSource::ROCSizeInY = pixColMAX/3
private

Definition at line 70 of file CTPPSPixelDQMSource.cc.

Referenced by bookHistograms(), and dqmBeginRun().

int CTPPSPixelDQMSource::RPindexValid[RPotsTotalNumber]
private

Definition at line 85 of file CTPPSPixelDQMSource.cc.

Referenced by analyze(), bookHistograms(), and dqmBeginRun().

int CTPPSPixelDQMSource::RPn_first = 3
staticprivate

Definition at line 57 of file CTPPSPixelDQMSource.cc.

Referenced by bookHistograms().

int CTPPSPixelDQMSource::RPn_last = 4
staticprivate

Definition at line 57 of file CTPPSPixelDQMSource.cc.

Referenced by bookHistograms().

int CTPPSPixelDQMSource::RPotsIDMAX =8
staticprivate

Definition at line 60 of file CTPPSPixelDQMSource.cc.

Referenced by dqmBeginRun().

int CTPPSPixelDQMSource::RPotsTotalNumber =NArms*NStationMAX*NRPotsMAX
staticprivate

Definition at line 83 of file CTPPSPixelDQMSource.cc.

Referenced by analyze().

int CTPPSPixelDQMSource::RPstatus[StationIDMAX][RPotsIDMAX]
private

Definition at line 106 of file CTPPSPixelDQMSource.cc.

Referenced by analyze(), bookHistograms(), and dqmBeginRun().

unsigned int CTPPSPixelDQMSource::rpStatusWord = 0x8008
private

Definition at line 105 of file CTPPSPixelDQMSource.cc.

Referenced by dqmBeginRun().

int CTPPSPixelDQMSource::StationIDMAX =4
staticprivate

Definition at line 59 of file CTPPSPixelDQMSource.cc.

Referenced by dqmBeginRun().

int CTPPSPixelDQMSource::StationStatus[StationIDMAX]
private

Definition at line 107 of file CTPPSPixelDQMSource.cc.

Referenced by analyze(), bookHistograms(), and dqmBeginRun().

CTPPSPixelIndices CTPPSPixelDQMSource::thePixIndices
private

Definition at line 66 of file CTPPSPixelDQMSource.cc.

Referenced by analyze(), and dqmBeginRun().

edm::EDGetTokenT< edm::DetSetVector<CTPPSPixelCluster> > CTPPSPixelDQMSource::tokenCluster
private

Definition at line 49 of file CTPPSPixelDQMSource.cc.

Referenced by analyze(), and CTPPSPixelDQMSource().

edm::EDGetTokenT<edm::DetSetVector<CTPPSPixelDigi> > CTPPSPixelDQMSource::tokenDigi
private

Definition at line 48 of file CTPPSPixelDQMSource.cc.

Referenced by analyze(), and CTPPSPixelDQMSource().

edm::EDGetTokenT< edm::DetSetVector<CTPPSPixelLocalTrack> > CTPPSPixelDQMSource::tokenTrack
private

Definition at line 50 of file CTPPSPixelDQMSource.cc.

Referenced by analyze(), and CTPPSPixelDQMSource().

int CTPPSPixelDQMSource::TrackFitDimension = 4
private

Definition at line 72 of file CTPPSPixelDQMSource.cc.

Referenced by analyze(), and dqmBeginRun().

unsigned int CTPPSPixelDQMSource::verbosity
private

Definition at line 46 of file CTPPSPixelDQMSource.cc.

Referenced by analyze(), dqmBeginRun(), and endRun().

float CTPPSPixelDQMSource::x0_MAX
private

Definition at line 142 of file CTPPSPixelDQMSource.cc.

Referenced by analyze(), dqmBeginRun(), and endRun().

float CTPPSPixelDQMSource::x0_MIN
private

Definition at line 142 of file CTPPSPixelDQMSource.cc.

Referenced by analyze(), dqmBeginRun(), and endRun().

float CTPPSPixelDQMSource::y0_MAX
private

Definition at line 142 of file CTPPSPixelDQMSource.cc.

Referenced by analyze(), dqmBeginRun(), and endRun().

float CTPPSPixelDQMSource::y0_MIN
private

Definition at line 142 of file CTPPSPixelDQMSource.cc.

Referenced by analyze(), dqmBeginRun(), and endRun().