CMS 3D CMS Logo

AlcaPCCEventProducer.cc
Go to the documentation of this file.
1 
7 // C++ standard
8 #include <string>
9 
10 // CMS
28 #include "TMath.h"
29 
30 //The class
32 public:
33  explicit AlcaPCCEventProducer(const edm::ParameterSet&);
34  ~AlcaPCCEventProducer() override = default;
35  void produce(edm::StreamID id, edm::Event& e, edm::EventSetup const& c) const final;
36  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
37 
38 private:
40  const std::string trigstring_; //specifies the trigger Rand or ZeroBias
41  const bool savePerROCInfo_; // save per ROC data (important for the special fills)
43 
44  static constexpr int rowsperroc = 52;
45  static constexpr int colsperroc = 80;
46  static constexpr int nROCcolumns = 8;
47 };
48 
49 //--------------------------------------------------------------------------------------------------
51  : pixelClusterLabel_(iConfig.getParameter<edm::InputTag>("pixelClusterLabel")),
52  trigstring_(iConfig.getUntrackedParameter<std::string>("trigstring", "alcaPCCEvent")),
53  savePerROCInfo_(iConfig.getParameter<bool>("savePerROCInfo")),
54  pixelToken_(consumes<edmNew::DetSetVector<SiPixelCluster> >(pixelClusterLabel_)) {
55  produces<reco::PixelClusterCountsInEvent, edm::Transition::Event>(trigstring_);
56 }
57 
58 //--------------------------------------------------------------------------------------------------
60  std::unique_ptr<reco::PixelClusterCountsInEvent> thePCCob = std::make_unique<reco::PixelClusterCountsInEvent>();
61  unsigned int bx = iEvent.bunchCrossing();
62 
63  //Looping over the clusters and adding the counts up
65  iEvent.getByToken(pixelToken_, hClusterColl);
66 
67  const edmNew::DetSetVector<SiPixelCluster>& clustColl = *(hClusterColl.product());
68  // ----------------------------------------------------------------------
69  // -- Clusters without tracks
70  for (auto const& mod : clustColl) {
71  if (mod.empty()) {
72  continue;
73  }
74  DetId detId = mod.id();
75 
76  if (savePerROCInfo_) {
77  // Iterate over Clusters in module to fill per ROC histogram
78  for (auto const& cluster : mod) {
79  for (int i = 0; i < cluster.size(); ++i) {
80  const auto pix = cluster.pixel(i);
81  // TODO: add roc threshold to config if(di.adc > fRocThreshold_) {
82  if (pix.adc > 0) {
83  int irow = pix.x / rowsperroc; /* constant column direction is along x-axis */
84  int icol = pix.y / colsperroc; /* constant row direction is along y-axis */
85  /* generate the folling roc index that is going to map with ROC id as
86  8 9 10 11 12 13 14 15
87  0 1 2 3 4 5 6 7 */
88  int key = icol + irow * nROCcolumns;
89  thePCCob->incrementRoc(((detId << 7) + key), 1);
90  }
91  }
92  }
93  }
94 
95  int nCluster = mod.size();
96  thePCCob->increment(detId(), nCluster);
97  thePCCob->setbxID(bx);
98  }
99 
100  iEvent.put(std::move(thePCCob), std::string(trigstring_));
101 }
102 
103 //--------------------------------------------------------------------------------------------------
105  edm::ParameterSetDescription evtParamDesc;
106  evtParamDesc.add<edm::InputTag>("pixelClusterLabel", edm::InputTag("siPixelClustersForLumi"));
107  evtParamDesc.addUntracked<std::string>("trigstring", "alcaPCCEvent");
108  evtParamDesc.add<bool>("savePerROCInfo", true);
109  descriptions.add("alcaPCCEventProducer", evtParamDesc);
110 }
111 
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
ParameterDescriptionBase * addUntracked(U const &iLabel, T const &value)
const edm::InputTag pixelClusterLabel_
AlcaPCCEventProducer(const edm::ParameterSet &)
T const * product() const
Definition: Handle.h:70
void produce(edm::StreamID id, edm::Event &e, edm::EventSetup const &c) const final
static constexpr int colsperroc
int iEvent
Definition: GenABIO.cc:224
const edm::EDGetTokenT< edmNew::DetSetVector< SiPixelCluster > > pixelToken_
~AlcaPCCEventProducer() override=default
static constexpr int nROCcolumns
static constexpr int rowsperroc
key
prepare the HTCondor submission files and eventually submit them
const std::string trigstring_
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
ParameterDescriptionBase * add(U const &iLabel, T const &value)
Definition: DetId.h:17
void add(std::string const &label, ParameterSetDescription const &psetDescription)
Pixel cluster – collection of neighboring pixels above threshold.
HLT enums.
T mod(const T &a, const T &b)
Definition: ecalDccMap.h:4
def move(src, dest)
Definition: eostools.py:511