CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
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
 

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

MonitorElementh2AllPlanesActive
 
MonitorElementh2CluSize [NArms][NStationMAX]
 
MonitorElementh2Efficiency [RPotsTotalNumber][NplaneMAX]
 
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
 
bool isPlanePlotsTurnedOff [NArms][NStationMAX][NRPotsMAX][NplaneMAX] = {}
 
const float mapXmax = 30. * TMath::Cos(18.4 / 180. * TMath::Pi())
 
const float mapXmin = 0. * TMath::Cos(18.4 / 180. * TMath::Pi())
 
long int nEvents = 0
 
bool offlinePlots = true
 
bool onlinePlots = true
 
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 constexpr int ADCMax = 256
 
static constexpr int ClusMultMAX = 10
 
static constexpr int ClusterSizeMax = 9
 
static constexpr int hitMultMAX = 50
 
static constexpr int mapXbins = 200
 
static constexpr int mapYbins = 240
 
static constexpr float mapYmax = 8.
 
static constexpr float mapYmin = -16.
 
static constexpr int NArms = 2
 
static constexpr int NLocalTracksMAX = 20
 
static constexpr int NPlaneBins = NplaneMAX * NRPotBinsInStation
 
static constexpr int NplaneMAX = 6
 
static constexpr int NROCsMAX = 6
 
static constexpr int NRPglobalBins = 4
 
static constexpr int NRPotBinsInStation = RPn_last - RPn_first
 
static constexpr int NRPotsMAX = 6
 
static constexpr int NStationMAX = 3
 
static constexpr int RPn_first = 3
 
static constexpr int RPn_last = 4
 
static constexpr int RPotsIDMAX = 8
 
static constexpr int RPotsTotalNumber = NArms * NStationMAX * NRPotsMAX
 
static constexpr int StationIDMAX = 4
 

Detailed Description

Definition at line 32 of file CTPPSPixelDQMSource.cc.

Constructor & Destructor Documentation

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

Definition at line 162 of file CTPPSPixelDQMSource.cc.

References submitPVResolutionJobs::count, Exception, edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), isPlanePlotsTurnedOff, NArms, NplaneMAX, NRPotsMAX, NStationMAX, offlinePlots, onlinePlots, alignCSCRings::s, relativeConstraints::station, tokenCluster, tokenDigi, tokenTrack, and verbosity.

163  : verbosity(ps.getUntrackedParameter<unsigned int>("verbosity", 0)),
164  rpStatusWord(ps.getUntrackedParameter<unsigned int>("RPStatusWord", 0x8008)) {
165  tokenDigi = consumes<DetSetVector<CTPPSPixelDigi>>(ps.getParameter<edm::InputTag>("tagRPixDigi"));
166  tokenCluster = consumes<DetSetVector<CTPPSPixelCluster>>(ps.getParameter<edm::InputTag>("tagRPixCluster"));
167  tokenTrack = consumes<DetSetVector<CTPPSPixelLocalTrack>>(ps.getParameter<edm::InputTag>("tagRPixLTrack"));
168  offlinePlots = ps.getUntrackedParameter<bool>("offlinePlots", true);
169  onlinePlots = ps.getUntrackedParameter<bool>("onlinePlots", true);
170 
171  vector<string> disabledPlanePlotsVec =
172  ps.getUntrackedParameter<vector<string>>("turnOffPlanePlots", vector<string>());
173 
174  // Parse the strings in disabledPlanePlotsVec and set the flags in
175  // isPlanePlotsTurnedOff
176  for (auto s : disabledPlanePlotsVec) {
177  // Check that the format is <arm>_<station>_<RP>_<Plane>
178  if (count(s.begin(), s.end(), '_') != 3)
179  throw cms::Exception("RPixPlaneCombinatoryTracking") << "Invalid string in turnOffPlanePlots: " << s;
180  else {
181  vector<string> armStationRpPlane;
182  size_t pos = 0;
183  while ((pos = s.find('_')) != string::npos) {
184  armStationRpPlane.push_back(s.substr(0, pos));
185  s.erase(0, pos + 1);
186  }
187  armStationRpPlane.push_back(s);
188 
189  int arm = stoi(armStationRpPlane.at(0));
190  int station = stoi(armStationRpPlane.at(1));
191  int rp = stoi(armStationRpPlane.at(2));
192  int plane = stoi(armStationRpPlane.at(3));
193 
194  if (arm < NArms && station < NStationMAX && rp < NRPotsMAX && plane < NplaneMAX) {
195  if (verbosity)
196  LogPrint("CTPPSPixelDQMSource")
197  << "Shutting off plots for: Arm " << arm << " Station " << station << " Rp " << rp << " Plane " << plane;
198  isPlanePlotsTurnedOff[arm][station][rp][plane] = true;
199  } else {
200  throw cms::Exception("RPixPlaneCombinatoryTracking") << "Invalid string in turnOffPlanePlots: " << s;
201  }
202  }
203  }
204 }
T getUntrackedParameter(std::string const &, T const &) const
edm::EDGetTokenT< edm::DetSetVector< CTPPSPixelCluster > > tokenCluster
edm::EDGetTokenT< edm::DetSetVector< CTPPSPixelDigi > > tokenDigi
static constexpr int NRPotsMAX
edm::EDGetTokenT< edm::DetSetVector< CTPPSPixelLocalTrack > > tokenTrack
Log< level::Warning, true > LogPrint
static constexpr int NStationMAX
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
static constexpr int NArms
static constexpr int NplaneMAX
bool isPlanePlotsTurnedOff[NArms][NStationMAX][NRPotsMAX][NplaneMAX]
CTPPSPixelDQMSource::~CTPPSPixelDQMSource ( )
override

Definition at line 208 of file CTPPSPixelDQMSource.cc.

208 {}

Member Function Documentation

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

Definition at line 474 of file CTPPSPixelDQMSource.cc.

References gpuClustering::adc, CTPPSDetId::arm(), edm::EventBase::bunchCrossing(), ClusterSizeMax, cuy::col, dqm::impl::MonitorElement::Fill(), HcalObjRepresent::Fill(), getDet(), getPixPlane(), getPlaneIndex(), getRPglobalBin(), getRPindex(), h2AllPlanesActive, h2CluSize, h2Efficiency, h2HitsMultipl, h2HitsMultROC, h2trackXY0, h2xyHits, hBX, hBXshort, hHitsMult, hitMultMAX, HitsMultPlane, HitsMultROC, hp2HitsMultROC_LS, hp2xyADC, hpixLTrack, hpRPactive, hRPotActivBX, hRPotActivBXall, hRPotActivBXroc, hRPotActivPlanes, htrackHits, htrackMult, mps_fire::i, isPlanePlotsTurnedOff, edm::HandleBase::isValid(), SiStripPI::max, nEvents, cms::cuda::nh, NLocalTracksMAX, np, NplaneMAX, NROCsMAX, NRPotsMAX, NStationMAX, offlinePlots, onlinePlots, AlCaHLTBitMon_ParallelJobs::p, prIndex(), alignCSCRings::r, PixelMapPlotter::roc, 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.

474  {
475  ++nEvents;
476  int lumiId = event.getLuminosityBlock().id().luminosityBlock();
477  if (lumiId < 0)
478  lumiId = 0;
479 
480  int RPactivity[RPotsTotalNumber], RPdigiSize[RPotsTotalNumber];
481  int pixRPTracks[RPotsTotalNumber];
482 
483  for (int rp = 0; rp < RPotsTotalNumber; rp++) {
484  RPactivity[rp] = RPdigiSize[rp] = pixRPTracks[rp] = 0;
485  }
486 
487  for (int ind = 0; ind < RPotsTotalNumber; ind++) {
488  for (int p = 0; p < NplaneMAX; p++) {
489  HitsMultPlane[ind][p] = 0;
490  }
491  }
492  for (int ind = 0; ind < RPotsTotalNumber * NplaneMAX; ind++) {
493  for (int roc = 0; roc < NROCsMAX; roc++) {
494  HitsMultROC[ind][roc] = 0;
495  }
496  }
498  event.getByToken(tokenDigi, pixDigi);
499 
501  event.getByToken(tokenCluster, pixClus);
502 
504  event.getByToken(tokenTrack, pixTrack);
505 
506  if (onlinePlots) {
507  hBX->Fill(event.bunchCrossing());
508  hBXshort->Fill(event.bunchCrossing());
509  }
510 
511  if (pixTrack.isValid()) {
512  for (const auto &ds_tr : *pixTrack) {
513  int idet = getDet(ds_tr.id);
514  if (idet != DetId::VeryForward) {
515  if (verbosity > 1)
516  LogPrint("CTPPSPixelDQMSource") << "not CTPPS: ds_tr.id" << ds_tr.id;
517  continue;
518  }
519  CTPPSDetId theId(ds_tr.id);
520  int arm = theId.arm() & 0x1;
521  int station = theId.station() & 0x3;
522  int rpot = theId.rp() & 0x7;
523  int rpInd = getRPindex(arm, station, rpot);
524 
525  for (DetSet<CTPPSPixelLocalTrack>::const_iterator dit = ds_tr.begin(); dit != ds_tr.end(); ++dit) {
526  ++pixRPTracks[rpInd];
527  int nh_tr = (dit->ndf() + TrackFitDimension) / 2;
528  if (onlinePlots) {
529  for (int i = 0; i <= NplaneMAX; i++) {
530  if (i == nh_tr)
531  htrackHits[rpInd]->Fill(nh_tr, 1.);
532  else
533  htrackHits[rpInd]->Fill(i, 0.);
534  }
535  }
536  float x0 = dit->x0();
537  float y0 = dit->y0();
538  h2trackXY0[rpInd]->Fill(x0, y0);
539 
540  if (x0_MAX < x0)
541  x0_MAX = x0;
542  if (y0_MAX < y0)
543  y0_MAX = y0;
544  if (x0_MIN > x0)
545  x0_MIN = x0;
546  if (y0_MIN > y0)
547  y0_MIN = y0;
548 
549  if (offlinePlots) {
550  edm::DetSetVector<CTPPSPixelFittedRecHit> fittedHits = dit->hits();
551 
552  std::map<int, int> numberOfPointPerPlaneEff;
553  for (const auto &ds_frh : fittedHits) {
554  int plane = getPixPlane(ds_frh.id);
555  for (DetSet<CTPPSPixelFittedRecHit>::const_iterator frh_it = ds_frh.begin(); frh_it != ds_frh.end();
556  ++frh_it) { // there should always be only one hit in each
557  // vector
558  if (frh_it != ds_frh.begin())
559  if (verbosity > 1)
560  LogPrint("CTPPSPixelDQMSource") << "More than one FittedRecHit found in plane " << plane;
561  if (frh_it->isRealHit())
562  for (int p = 0; p < NplaneMAX; p++) {
563  if (p != plane)
564  numberOfPointPerPlaneEff[p]++;
565  }
566  }
567  }
568 
569  if (verbosity > 1)
570  for (auto planeAndHitsOnOthers : numberOfPointPerPlaneEff) {
571  LogPrint("CTPPSPixelDQMSource")
572  << "For plane " << planeAndHitsOnOthers.first << ", " << planeAndHitsOnOthers.second
573  << " hits on other planes were found" << endl;
574  }
575 
576  for (const auto &ds_frh : fittedHits) {
577  int plane = getPixPlane(ds_frh.id);
578  if (isPlanePlotsTurnedOff[arm][station][rpot][plane])
579  continue;
580  for (DetSet<CTPPSPixelFittedRecHit>::const_iterator frh_it = ds_frh.begin(); frh_it != ds_frh.end();
581  ++frh_it) {
582  float frhX0 = frh_it->globalCoordinates().x() + frh_it->xResidual();
583  float frhY0 = frh_it->globalCoordinates().y() + frh_it->yResidual();
584  if (numberOfPointPerPlaneEff[plane] >= 3) {
585  if (frh_it->isRealHit())
586  h2Efficiency[rpInd][plane]->Fill(frhX0, frhY0, 1);
587  else
588  h2Efficiency[rpInd][plane]->Fill(frhX0, frhY0, 0);
589  }
590  }
591  }
592  }
593  }
594  }
595  } // end if(pixTrack.isValid())
596 
597  bool valid = false;
598  valid |= pixDigi.isValid();
599  // valid |= Clus.isValid();
600 
601  if (!valid && verbosity)
602  LogPrint("CTPPSPixelDQMSource") << "No valid data in Event " << nEvents;
603 
604  if (pixDigi.isValid()) {
605  for (const auto &ds_digi : *pixDigi) {
606  int idet = getDet(ds_digi.id);
607  if (idet != DetId::VeryForward) {
608  if (verbosity > 1)
609  LogPrint("CTPPSPixelDQMSource") << "not CTPPS: ds_digi.id" << ds_digi.id;
610  continue;
611  }
612  // int subdet = getSubdet(ds_digi.id);
613 
614  int plane = getPixPlane(ds_digi.id);
615 
616  CTPPSDetId theId(ds_digi.id);
617  int arm = theId.arm() & 0x1;
618  int station = theId.station() & 0x3;
619  int rpot = theId.rp() & 0x7;
620  int rpInd = getRPindex(arm, station, rpot);
621  RPactivity[rpInd] = 1;
622  ++RPdigiSize[rpInd];
623 
624  if (StationStatus[station] && RPstatus[station][rpot]) {
625  if (onlinePlots) {
626  h2HitsMultipl[arm][station]->Fill(prIndex(rpot, plane), ds_digi.data.size());
627  h2AllPlanesActive->Fill(plane, getRPglobalBin(arm, station));
628  }
629  int index = getRPindex(arm, station, rpot);
630  HitsMultPlane[index][plane] += ds_digi.data.size();
631  if (RPindexValid[index]) {
632  int nh = ds_digi.data.size();
633  if (nh > hitMultMAX)
634  nh = hitMultMAX;
635  if (!isPlanePlotsTurnedOff[arm][station][rpot][plane])
636  if (onlinePlots)
637  hHitsMult[index][plane]->Fill(nh);
638  }
639  int rocHistIndex = getPlaneIndex(arm, station, rpot, plane);
640 
641  for (DetSet<CTPPSPixelDigi>::const_iterator dit = ds_digi.begin(); dit != ds_digi.end(); ++dit) {
642  int row = dit->row();
643  int col = dit->column();
644  int adc = dit->adc();
645 
646  if (RPindexValid[index]) {
647  if (!isPlanePlotsTurnedOff[arm][station][rpot][plane]) {
648  if (onlinePlots)
649  h2xyHits[index][plane]->Fill(col, row);
650  hp2xyADC[index][plane]->Fill(col, row, adc);
651  }
652  int colROC, rowROC;
653  int trocId;
654  if (!thePixIndices.transformToROC(col, row, trocId, colROC, rowROC)) {
655  if (trocId >= 0 && trocId < NROCsMAX) {
656  ++HitsMultROC[rocHistIndex][trocId];
657  }
658  }
659  } // end if(RPindexValid[index]) {
660  }
661  } // end if(StationStatus[station]) {
662  } // end for(const auto &ds_digi : *pixDigi)
663  } // if(pixDigi.isValid()) {
664 
665  if (pixClus.isValid() && onlinePlots)
666  for (const auto &ds : *pixClus) {
667  int idet = getDet(ds.id);
668  if (idet != DetId::VeryForward && verbosity > 1) {
669  LogPrint("CTPPSPixelDQMSource") << "not CTPPS: cluster.id" << ds.id;
670  continue;
671  }
672 
673  CTPPSDetId theId(ds.id);
674  int plane = getPixPlane(ds.id);
675  int arm = theId.arm() & 0x1;
676  int station = theId.station() & 0x3;
677  int rpot = theId.rp() & 0x7;
678 
679  if ((StationStatus[station] == 0) || (RPstatus[station][rpot] == 0))
680  continue;
681 
682  for (const auto &p : ds) {
683  int clusize = p.size();
684 
685  if (clusize > ClusterSizeMax)
686  clusize = ClusterSizeMax;
687 
688  h2CluSize[arm][station]->Fill(prIndex(rpot, plane), clusize);
689  }
690  } // end if(pixClus.isValid()) for(const auto &ds : *pixClus)
691 
692  bool allRPactivity = false;
693  for (int rp = 0; rp < RPotsTotalNumber; rp++)
694  if (RPactivity[rp] > 0)
695  allRPactivity = true;
696  for (int arm = 0; arm < 2; arm++) {
697  for (int stn = 0; stn < NStationMAX; stn++) {
698  for (int rp = 0; rp < NRPotsMAX; rp++) {
699  int index = getRPindex(arm, stn, rp);
700  if (RPindexValid[index] == 0)
701  continue;
702 
703  if (onlinePlots)
704  hpRPactive->Fill(getRPglobalBin(arm, stn), RPactivity[index]);
705  // if(RPactivity[index]==0) continue;
706  if (!allRPactivity)
707  continue;
708  if (onlinePlots)
709  hpixLTrack->Fill(getRPglobalBin(arm, stn), pixRPTracks[index]);
710  int ntr = pixRPTracks[index];
711  if (ntr > NLocalTracksMAX)
712  ntr = NLocalTracksMAX;
713  for (int i = 0; i <= NLocalTracksMAX; i++) {
714  if (i == ntr)
715  htrackMult[index]->Fill(ntr, 1.);
716  else
717  htrackMult[index]->Fill(i, 0.);
718  }
719 
720  int np = 0;
721  for (int p = 0; p < NplaneMAX; p++)
722  if (HitsMultPlane[index][p] > 0)
723  np++;
724  for (int p = 0; p <= NplaneMAX; p++) {
725  if (p == np)
726  hRPotActivPlanes[index]->Fill(p, 1.);
727  else
728  hRPotActivPlanes[index]->Fill(p, 0.);
729  }
730  if (onlinePlots) {
731  if (np >= 5)
732  hRPotActivBX[index]->Fill(event.bunchCrossing());
733  hRPotActivBXall[index]->Fill(event.bunchCrossing(), float(RPdigiSize[index]));
734  }
735 
736  int planesFiredAtROC[NROCsMAX]; // how many planes registered hits on ROC r
737  for (int r = 0; r < NROCsMAX; r++)
738  planesFiredAtROC[r] = 0;
739  for (int p = 0; p < NplaneMAX; p++) {
740  int indp = getPlaneIndex(arm, stn, rp, p);
741  for (int r = 0; r < NROCsMAX; r++)
742  if (HitsMultROC[indp][r] > 0)
743  ++planesFiredAtROC[r];
744  for (int r = 0; r < NROCsMAX; r++) {
745  if (onlinePlots)
746  h2HitsMultROC[index]->Fill(p, r, HitsMultROC[indp][r]);
747  hp2HitsMultROC_LS[index]->Fill(lumiId, p * NROCsMAX + r, HitsMultROC[indp][r]);
748  }
749  }
750  int max = 0;
751  for (int r = 0; r < NROCsMAX; r++)
752  if (max < planesFiredAtROC[r])
753  max = planesFiredAtROC[r];
754  if (max >= 4 && onlinePlots) // fill only if there are at least 4 aligned ROCs firing
755  hRPotActivBXroc[index]->Fill(event.bunchCrossing());
756  } // end for(int rp=0; rp<NRPotsMAX; rp++) {
757  }
758  } // end for(int arm=0; arm<2; arm++) {
759 
760  if ((nEvents % 100))
761  return;
762  if (verbosity)
763  LogPrint("CTPPSPixelDQMSource") << "analyze event " << nEvents;
764 }
int StationStatus[StationIDMAX]
MonitorElement * hHitsMult[RPotsTotalNumber][NplaneMAX]
int HitsMultPlane[RPotsTotalNumber][NplaneMAX]
MonitorElement * hp2xyADC[RPotsTotalNumber][NplaneMAX]
MonitorElement * hpRPactive
MonitorElement * hp2HitsMultROC_LS[RPotsTotalNumber]
edm::EDGetTokenT< edm::DetSetVector< CTPPSPixelCluster > > tokenCluster
MonitorElement * h2CluSize[NArms][NStationMAX]
edm::EDGetTokenT< edm::DetSetVector< CTPPSPixelDigi > > tokenDigi
MonitorElement * htrackMult[RPotsTotalNumber]
MonitorElement * h2HitsMultROC[RPotsTotalNumber]
MonitorElement * h2trackXY0[RPotsTotalNumber]
MonitorElement * hpixLTrack
static constexpr int hitMultMAX
static constexpr int NRPotsMAX
CTPPSPixelIndices thePixIndices
static constexpr int RPotsTotalNumber
void Fill(long long x)
MonitorElement * htrackHits[RPotsTotalNumber]
edm::EDGetTokenT< edm::DetSetVector< CTPPSPixelLocalTrack > > tokenTrack
int prIndex(int rp, int plane)
int RPindexValid[RPotsTotalNumber]
MonitorElement * h2Efficiency[RPotsTotalNumber][NplaneMAX]
int np
Definition: AMPTWrapper.h:43
void Fill(HcalDetId &id, double val, std::vector< TH2F > &depth)
MonitorElement * hRPotActivBX[RPotsTotalNumber]
static constexpr int NLocalTracksMAX
static constexpr int ClusterSizeMax
uint32_t nh
bool isValid() const
Definition: HandleBase.h:70
int getRPindex(int arm, int station, int rp)
Log< level::Warning, true > LogPrint
MonitorElement * hRPotActivBXroc[RPotsTotalNumber]
MonitorElement * hRPotActivBXall[RPotsTotalNumber]
static constexpr int NStationMAX
static constexpr int NROCsMAX
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
MonitorElement * h2AllPlanesActive
int getRPglobalBin(int arm, int stn)
int HitsMultROC[RPotsTotalNumber *NplaneMAX][NROCsMAX]
static constexpr int NplaneMAX
collection_type::const_iterator const_iterator
Definition: DetSet.h:31
int getPlaneIndex(int arm, int station, int rp, int plane)
int col
Definition: cuy.py:1009
MonitorElement * h2xyHits[RPotsTotalNumber][NplaneMAX]
MonitorElement * hRPotActivPlanes[RPotsTotalNumber]
int RPstatus[StationIDMAX][RPotsIDMAX]
MonitorElement * hBXshort
bool isPlanePlotsTurnedOff[NArms][NStationMAX][NRPotsMAX][NplaneMAX]
MonitorElement * h2HitsMultipl[NArms][NStationMAX]
uint16_t *__restrict__ uint16_t const *__restrict__ adc
void CTPPSPixelDQMSource::bookHistograms ( DQMStore::IBooker ibooker,
edm::Run const &  ,
edm::EventSetup const &   
)
overrideprotected

Definition at line 247 of file CTPPSPixelDQMSource.cc.

References CTPPSDetId::armName(), dqm::implementation::IBooker::book1D(), dqm::implementation::IBooker::book1DD(), dqm::implementation::IBooker::book2D(), dqm::implementation::IBooker::book2DD(), dqm::implementation::IBooker::bookProfile(), dqm::implementation::IBooker::bookProfile2D(), dqm::implementation::NavigatorBase::cd(), ClusterSizeMax, getRPglobalBin(), getRPindex(), dqm::legacy::MonitorElement::getTH2D(), dqm::legacy::MonitorElement::getTH2F(), dqm::legacy::MonitorElement::getTProfile(), dqm::legacy::MonitorElement::getTProfile2D(), h2AllPlanesActive, h2CluSize, h2Efficiency, h2HitsMultipl, h2HitsMultROC, h2trackXY0, h2xyHits, hBX, hBXshort, hHitsMult, hitMultMAX, hp2HitsMultROC_LS, hp2xyADC, hpixLTrack, hpRPactive, hRPotActivBX, hRPotActivBXall, hRPotActivBXroc, hRPotActivPlanes, htrackHits, htrackMult, isPlanePlotsTurnedOff, mapXbins, mapXmax, mapXmin, mapYbins, mapYmax, mapYmin, hlt_dqm_clientPB-live_cfg::nbins, CTPPSDetId::nFull, NLocalTracksMAX, CTPPSDetId::nPath, NPlaneBins, NplaneMAX, NROCsMAX, NRPglobalBins, CTPPSDetId::nShort, NStationMAX, offlinePlots, onlinePlots, AlCaHLTBitMon_ParallelJobs::p, alignCSCRings::r, CTPPSDetId::rpId(), RPindexValid, RPn_first, RPn_last, RPstatus, alignCSCRings::s, sd, CTPPSDetId::sdTrackingPixel, dqm::implementation::NavigatorBase::setCurrentFolder(), CTPPSDetId::setRP(), CTPPSDetId::setStation(), CTPPSDetId::stationId(), StationStatus, and AlCaHLTBitMon_QueryRunRegistry::string.

247  {
248  ibooker.cd();
249  ibooker.setCurrentFolder("CTPPS/TrackingPixel");
250  char s[50];
251  string armTitleShort, stnTitleShort;
252 
253  TAxis *yah1st = nullptr;
254  TAxis *xaRPact = nullptr;
255  TAxis *xah1trk = nullptr;
256  if (onlinePlots) {
257  hBX = ibooker.book1D("events per BX", "ctpps_pixel;Event.BX", 4002, -1.5, 4000. + 0.5);
258  hBXshort = ibooker.book1D("events per BX(short)", "ctpps_pixel;Event.BX", 102, -1.5, 100. + 0.5);
259 
260  string str1st = "Pixel planes activity";
261  h2AllPlanesActive = ibooker.book2DD(
262  str1st, str1st + "(digi task);Plane #", NplaneMAX, 0, NplaneMAX, NRPglobalBins, 0.5, NRPglobalBins + 0.5);
263  TH2D *h1st = h2AllPlanesActive->getTH2D();
264  h1st->SetOption("colz");
265  yah1st = h1st->GetYaxis();
266 
267  string str2 = "Pixel RP active";
268  hpRPactive = ibooker.bookProfile(
269  str2, str2 + " per event(digi task)", NRPglobalBins, 0.5, NRPglobalBins + 0.5, -0.1, 1.1, "");
270  xaRPact = hpRPactive->getTProfile()->GetXaxis();
271  hpRPactive->getTProfile()->SetOption("hist");
272  hpRPactive->getTProfile()->SetMinimum(0.);
273  hpRPactive->getTProfile()->SetMaximum(1.1);
274 
275  str2 = "Pixel Local Tracks";
276  hpixLTrack = ibooker.bookProfile(
277  str2, str2 + " per event", NRPglobalBins, 0.5, NRPglobalBins + 0.5, -0.1, NLocalTracksMAX, "");
278 
279  xah1trk = hpixLTrack->getTProfile()->GetXaxis();
280  hpixLTrack->getTProfile()->GetYaxis()->SetTitle("average number of tracks per event");
281  hpixLTrack->getTProfile()->SetOption("hist");
282  }
283 
284  for (int arm = 0; arm < 2; arm++) {
286  string sd, armTitle;
287  ID.armName(sd, CTPPSDetId::nPath);
288  ID.armName(armTitle, CTPPSDetId::nFull);
289  ID.armName(armTitleShort, CTPPSDetId::nShort);
290 
291  ibooker.setCurrentFolder(sd);
292 
293  for (int stn = 0; stn < NStationMAX; stn++) {
294  if (StationStatus[stn] == 0)
295  continue;
296  ID.setStation(stn);
297  string stnd, stnTitle;
298 
299  CTPPSDetId(ID.stationId()).stationName(stnd, CTPPSDetId::nPath);
300  CTPPSDetId(ID.stationId()).stationName(stnTitle, CTPPSDetId::nFull);
301  CTPPSDetId(ID.stationId()).stationName(stnTitleShort, CTPPSDetId::nShort);
302 
303  ibooker.setCurrentFolder(stnd);
304  //--------- RPots ---
305  int pixBinW = 4;
306  for (int rp = RPn_first; rp < RPn_last; rp++) { // only installed pixel pots
307  ID.setRP(rp);
308  string rpd, rpTitle;
309  CTPPSDetId(ID.rpId()).rpName(rpTitle, CTPPSDetId::nShort);
310  string rpBinName = armTitleShort + "_" + stnTitleShort + "_" + rpTitle;
311  if (onlinePlots) {
312  yah1st->SetBinLabel(getRPglobalBin(arm, stn), rpBinName.c_str());
313  xah1trk->SetBinLabel(getRPglobalBin(arm, stn), rpBinName.c_str());
314  xaRPact->SetBinLabel(getRPglobalBin(arm, stn), rpBinName.c_str());
315  }
316  if (RPstatus[stn][rp] == 0)
317  continue;
318  int indexP = getRPindex(arm, stn, rp);
319  RPindexValid[indexP] = 1;
320 
321  CTPPSDetId(ID.rpId()).rpName(rpTitle, CTPPSDetId::nFull);
322  CTPPSDetId(ID.rpId()).rpName(rpd, CTPPSDetId::nPath);
323 
324  ibooker.setCurrentFolder(rpd);
325 
326  const float x0Maximum = 70.;
327  const float y0Maximum = 15.;
328  string st = "track intercept point";
329  string st2 = ": " + stnTitle;
330  h2trackXY0[indexP] = ibooker.book2D(
331  st, st + st2 + ";x0;y0", int(x0Maximum) * 2, 0., x0Maximum, int(y0Maximum) * 4, -y0Maximum, y0Maximum);
332  h2trackXY0[indexP]->getTH2F()->SetOption("colz");
333 
334  st = "number of tracks per event";
335  htrackMult[indexP] = ibooker.bookProfile(st,
336  rpTitle + ";number of tracks",
337  NLocalTracksMAX + 1,
338  -0.5,
339  NLocalTracksMAX + 0.5,
340  -0.5,
341  NLocalTracksMAX + 0.5,
342  "");
343  htrackMult[indexP]->getTProfile()->SetOption("hist");
344 
345  hRPotActivPlanes[indexP] = ibooker.bookProfile("number of fired planes per event",
346  rpTitle + ";nPlanes;Probability",
347  NplaneMAX + 1,
348  -0.5,
349  NplaneMAX + 0.5,
350  -0.5,
351  NplaneMAX + 0.5,
352  "");
353  hRPotActivPlanes[indexP]->getTProfile()->SetOption("hist");
354 
355  hp2HitsMultROC_LS[indexP] = ibooker.bookProfile2D("ROCs hits multiplicity per event vs LS",
356  rpTitle + ";LumiSection;Plane#___ROC#",
357  1000,
358  0.,
359  1000.,
361  0.,
362  double(NplaneMAX * NROCsMAX),
363  0.,
364  ROCSizeInX *ROCSizeInY,
365  "");
366  hp2HitsMultROC_LS[indexP]->getTProfile2D()->SetOption("colz");
367  hp2HitsMultROC_LS[indexP]->getTProfile2D()->SetMinimum(1.0e-10);
368  hp2HitsMultROC_LS[indexP]->getTProfile2D()->SetCanExtend(TProfile2D::kXaxis);
369  TAxis *yahp2 = hp2HitsMultROC_LS[indexP]->getTProfile2D()->GetYaxis();
370  for (int p = 0; p < NplaneMAX; p++) {
371  sprintf(s, "plane%d_0", p);
372  yahp2->SetBinLabel(p * NplaneMAX + 1, s);
373  for (int r = 1; r < NROCsMAX; r++) {
374  sprintf(s, " %d_%d", p, r);
375  yahp2->SetBinLabel(p * NplaneMAX + r + 1, s);
376  }
377  }
378 
379  if (onlinePlots) {
380  string st3 = ";PlaneIndex(=pixelPot*PlaneMAX + plane)";
381 
382  st = "hit multiplicity in planes";
383  h2HitsMultipl[arm][stn] = ibooker.book2DD(
384  st, st + st2 + st3 + ";multiplicity", NPlaneBins, 0, NPlaneBins, hitMultMAX, 0, hitMultMAX);
385  h2HitsMultipl[arm][stn]->getTH2D()->SetOption("colz");
386 
387  st = "cluster size in planes";
388  h2CluSize[arm][stn] = ibooker.book2D(st,
389  st + st2 + st3 + ";Cluster size",
390  NPlaneBins,
391  0,
392  NPlaneBins,
393  ClusterSizeMax + 1,
394  0,
395  ClusterSizeMax + 1);
396  h2CluSize[arm][stn]->getTH2F()->SetOption("colz");
397 
398  st = "number of hits per track";
399  htrackHits[indexP] = ibooker.bookProfile(st, rpTitle + ";number of hits", 5, 1.5, 6.5, -0.1, 1.1, "");
400  htrackHits[indexP]->getTProfile()->SetOption("hist");
401 
402  h2HitsMultROC[indexP] = ibooker.bookProfile2D("ROCs hits multiplicity per event",
403  rpTitle + ";plane # ;ROC #",
404  NplaneMAX,
405  -0.5,
406  NplaneMAX - 0.5,
407  NROCsMAX,
408  -0.5,
409  NROCsMAX - 0.5,
410  0.,
411  ROCSizeInX * ROCSizeInY,
412  "");
413  h2HitsMultROC[indexP]->getTProfile2D()->SetOption("colztext");
414  h2HitsMultROC[indexP]->getTProfile2D()->SetMinimum(1.e-10);
415 
416  ibooker.setCurrentFolder(rpd + "/latency");
417  hRPotActivBX[indexP] =
418  ibooker.book1D("5 fired planes per BX", rpTitle + ";Event.BX", 4002, -1.5, 4000. + 0.5);
419 
420  hRPotActivBXroc[indexP] =
421  ibooker.book1D("4 fired ROCs per BX", rpTitle + ";Event.BX", 4002, -1.5, 4000. + 0.5);
422 
423  hRPotActivBXall[indexP] = ibooker.book1D("hits per BX", rpTitle + ";Event.BX", 4002, -1.5, 4000. + 0.5);
424  }
425  int nbins = defaultDetSizeInX / pixBinW;
426 
427  for (int p = 0; p < NplaneMAX; p++) {
428  if (isPlanePlotsTurnedOff[arm][stn][rp][p])
429  continue;
430  sprintf(s, "plane_%d", p);
431  string pd = rpd + "/" + string(s);
432  ibooker.setCurrentFolder(pd);
433  string st1 = ": " + rpTitle + "_" + string(s);
434 
435  st = "adc average value";
436  hp2xyADC[indexP][p] = ibooker.bookProfile2D(
437  st, st1 + ";pix col;pix row", nbins, 0, defaultDetSizeInX, nbins, 0, defaultDetSizeInX, 0., 512., "");
438  hp2xyADC[indexP][p]->getTProfile2D()->SetOption("colz");
439 
440  if (onlinePlots) {
441  st = "hits position";
442  h2xyHits[indexP][p] = ibooker.book2DD(st,
443  st1 + ";pix col;pix row",
444  defaultDetSizeInX,
445  0,
446  defaultDetSizeInX,
447  defaultDetSizeInX,
448  0,
449  defaultDetSizeInX);
450  h2xyHits[indexP][p]->getTH2D()->SetOption("colz");
451 
452  st = "hits multiplicity";
453  hHitsMult[indexP][p] =
454  ibooker.book1DD(st, st1 + ";number of hits;N / 1 hit", hitMultMAX + 1, -0.5, hitMultMAX + 0.5);
455  }
456 
457  if (offlinePlots) {
458  st = "plane efficiency";
459  h2Efficiency[indexP][p] = ibooker.bookProfile2D(
460  st, st1 + ";x0;y0", mapXbins, mapXmin, mapXmax, mapYbins, mapYmin, mapYmax, 0, 1, "");
461  h2Efficiency[indexP][p]->getTProfile2D()->SetOption("colz");
462  }
463  } // end of for(int p=0; p<NplaneMAX;..
464 
465  } // end for(int rp=0; rp<NRPotsMAX;...
466  } // end of for(int stn=0; stn<
467  } // end of for(int arm=0; arm<2;...
468 
469  return;
470 }
int StationStatus[StationIDMAX]
virtual TH2D * getTH2D() const
MonitorElement * hHitsMult[RPotsTotalNumber][NplaneMAX]
MonitorElement * hp2xyADC[RPotsTotalNumber][NplaneMAX]
MonitorElement * hpRPactive
MonitorElement * hp2HitsMultROC_LS[RPotsTotalNumber]
static constexpr int mapYbins
virtual TH2F * getTH2F() const
MonitorElement * bookProfile2D(TString const &name, TString const &title, int nchX, double lowX, double highX, int nchY, double lowY, double highY, double lowZ, double highZ, char const *option="s", FUNC onbooking=NOOP())
Definition: DQMStore.h:399
MonitorElement * h2CluSize[NArms][NStationMAX]
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:32
MonitorElement * htrackMult[RPotsTotalNumber]
MonitorElement * h2HitsMultROC[RPotsTotalNumber]
uint32_t ID
Definition: Definitions.h:24
MonitorElement * h2trackXY0[RPotsTotalNumber]
MonitorElement * hpixLTrack
static constexpr int RPn_last
static constexpr int hitMultMAX
static constexpr int NPlaneBins
static constexpr int mapXbins
static constexpr float mapYmax
MonitorElement * htrackHits[RPotsTotalNumber]
MonitorElement * book1DD(TString const &name, TString const &title, int nchX, double lowX, double highX, FUNC onbooking=NOOP())
Definition: DQMStore.h:155
int RPindexValid[RPotsTotalNumber]
MonitorElement * h2Efficiency[RPotsTotalNumber][NplaneMAX]
MonitorElement * bookProfile(TString const &name, TString const &title, int nchX, double lowX, double highX, int, double lowY, double highY, char const *option="s", FUNC onbooking=NOOP())
Definition: DQMStore.h:322
MonitorElement * hRPotActivBX[RPotsTotalNumber]
static constexpr int NLocalTracksMAX
static constexpr int ClusterSizeMax
int getRPindex(int arm, int station, int rp)
virtual TProfile2D * getTProfile2D() const
MonitorElement * hRPotActivBXroc[RPotsTotalNumber]
MonitorElement * hRPotActivBXall[RPotsTotalNumber]
static constexpr int NRPglobalBins
double sd
virtual TProfile * getTProfile() const
static constexpr int NStationMAX
static constexpr int NROCsMAX
MonitorElement * book2D(TString const &name, TString const &title, int nchX, double lowX, double highX, int nchY, double lowY, double highY, FUNC onbooking=NOOP())
Definition: DQMStore.h:177
MonitorElement * book2DD(TString const &name, TString const &title, int nchX, double lowX, double highX, int nchY, double lowY, double highY, FUNC onbooking=NOOP())
Definition: DQMStore.h:261
Base class for CTPPS detector IDs.
Definition: CTPPSDetId.h:32
MonitorElement * h2AllPlanesActive
static constexpr int RPn_first
int getRPglobalBin(int arm, int stn)
static constexpr int NplaneMAX
MonitorElement * book1D(TString const &name, TString const &title, int const nchX, double const lowX, double const highX, FUNC onbooking=NOOP())
Definition: DQMStore.h:98
MonitorElement * h2xyHits[RPotsTotalNumber][NplaneMAX]
MonitorElement * hRPotActivPlanes[RPotsTotalNumber]
int RPstatus[StationIDMAX][RPotsIDMAX]
MonitorElement * hBXshort
bool isPlanePlotsTurnedOff[NArms][NStationMAX][NRPotsMAX][NplaneMAX]
MonitorElement * h2HitsMultipl[NArms][NStationMAX]
static constexpr float mapYmin
void CTPPSPixelDQMSource::dqmBeginRun ( edm::Run const &  run,
edm::EventSetup const &   
)
overrideprotected

Definition at line 212 of file CTPPSPixelDQMSource.cc.

References CTPPSPixelLocalTrack::dimension, nEvents, NRPotsMAX, NStationMAX, RPindexValid, RPotsIDMAX, RPstatus, rpStatusWord, StationIDMAX, StationStatus, TrackFitDimension, verbosity, x0_MAX, x0_MIN, y0_MAX, and y0_MIN.

212  {
213  if (verbosity)
214  LogPrint("CTPPSPixelDQMSource") << "RPstatusWord= " << rpStatusWord;
215  nEvents = 0;
216 
217  CTPPSPixelLocalTrack thePixelLocalTrack;
218  TrackFitDimension = thePixelLocalTrack.dimension;
219 
220  for (int stn = 0; stn < StationIDMAX; stn++) {
221  StationStatus[stn] = 0;
222  for (int rp = 0; rp < RPotsIDMAX; rp++)
223  RPstatus[stn][rp] = 0;
224  }
225 
226  unsigned int rpSts = rpStatusWord << 1;
227  for (int stn = 0; stn < NStationMAX; stn++) {
228  int stns = 0;
229  for (int rp = 0; rp < NRPotsMAX; rp++) {
230  rpSts = (rpSts >> 1);
231  RPstatus[stn][rp] = rpSts & 1;
232  if (RPstatus[stn][rp] > 0)
233  stns = 1;
234  }
235  StationStatus[stn] = stns;
236  }
237 
238  for (int ind = 0; ind < 2 * 3 * NRPotsMAX; ind++)
239  RPindexValid[ind] = 0;
240 
241  x0_MIN = y0_MIN = 1.0e06;
242  x0_MAX = y0_MAX = -1.0e06;
243 }
int StationStatus[StationIDMAX]
static constexpr int RPotsIDMAX
static constexpr int dimension
static constexpr int NRPotsMAX
int RPindexValid[RPotsTotalNumber]
Log< level::Warning, true > LogPrint
static constexpr int StationIDMAX
static constexpr int NStationMAX
int RPstatus[StationIDMAX][RPotsIDMAX]
int CTPPSPixelDQMSource::getDet ( int  id)
inlineprivate

Definition at line 148 of file CTPPSPixelDQMSource.cc.

References DetId::kDetOffset.

Referenced by analyze().

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

Definition at line 149 of file CTPPSPixelDQMSource.cc.

Referenced by analyze().

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

Definition at line 125 of file CTPPSPixelDQMSource.cc.

References getRPindex(), IndexNotValid, and NplaneMAX.

Referenced by analyze().

125  {
126  if (plane < 0 || plane >= NplaneMAX)
127  return (IndexNotValid);
128  int rc = getRPindex(arm, station, rp);
129  if (rc == IndexNotValid)
130  return (IndexNotValid);
131  return (rc * NplaneMAX + plane);
132  }
int getRPindex(int arm, int station, int rp)
static constexpr int NplaneMAX
int CTPPSPixelDQMSource::getRPglobalBin ( int  arm,
int  stn 
)
inlineprivate

Definition at line 138 of file CTPPSPixelDQMSource.cc.

References NStationMAX.

Referenced by analyze(), and bookHistograms().

138  {
139  static constexpr int stationBinOrder[NStationMAX] = {0, 4, 1};
140  return (arm * 2 + stationBinOrder[stn] + 1);
141  }
static constexpr int NStationMAX
int CTPPSPixelDQMSource::getRPindex ( int  arm,
int  station,
int  rp 
)
inlineprivate

Definition at line 116 of file CTPPSPixelDQMSource.cc.

References IndexNotValid, NRPotsMAX, NStationMAX, and relativeConstraints::station.

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

116  {
117  if (arm < 0 || station < 0 || rp < 0)
118  return (IndexNotValid);
119  if (arm > 1 || station >= NStationMAX || rp >= NRPotsMAX)
120  return (IndexNotValid);
121  int rc = (arm * NStationMAX + station) * NRPotsMAX + rp;
122  return (rc);
123  }
static constexpr int NRPotsMAX
static constexpr int NStationMAX
int CTPPSPixelDQMSource::getRPInStationBin ( int  rp)
inlineprivate

Definition at line 134 of file CTPPSPixelDQMSource.cc.

References RPn_first.

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

Definition at line 143 of file CTPPSPixelDQMSource.cc.

References NplaneMAX, and RPn_first.

Referenced by analyze().

145  {
146  return ((rp - RPn_first) * NplaneMAX + plane);
147  }
static constexpr int RPn_first
static constexpr int NplaneMAX

Member Data Documentation

constexpr int CTPPSPixelDQMSource::ADCMax = 256
staticprivate

Definition at line 55 of file CTPPSPixelDQMSource.cc.

constexpr int CTPPSPixelDQMSource::ClusMultMAX = 10
staticprivate

Definition at line 60 of file CTPPSPixelDQMSource.cc.

constexpr int CTPPSPixelDQMSource::ClusterSizeMax = 9
staticprivate

Definition at line 61 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::h2Efficiency[RPotsTotalNumber][NplaneMAX]
private

Definition at line 97 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 98 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().

constexpr int CTPPSPixelDQMSource::hitMultMAX = 50
staticprivate

Definition at line 59 of file CTPPSPixelDQMSource.cc.

Referenced by analyze(), and bookHistograms().

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

Definition at line 102 of file CTPPSPixelDQMSource.cc.

Referenced by analyze().

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

Definition at line 101 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 99 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 100 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 114 of file CTPPSPixelDQMSource.cc.

Referenced by getPlaneIndex(), and getRPindex().

bool CTPPSPixelDQMSource::isPlanePlotsTurnedOff[NArms][NStationMAX][NRPotsMAX][NplaneMAX] = {}
private

Definition at line 109 of file CTPPSPixelDQMSource.cc.

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

constexpr int CTPPSPixelDQMSource::mapXbins = 200
staticprivate

Definition at line 63 of file CTPPSPixelDQMSource.cc.

Referenced by bookHistograms().

const float CTPPSPixelDQMSource::mapXmax = 30. * TMath::Cos(18.4 / 180. * TMath::Pi())
private

Definition at line 68 of file CTPPSPixelDQMSource.cc.

Referenced by bookHistograms().

const float CTPPSPixelDQMSource::mapXmin = 0. * TMath::Cos(18.4 / 180. * TMath::Pi())
private

Definition at line 67 of file CTPPSPixelDQMSource.cc.

Referenced by bookHistograms().

constexpr int CTPPSPixelDQMSource::mapYbins = 240
staticprivate

Definition at line 64 of file CTPPSPixelDQMSource.cc.

Referenced by bookHistograms().

constexpr float CTPPSPixelDQMSource::mapYmax = 8.
staticprivate

Definition at line 66 of file CTPPSPixelDQMSource.cc.

Referenced by bookHistograms().

constexpr float CTPPSPixelDQMSource::mapYmin = -16.
staticprivate

Definition at line 65 of file CTPPSPixelDQMSource.cc.

Referenced by bookHistograms().

constexpr int CTPPSPixelDQMSource::NArms = 2
staticprivate

Definition at line 49 of file CTPPSPixelDQMSource.cc.

Referenced by CTPPSPixelDQMSource().

long int CTPPSPixelDQMSource::nEvents = 0
private

Definition at line 44 of file CTPPSPixelDQMSource.cc.

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

constexpr int CTPPSPixelDQMSource::NLocalTracksMAX = 20
staticprivate

Definition at line 58 of file CTPPSPixelDQMSource.cc.

Referenced by analyze(), and bookHistograms().

constexpr int CTPPSPixelDQMSource::NPlaneBins = NplaneMAX * NRPotBinsInStation
staticprivate

Definition at line 75 of file CTPPSPixelDQMSource.cc.

Referenced by bookHistograms().

constexpr int CTPPSPixelDQMSource::NplaneMAX = 6
staticprivate
constexpr int CTPPSPixelDQMSource::NROCsMAX = 6
staticprivate

Definition at line 53 of file CTPPSPixelDQMSource.cc.

Referenced by analyze(), and bookHistograms().

constexpr int CTPPSPixelDQMSource::NRPglobalBins = 4
staticprivate

Definition at line 136 of file CTPPSPixelDQMSource.cc.

Referenced by bookHistograms().

constexpr int CTPPSPixelDQMSource::NRPotBinsInStation = RPn_last - RPn_first
staticprivate

Definition at line 74 of file CTPPSPixelDQMSource.cc.

constexpr int CTPPSPixelDQMSource::NRPotsMAX = 6
staticprivate

Definition at line 51 of file CTPPSPixelDQMSource.cc.

Referenced by analyze(), CTPPSPixelDQMSource(), dqmBeginRun(), and getRPindex().

constexpr int CTPPSPixelDQMSource::NStationMAX = 3
staticprivate
bool CTPPSPixelDQMSource::offlinePlots = true
private

Definition at line 105 of file CTPPSPixelDQMSource.cc.

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

bool CTPPSPixelDQMSource::onlinePlots = true
private

Definition at line 106 of file CTPPSPixelDQMSource.cc.

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

int CTPPSPixelDQMSource::RPindexValid[RPotsTotalNumber]
private

Definition at line 85 of file CTPPSPixelDQMSource.cc.

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

constexpr int CTPPSPixelDQMSource::RPn_first = 3
staticprivate

Definition at line 54 of file CTPPSPixelDQMSource.cc.

Referenced by bookHistograms(), getRPInStationBin(), and prIndex().

constexpr int CTPPSPixelDQMSource::RPn_last = 4
staticprivate

Definition at line 54 of file CTPPSPixelDQMSource.cc.

Referenced by bookHistograms().

constexpr int CTPPSPixelDQMSource::RPotsIDMAX = 8
staticprivate

Definition at line 57 of file CTPPSPixelDQMSource.cc.

Referenced by dqmBeginRun().

constexpr 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 112 of file CTPPSPixelDQMSource.cc.

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

unsigned int CTPPSPixelDQMSource::rpStatusWord = 0x8008
private

Definition at line 111 of file CTPPSPixelDQMSource.cc.

Referenced by dqmBeginRun().

constexpr int CTPPSPixelDQMSource::StationIDMAX = 4
staticprivate

Definition at line 56 of file CTPPSPixelDQMSource.cc.

Referenced by dqmBeginRun().

int CTPPSPixelDQMSource::StationStatus[StationIDMAX]
private

Definition at line 113 of file CTPPSPixelDQMSource.cc.

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

CTPPSPixelIndices CTPPSPixelDQMSource::thePixIndices
private

Definition at line 70 of file CTPPSPixelDQMSource.cc.

Referenced by analyze().

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

Definition at line 46 of file CTPPSPixelDQMSource.cc.

Referenced by analyze(), and CTPPSPixelDQMSource().

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

Definition at line 45 of file CTPPSPixelDQMSource.cc.

Referenced by analyze(), and CTPPSPixelDQMSource().

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

Definition at line 47 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 43 of file CTPPSPixelDQMSource.cc.

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

float CTPPSPixelDQMSource::x0_MAX
private

Definition at line 152 of file CTPPSPixelDQMSource.cc.

Referenced by analyze(), and dqmBeginRun().

float CTPPSPixelDQMSource::x0_MIN
private

Definition at line 152 of file CTPPSPixelDQMSource.cc.

Referenced by analyze(), and dqmBeginRun().

float CTPPSPixelDQMSource::y0_MAX
private

Definition at line 152 of file CTPPSPixelDQMSource.cc.

Referenced by analyze(), and dqmBeginRun().

float CTPPSPixelDQMSource::y0_MIN
private

Definition at line 152 of file CTPPSPixelDQMSource.cc.

Referenced by analyze(), and dqmBeginRun().