CMS 3D CMS Logo

SiStripMonitorHLT.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: SiStripMonitorCluster
4 // Class: SiStripMonitorHLT
5 //
6 // class SiStripMonitorHLT SiStripMonitorHLT.cc
7 // DQM/SiStripMonitorCluster/src/SiStripMonitorHLT.cc
8 #include <vector>
9 
10 #include <iostream>
11 #include <numeric>
12 
15 
17 
20 
22  HLTDirectory = "HLTResults";
24  conf_ = iConfig;
25 
27  consumes<int>(conf_.getParameter<std::string>("HLTProducer"));
29  consumes<uint>(conf_.getParameter<std::string>("HLTProducer"));
31  consumes<std::map<uint, std::vector<SiStripCluster> > >(
32  conf_.getParameter<std::string>("HLTProducer"));
33 }
34 
36  const edm::Run& run,
37  const edm::EventSetup& es) {
39  std::string HLTProducer = conf_.getParameter<std::string>("HLTProducer");
40  HLTDecision = ibooker.book1D(HLTProducer + "_HLTDecision",
41  HLTProducer + "HLTDecision", 2, -0.5, 1.5);
42  // all
44  "SumOfClusterCharges_all", "SumOfClusterCharges_all", 50, 0, 2000);
46  ibooker.book1D("ChargeOfEachClusterTIB_all", "ChargeOfEachClusterTIB_all",
47  400, -0.5, 400.5);
49  ibooker.book1D("ChargeOfEachClusterTOB_all", "ChargeOfEachClusterTOB_all",
50  400, -0.5, 400.5);
52  ibooker.book1D("ChargeOfEachClusterTEC_all", "ChargeOfEachClusterTEC_all",
53  400, -0.5, 400.5);
55  ibooker.book1D("NumberOfClustersAboveThreshold_all",
56  "NumberOfClustersAboveThreshold_all", 30, 30.5, 60.5);
57  // 31 = TIB2, 32 = TIB2, 33 = TIB3, 51 = TOB1, 52=TOB2, 60 = TEC
58  // accepted from HLT
60  "SumOfClusterCharges_hlt", "SumOfClusterCharges_hlt", 50, 0, 2000);
62  ibooker.book1D("ChargeOfEachClusterTIB_hlt", "ChargeOfEachClusterTIB_hlt",
63  400, -0.5, 400.5);
65  ibooker.book1D("ChargeOfEachClusterTOB_hlt", "ChargeOfEachClusterTOB_hlt",
66  400, -0.5, 400.5);
68  ibooker.book1D("ChargeOfEachClusterTEC_hlt", "ChargeOfEachClusterTEC_hlt",
69  400, -0.5, 400.5);
71  ibooker.book1D("NumberOfClustersAboveThreshold_hlt",
72  "NumberOfClustersAboveThreshold_hlt", 30, 30.5, 60.5);
73 }
74 
76  const edm::EventSetup& iSetup) {
77  // get from event
78  std::string HLTProducer = conf_.getParameter<std::string>("HLTProducer");
79  edm::Handle<int> filter_decision;
80  iEvent.getByToken(filerDecisionToken_, filter_decision); // filter decision
81  edm::Handle<uint> sum_of_clustch;
83  sum_of_clustch); // sum of cluster charges
84  // first element of pair: layer: TIB1, ...., TEC; second element: nr of
85  // clusters above threshold
87  clusters_in_subcomponents;
88  if (HLTProducer == "ClusterMTCCFilter")
89  iEvent.getByToken(clusterInSubComponentsToken_, clusters_in_subcomponents);
90 
91  // trigger decision
92  HLTDecision->Fill(*filter_decision);
93 
94  // sum of charges of clusters
95  SumOfClusterCharges_all->Fill(*sum_of_clustch);
96  if (*filter_decision) SumOfClusterCharges_hlt->Fill(*sum_of_clustch);
97 
98  // clusters in different layers
99  if (HLTProducer == "ClusterMTCCFilter") {
100  // loop over layers ("subcomponents")
101  for (std::map<uint, std::vector<SiStripCluster> >::const_iterator it =
102  clusters_in_subcomponents->begin();
103  it != clusters_in_subcomponents->end(); it++) {
104  int generalized_layer = it->first;
105  std::vector<SiStripCluster> theclusters = it->second;
107  generalized_layer,
108  theclusters.size()); // number of clusters in this generalized layer
109  if (*filter_decision)
110  NumberOfClustersAboveThreshold_hlt->Fill(generalized_layer,
111  theclusters.size());
112  // loop over clusters (and detids)
113  for (std::vector<SiStripCluster>::const_iterator icluster =
114  theclusters.begin();
115  icluster != theclusters.end(); icluster++) {
116  // calculate sum of amplitudes
117  unsigned int amplclus = 0;
118  for (auto ia = icluster->amplitudes().begin();
119  ia != icluster->amplitudes().end(); ia++) {
120  if ((*ia) > 0) amplclus += (*ia); // why should this be negative?
121  }
122  if (generalized_layer == 31 || generalized_layer == 32 ||
123  generalized_layer ==
124  33) { // you can also ask the detid here whether is TIB
125  ChargeOfEachClusterTIB_all->Fill(amplclus, 1.);
126  if (*filter_decision) ChargeOfEachClusterTIB_hlt->Fill(amplclus, 1.);
127  }
128  if (generalized_layer == 51 || generalized_layer == 52) {
129  ChargeOfEachClusterTOB_all->Fill(amplclus, 1.);
130  if (*filter_decision) ChargeOfEachClusterTOB_hlt->Fill(amplclus, 1.);
131  }
132  if (generalized_layer == 60) {
133  ChargeOfEachClusterTEC_all->Fill(amplclus, 1.);
134  if (*filter_decision) ChargeOfEachClusterTEC_hlt->Fill(amplclus, 1.);
135  }
136  }
137  }
138  }
139 }
T getParameter(std::string const &) const
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:579
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
MonitorElement * NumberOfClustersAboveThreshold_all
void Fill(long long x)
MonitorElement * SumOfClusterCharges_all
int iEvent
Definition: GenABIO.cc:230
MonitorElement * ChargeOfEachClusterTEC_hlt
void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:268
std::string HLTDirectory
MonitorElement * ChargeOfEachClusterTIB_all
MonitorElement * NumberOfClustersAboveThreshold_hlt
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:106
MonitorElement * HLTDecision
MonitorElement * ChargeOfEachClusterTEC_all
MonitorElement * SumOfClusterCharges_hlt
def uint(string)
MonitorElement * ChargeOfEachClusterTOB_hlt
SiStripMonitorHLT(const edm::ParameterSet &)
edm::EDGetTokenT< int > filerDecisionToken_
MonitorElement * ChargeOfEachClusterTIB_hlt
edm::EDGetTokenT< std::map< uint, std::vector< SiStripCluster > > > clusterInSubComponentsToken_
MonitorElement * ChargeOfEachClusterTOB_all
edm::ParameterSet conf_
void analyze(const edm::Event &, const edm::EventSetup &) override
Definition: Run.h:44
edm::EDGetTokenT< uint > sumOfClusterToken_