CMS 3D CMS Logo

CTPPSRandomDQMSource.cc
Go to the documentation of this file.
1 /******************************************
2  *
3  * This is a part of CTPPSDQM software.
4  * Authors:
5  * A. Bellora (Universita' e INFN Torino)
6  *
7  *******************************************/
8 
14 
17 
19 
23 
24 #include <string>
25 
26 //-----------------------------------------------------------------------------
27 
29 public:
31  ~CTPPSRandomDQMSource() override = default;
32  static void fillDescriptions(edm::ConfigurationDescriptions &descriptions);
33 
34 protected:
35  void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override;
36  void analyze(edm::Event const &e, edm::EventSetup const &eSetup) override;
37 
38 private:
40 
41  static constexpr int kNArms_ = 2;
42  static constexpr int kNStationMAX_ = 3; // in an arm
43  static constexpr int kNRPotsMAX_ = 6; // per station
44  static constexpr int kNplaneMAX_ = 6; // per RPot
45  static constexpr int kFirstRPn_ = 3, kLastRPn_ = 4;
46  static constexpr int kStationIDMAX_ = 4; // possible range of ID
47  static constexpr int kRPotsIDMAX_ = 8; // possible range of ID
48 
49  const std::string folderName_ = "PPSRANDOM/RandomPixel";
50 
51  unsigned int rpStatusWord_ = 0x8008; // 220_fr_hr(stn2rp3)+ 210_fr_hr
52  int rpStatus_[kStationIDMAX_][kRPotsIDMAX_]; // symmetric in both arms
53  int stationStatus_[kStationIDMAX_]; // symmetric in both arms
54  const int kIndexNotValid = 0;
55 
57 
59 
62 
63  int getRPindex(int arm, int station, int rp) const {
64  if (arm < 0 || station < 0 || rp < 0)
65  return (kIndexNotValid);
66  if (arm > 1 || station >= kNStationMAX_ || rp >= kNRPotsMAX_)
67  return (kIndexNotValid);
68  int rc = (arm * kNStationMAX_ + station) * kNRPotsMAX_ + rp;
69  return (rc);
70  }
71 };
72 
73 //-------------------------------------------------------------------------------
74 
76  : tokenDigi_(consumes<edm::DetSetVector<CTPPSPixelDigi>>(ps.getParameter<edm::InputTag>("tagRPixDigi"))),
77  folderName_(ps.getUntrackedParameter<std::string>("folderName", "PPSRANDOM/RandomPixel")),
78  rpStatusWord_(ps.getUntrackedParameter<unsigned int>("RPStatusWord", 0x8008)) {
79  for (int stn = 0; stn < kStationIDMAX_; stn++) {
80  stationStatus_[stn] = 0;
81  for (int rp = 0; rp < kRPotsIDMAX_; rp++)
82  rpStatus_[stn][rp] = 0;
83  }
84 
85  unsigned int rpSts = rpStatusWord_ << 1;
86  for (int stn = 0; stn < kNStationMAX_; stn++) {
87  int stns = 0;
88  for (int rp = 0; rp < kNRPotsMAX_; rp++) {
89  rpSts = (rpSts >> 1);
90  rpStatus_[stn][rp] = rpSts & 1;
91  if (rpStatus_[stn][rp] > 0)
92  stns = 1;
93  }
94  stationStatus_[stn] = stns;
95  }
96 
97  for (int index = 0; index < 2 * 3 * kNRPotsMAX_; index++)
98  RPindexValid_[index] = 0;
99 }
100 
101 //--------------------------------------------------------------------------
102 
104  ibooker.cd();
105  ibooker.setCurrentFolder(folderName_);
106 
107  hBX_ = ibooker.book1D("events per BX", "ctpps_pixel;Event.BX", 4002, -1.5, 4000. + 0.5);
108 
109  for (int arm = 0; arm < kNArms_; arm++) {
111  std::string sd;
112  ID.armName(sd, CTPPSDetId::nShort);
113  sd = folderName_ + "/sector " + sd;
114 
115  ibooker.setCurrentFolder(sd);
116 
117  for (int stn = 0; stn < kNStationMAX_; stn++) {
118  if (stationStatus_[stn] == 0)
119  continue;
120  ID.setStation(stn);
121  std::string stnd;
122  CTPPSDetId(ID.stationId()).stationName(stnd, CTPPSDetId::nShort);
123  stnd = sd + "/station " + stnd;
124 
125  ibooker.setCurrentFolder(stnd);
126 
127  for (int rp = kFirstRPn_; rp < kLastRPn_; rp++) { // only installed pixel pots
128  ID.setRP(rp);
129  std::string rpd, rpTitle;
130  CTPPSDetId(ID.rpId()).rpName(rpTitle, CTPPSDetId::nFull);
131  CTPPSDetId(ID.rpId()).rpName(rpd, CTPPSDetId::nShort);
132  rpd = stnd + "/" + rpd;
133 
134  ibooker.setCurrentFolder(rpd);
135 
136  int indexP = getRPindex(arm, stn, rp);
137  RPindexValid_[indexP] = 1;
138 
139  h2HitsVsBXRandoms_[indexP] = ibooker.book2D("Digi per plane per BX - random triggers",
140  rpTitle + ";Event.BX;Plane",
141  4002,
142  -1.5,
143  4000. + 0.5,
144  kNplaneMAX_,
145  0,
146  kNplaneMAX_);
147  h2HitsVsBXRandoms_[indexP]->getTH2F()->SetOption("colz");
148 
149  } // end for(int rp=0; rp<kNRPotsMAX_;...
150  } // end of for(int stn=0; stn<
151  } // end of for(int arm=0; arm<2;...
152 
153  return;
154 }
155 
156 //-------------------------------------------------------------------------------
157 
159  auto const pixDigi = event.getHandle(tokenDigi_);
160 
161  if (!pixDigi.isValid())
162  return;
163 
164  hBX_->Fill(event.bunchCrossing());
165 
166  for (int arm = 0; arm < 2; arm++) {
167  for (int stn = 0; stn < kNStationMAX_; stn++) {
168  if (!stationStatus_[stn])
169  continue;
170  for (int rp = 0; rp < kNRPotsMAX_; rp++) {
171  if (!rpStatus_[stn][rp])
172  continue;
173  int index = getRPindex(arm, stn, rp);
174  if (RPindexValid_[index] == 0)
175  continue;
176 
177  for (int p = 0; p < kNplaneMAX_; p++) {
178  CTPPSPixelDetId planeId(arm, stn, rp, p);
179  auto pix_d = pixDigi->find(planeId.rawId());
180  if (pix_d != pixDigi->end()) {
181  int n_digis = pix_d->size();
182  h2HitsVsBXRandoms_[index]->Fill(event.bunchCrossing(), p, n_digis);
183  }
184  }
185  } // end for (int rp=0; rp<kNRPotsMAX_; rp++) {
186  } // end for (int stn = 0; stn < kNStationMAX_; stn++) {
187  } // end for (int arm=0; arm<2; arm++) {
188 }
189 
190 //---------------------------------------------------------------------------
191 
194  desc.add<edm::InputTag>("tagRPixDigi", edm::InputTag("ctppsPixelDigisAlCaRecoProducer", ""));
195  desc.addUntracked<std::string>("folderName", "PPSRANDOM/RandomPixel");
196  desc.addUntracked<unsigned int>("RPStatusWord", 0x8008);
197  descriptions.add("ctppsRandomDQMSource", desc);
198 }
199 
200 //---------------------------------------------------------------------------
int rpStatus_[kStationIDMAX_][kRPotsIDMAX_]
static constexpr int kNStationMAX_
static constexpr int kFirstRPn_
static constexpr int kLastRPn_
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
int RPindexValid_[kRPotsTotalNumber_]
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:36
static constexpr int kNplaneMAX_
uint32_t ID
Definition: Definitions.h:24
static constexpr int kRPotsIDMAX_
CTPPSRandomDQMSource(const edm::ParameterSet &ps)
edm::EDGetTokenT< edm::DetSetVector< CTPPSPixelDigi > > const tokenDigi_
static constexpr int kStationIDMAX_
void Fill(long long x)
static constexpr int kNRPotsMAX_
int getRPindex(int arm, int station, int rp) const
static constexpr int kNArms_
void analyze(edm::Event const &e, edm::EventSetup const &eSetup) override
~CTPPSRandomDQMSource() override=default
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
MonitorElement * h2HitsVsBXRandoms_[kRPotsTotalNumber_]
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
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:212
void add(std::string const &label, ParameterSetDescription const &psetDescription)
Base class for CTPPS detector IDs.
Definition: CTPPSDetId.h:32
HLT enums.
int stationStatus_[kStationIDMAX_]
const std::string folderName_
static constexpr int kRPotsTotalNumber_
MonitorElement * book1D(TString const &name, TString const &title, int const nchX, double const lowX, double const highX, FUNC onbooking=NOOP())
Definition: DQMStore.h:98
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
Definition: event.py:1
Definition: Run.h:45