CMS 3D CMS Logo

ESRecoSummary.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: ESRecoSummary
4 // Class: ESRecoSummary
5 // Original Author: Martina Malberti
6 //
7 // system include files
8 #include <memory>
9 #include <iostream>
10 #include <cmath>
11 
12 // user include files
14 
16 
19 
21 
28 
38 
39 //
40 // constructors and destructor
41 //
43 {
44  prefixME_ = ps.getUntrackedParameter<std::string>("prefixME", "");
45 
46  //now do what ever initialization is needed
47  esRecHitCollection_ = consumes<ESRecHitCollection>(ps.getParameter<edm::InputTag>("recHitCollection_ES"));
48  esClusterCollectionX_ = consumes<reco::PreshowerClusterCollection>(ps.getParameter<edm::InputTag>("ClusterCollectionX_ES"));
49  esClusterCollectionY_ = consumes<reco::PreshowerClusterCollection>(ps.getParameter<edm::InputTag>("ClusterCollectionY_ES"));
50 
51  superClusterCollection_EE_ = consumes<reco::SuperClusterCollection>(ps.getParameter<edm::InputTag>("superClusterCollection_EE"));
52 }
53 
54 void
56 {
57  // Monitor Elements (ex THXD)
58  iBooker.setCurrentFolder(prefixME_ + "/ESRecoSummary"); // to organise the histos in folders
59 
60 
61  // Preshower ----------------------------------------------
62  h_recHits_ES_energyMax = iBooker.book1D("recHits_ES_energyMax","recHits_ES_energyMax",200,0.,0.01);
63  h_recHits_ES_time = iBooker.book1D("recHits_ES_time","recHits_ES_time",200,-100.,100.);
64 
65  h_esClusters_energy_plane1 = iBooker.book1D("esClusters_energy_plane1","esClusters_energy_plane1",200,0.,0.01);
66  h_esClusters_energy_plane2 = iBooker.book1D("esClusters_energy_plane2","esClusters_energy_plane2",200,0.,0.01);
67  h_esClusters_energy_ratio = iBooker.book1D("esClusters_energy_ratio","esClusters_energy_ratio",200,0.,20.);
68 }
69 
70 //
71 // member functions
72 //
73 
74 // ------------ method called to for each event ------------
76 {
77  //Preshower RecHits
79  ev.getByToken (esRecHitCollection_, recHitsES) ;
80  const ESRecHitCollection* thePreShowerRecHits = recHitsES.product () ;
81 
82  if ( ! recHitsES.isValid() ) {
83  std::cerr << "ESRecoSummary::analyze --> recHitsES not found" << std::endl;
84  }
85 
86  float maxRecHitEnergyES = -999.;
87 
88  for (ESRecHitCollection::const_iterator esItr = thePreShowerRecHits->begin(); esItr != thePreShowerRecHits->end(); ++esItr)
89  {
90 
91  h_recHits_ES_time -> Fill(esItr->time());
92  if (esItr -> energy() > maxRecHitEnergyES ) maxRecHitEnergyES = esItr -> energy() ;
93 
94  } // end loop over ES rec Hits
95 
96  h_recHits_ES_energyMax -> Fill(maxRecHitEnergyES );
97 
98  // ES clusters in X plane
100  ev.getByToken( esClusterCollectionX_, esClustersX);
101  const reco::PreshowerClusterCollection *ESclustersX = esClustersX.product();
102 
103  // ES clusters in Y plane
105  ev.getByToken( esClusterCollectionY_, esClustersY);
106  const reco::PreshowerClusterCollection *ESclustersY = esClustersY.product();
107 
108 
109  // ... endcap
111  ev.getByToken( superClusterCollection_EE_, superClusters_EE_h );
112  const reco::SuperClusterCollection* theEndcapSuperClusters = superClusters_EE_h.product () ;
113  if ( ! superClusters_EE_h.isValid() ) {
114  std::cerr << "EcalRecHitSummary::analyze --> superClusters_EE_h not found" << std::endl;
115  }
116 
117  // loop over all super clusters
118  for (reco::SuperClusterCollection::const_iterator itSC = theEndcapSuperClusters->begin();
119  itSC != theEndcapSuperClusters->end(); ++itSC ) {
120 
121  if ( fabs(itSC->eta()) < 1.65 || fabs(itSC->eta()) > 2.6 ) continue;
122 
123  float ESenergyPlane1 = 0.;
124  float ESenergyPlane2 = 0.;
125 
126 
127  // Loop over all ECAL Basic clusters in the supercluster
128  for (reco::CaloCluster_iterator ecalBasicCluster = itSC->clustersBegin(); ecalBasicCluster!= itSC->clustersEnd();
129  ecalBasicCluster++) {
130  const reco::CaloClusterPtr ecalBasicClusterPtr = *(ecalBasicCluster);
131 
132  for (reco::PreshowerClusterCollection::const_iterator iESClus = ESclustersX->begin(); iESClus != ESclustersX->end();
133  ++iESClus) {
134  const reco::CaloClusterPtr preshBasicCluster = iESClus->basicCluster();
135  const reco::PreshowerCluster *esCluster = &*iESClus;
136  if (preshBasicCluster == ecalBasicClusterPtr) {
137  ESenergyPlane1 += esCluster->energy();
138  }
139  } // end of x loop
140 
141  for (reco::PreshowerClusterCollection::const_iterator iESClus = ESclustersY->begin(); iESClus != ESclustersY->end();
142  ++iESClus) {
143  const reco::CaloClusterPtr preshBasicCluster = iESClus->basicCluster();
144  const reco::PreshowerCluster *esCluster = &*iESClus;
145  if (preshBasicCluster == ecalBasicClusterPtr) {
146  ESenergyPlane2 += esCluster->energy();
147  }
148  } // end of y loop
149  } // end loop over all basic clusters in the supercluster
150 
151  //cout<<"DQM : "<<ESenergyPlane1<<" "<<ESenergyPlane2<<endl;
152  h_esClusters_energy_plane1->Fill(ESenergyPlane1);
153  h_esClusters_energy_plane2->Fill(ESenergyPlane2);
154  if (ESenergyPlane1 > 0 && ESenergyPlane2 > 0) h_esClusters_energy_ratio -> Fill(ESenergyPlane1/ESenergyPlane2);
155 
156  }// end loop over superclusters
157 
158 }
159 
160 //define this as a plug-in
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:579
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
std::vector< EcalRecHit >::const_iterator const_iterator
bool ev
MonitorElement * h_esClusters_energy_plane1
Definition: ESRecoSummary.h:30
MonitorElement * h_recHits_ES_time
Definition: ESRecoSummary.h:28
edm::EDGetTokenT< ESRecHitCollection > esRecHitCollection_
Definition: ESRecoSummary.h:39
MonitorElement * h_esClusters_energy_ratio
Definition: ESRecoSummary.h:32
void Fill(long long x)
ESRecoSummary(const edm::ParameterSet &)
edm::EDGetTokenT< reco::PreshowerClusterCollection > esClusterCollectionY_
Definition: ESRecoSummary.h:41
std::string prefixME_
Definition: ESRecoSummary.h:24
void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:268
std::vector< SuperCluster > SuperClusterCollection
collection of SuperCluser objectr
MonitorElement * h_recHits_ES_energyMax
Definition: ESRecoSummary.h:27
void Fill(HcalDetId &id, double val, std::vector< TH2F > &depth)
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:106
std::vector< PreshowerCluster > PreshowerClusterCollection
collection of PreshowerCluster objects
double energy() const
cluster energy
Definition: CaloCluster.h:126
bool isValid() const
Definition: HandleBase.h:74
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
const_iterator end() const
T const * product() const
Definition: Handle.h:81
void analyze(const edm::Event &, const edm::EventSetup &) override
edm::EDGetTokenT< reco::PreshowerClusterCollection > esClusterCollectionX_
Definition: ESRecoSummary.h:40
MonitorElement * h_esClusters_energy_plane2
Definition: ESRecoSummary.h:31
edm::EDGetTokenT< reco::SuperClusterCollection > superClusterCollection_EE_
Definition: ESRecoSummary.h:38
const_iterator begin() const
Definition: Run.h:44