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 DQM/SiStripMonitorCluster/src/SiStripMonitorHLT.cc
7 #include <vector>
8 
9 #include <numeric>
10 #include <iostream>
11 
14 
15 
17 
20 
21 
23 {
24  HLTDirectory="HLTResults";
26  conf_ = iConfig;
27 
28  filerDecisionToken_ = consumes<int>(conf_.getParameter<std::string>("HLTProducer") );
29  sumOfClusterToken_ = consumes<uint>(conf_.getParameter<std::string>("HLTProducer") );
30  clusterInSubComponentsToken_ = consumes<std::map<uint,std::vector<SiStripCluster> > >(conf_.getParameter<std::string>("HLTProducer") );
31 }
32 
34 {
36  std::string HLTProducer = conf_.getParameter<std::string>("HLTProducer");
37  HLTDecision = ibooker.book1D(HLTProducer+"_HLTDecision", HLTProducer+"HLTDecision", 2, -0.5, 1.5);
38  // all
39  SumOfClusterCharges_all = ibooker.book1D("SumOfClusterCharges_all", "SumOfClusterCharges_all", 50, 0, 2000);
40  ChargeOfEachClusterTIB_all = ibooker.book1D("ChargeOfEachClusterTIB_all", "ChargeOfEachClusterTIB_all", 400, -0.5, 400.5);
41  ChargeOfEachClusterTOB_all = ibooker.book1D("ChargeOfEachClusterTOB_all", "ChargeOfEachClusterTOB_all", 400, -0.5, 400.5);
42  ChargeOfEachClusterTEC_all = ibooker.book1D("ChargeOfEachClusterTEC_all", "ChargeOfEachClusterTEC_all", 400, -0.5, 400.5);
43  NumberOfClustersAboveThreshold_all = ibooker.book1D("NumberOfClustersAboveThreshold_all", "NumberOfClustersAboveThreshold_all", 30, 30.5, 60.5);
44  // 31 = TIB2, 32 = TIB2, 33 = TIB3, 51 = TOB1, 52=TOB2, 60 = TEC
45  // accepted from HLT
46  SumOfClusterCharges_hlt = ibooker.book1D("SumOfClusterCharges_hlt", "SumOfClusterCharges_hlt", 50, 0, 2000);
47  ChargeOfEachClusterTIB_hlt = ibooker.book1D("ChargeOfEachClusterTIB_hlt", "ChargeOfEachClusterTIB_hlt", 400, -0.5, 400.5);
48  ChargeOfEachClusterTOB_hlt = ibooker.book1D("ChargeOfEachClusterTOB_hlt", "ChargeOfEachClusterTOB_hlt", 400, -0.5, 400.5);
49  ChargeOfEachClusterTEC_hlt = ibooker.book1D("ChargeOfEachClusterTEC_hlt", "ChargeOfEachClusterTEC_hlt", 400, -0.5, 400.5);
50  NumberOfClustersAboveThreshold_hlt = ibooker.book1D("NumberOfClustersAboveThreshold_hlt", "NumberOfClustersAboveThreshold_hlt", 30, 30.5, 60.5);
51 
52 }
53 
54 
56 {
57 
58  // get from event
59  std::string HLTProducer = conf_.getParameter<std::string>("HLTProducer");
60  edm::Handle<int> filter_decision; iEvent.getByToken(filerDecisionToken_,filter_decision); // filter decision
61  edm::Handle<uint> sum_of_clustch; iEvent.getByToken(sumOfClusterToken_, sum_of_clustch); // sum of cluster charges
62  // first element of pair: layer: TIB1, ...., TEC; second element: nr of clusters above threshold
63  edm::Handle<std::map<uint,std::vector<SiStripCluster> > > clusters_in_subcomponents;
64  if(HLTProducer=="ClusterMTCCFilter") iEvent.getByToken(clusterInSubComponentsToken_,clusters_in_subcomponents);
65 
66  // trigger decision
67  HLTDecision->Fill(*filter_decision);
68 
69  // sum of charges of clusters
70  SumOfClusterCharges_all->Fill(*sum_of_clustch);
71  if(*filter_decision) SumOfClusterCharges_hlt->Fill(*sum_of_clustch);
72 
73  //clusters in different layers
74  if(HLTProducer=="ClusterMTCCFilter"){
75  // loop over layers ("subcomponents")
76  for(std::map<uint,std::vector<SiStripCluster> >::const_iterator it = clusters_in_subcomponents->begin(); it != clusters_in_subcomponents->end(); it++){
77  int generalized_layer = it->first;
78  std::vector<SiStripCluster> theclusters = it->second;
79  NumberOfClustersAboveThreshold_all->Fill( generalized_layer, theclusters.size() ); // number of clusters in this generalized layer
80  if(*filter_decision) NumberOfClustersAboveThreshold_hlt->Fill( generalized_layer, theclusters.size() );
81  //loop over clusters (and detids)
82  for(std::vector<SiStripCluster>::const_iterator icluster = theclusters.begin(); icluster != theclusters.end(); icluster++){
83  // calculate sum of amplitudes
84  unsigned int amplclus=0;
85  for(auto ia=icluster->amplitudes().begin(); ia!=icluster->amplitudes().end(); ia++) {
86  if ((*ia)>0) amplclus+=(*ia); // why should this be negative?
87  }
88  if(generalized_layer==31 || generalized_layer==32 || generalized_layer==33){ // you can also ask the detid here whether is TIB
89  ChargeOfEachClusterTIB_all->Fill(amplclus,1.);
90  if(*filter_decision) ChargeOfEachClusterTIB_hlt->Fill(amplclus,1.);
91  }
92  if(generalized_layer==51 || generalized_layer==52){
93  ChargeOfEachClusterTOB_all->Fill(amplclus,1.);
94  if(*filter_decision) ChargeOfEachClusterTOB_hlt->Fill(amplclus,1.);
95  }
96  if(generalized_layer==60 ){
97  ChargeOfEachClusterTEC_all->Fill(amplclus,1.);
98  if(*filter_decision) ChargeOfEachClusterTEC_hlt->Fill(amplclus,1.);
99  }
100  }
101  }
102  }
103 }
T getParameter(std::string const &) const
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:508
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
std::string HLTDirectory
MonitorElement * ChargeOfEachClusterTIB_all
MonitorElement * NumberOfClustersAboveThreshold_hlt
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:115
MonitorElement * HLTDecision
MonitorElement * ChargeOfEachClusterTEC_all
MonitorElement * SumOfClusterCharges_hlt
edm::EDGetTokenT< std::map< uint, std::vector< SiStripCluster > > > clusterInSubComponentsToken_
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:277
MonitorElement * ChargeOfEachClusterTOB_hlt
SiStripMonitorHLT(const edm::ParameterSet &)
edm::EDGetTokenT< int > filerDecisionToken_
MonitorElement * ChargeOfEachClusterTIB_hlt
MonitorElement * ChargeOfEachClusterTOB_all
edm::ParameterSet conf_
void analyze(const edm::Event &, const edm::EventSetup &) override
Definition: Run.h:43
edm::EDGetTokenT< uint > sumOfClusterToken_