CMS 3D CMS Logo

AlcaPCCEventProducer.cc
Go to the documentation of this file.
1 
10 // C++ standard
11 #include <string>
12 // CMS
30 #include "TMath.h"
31 //The class
33 public:
34  explicit AlcaPCCEventProducer(const edm::ParameterSet&);
35  ~AlcaPCCEventProducer() override;
36  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
37 
38 private:
39  void produce(edm::Event& iEvent, const edm::EventSetup& iSetup) override;
40 
43 
44  std::string trigstring_; //specifies the trigger Rand or ZeroBias
45  int countEvt_; //counter
46  int countLumi_; //counter
47 
48  std::unique_ptr<reco::PixelClusterCountsInEvent> thePCCob;
49 };
50 
51 //--------------------------------------------------------------------------------------------------
53  fPixelClusterLabel = iConfig.getParameter<edm::InputTag>("pixelClusterLabel");
54  trigstring_ = iConfig.getUntrackedParameter<std::string>("trigstring", "alcaPCCEvent");
55  produces<reco::PixelClusterCountsInEvent, edm::Transition::Event>(trigstring_);
56  pixelToken = consumes<edmNew::DetSetVector<SiPixelCluster> >(fPixelClusterLabel);
57 }
58 
59 //--------------------------------------------------------------------------------------------------
61 
62 //--------------------------------------------------------------------------------------------------
64  countEvt_++;
65  thePCCob = std::make_unique<reco::PixelClusterCountsInEvent>();
66 
67  unsigned int bx = iEvent.bunchCrossing();
68 
69  //Looping over the clusters and adding the counts up
71  iEvent.getByToken(pixelToken, hClusterColl);
72 
73  const edmNew::DetSetVector<SiPixelCluster>& clustColl = *(hClusterColl.product());
74  // ----------------------------------------------------------------------
75  // -- Clusters without tracks
76  for (auto const& mod : clustColl) {
77  if (mod.empty()) {
78  continue;
79  }
80  DetId detId = mod.id();
81 
82  //--The following will be used when we make a theshold for the clusters.
83  //--Keeping this for features that may be implemented later.
84  // -- clusters on this det
85  //edmNew::DetSet<SiPixelCluster>::const_iterator di;
86  //int nClusterCount=0;
87  //for (di = mod.begin(); di != mod.end(); ++di) {
88  // nClusterCount++;
89  //}
90  int nCluster = mod.size();
91  thePCCob->increment(detId(), nCluster);
92  thePCCob->setbxID(bx);
93  }
94 
96 }
97 
98 //--------------------------------------------------------------------------------------------------
100  edm::ParameterSetDescription evtParamDesc;
101  evtParamDesc.add<edm::InputTag>("pixelClusterLabel", edm::InputTag("siPixelClustersForLumi"));
102  evtParamDesc.addUntracked<std::string>("trigstring", "alcaPCCEvent");
103  descriptions.add("alcaPCCEventProducer", evtParamDesc);
104 }
105 
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
ParameterDescriptionBase * addUntracked(U const &iLabel, T const &value)
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
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
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