CMS 3D CMS Logo

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