00001
00002
00003
00004
00005
00006
00007
00008 #include <memory>
00009
00010
00011 #include "FWCore/Framework/interface/Frameworkfwd.h"
00012 #include "FWCore/Framework/interface/EDAnalyzer.h"
00013
00014 #include "FWCore/Framework/interface/Event.h"
00015 #include "FWCore/Framework/interface/MakerMacros.h"
00016 #include "FWCore/Common/interface/EventBase.h"
00017 #include "DataFormats/Provenance/interface/EventAuxiliary.h"
00018
00019 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00020 #include "FWCore/ServiceRegistry/interface/Service.h"
00021 #include "CommonTools/UtilAlgos/interface/TFileService.h"
00022
00023
00024 #include "Geometry/CaloTopology/interface/CaloTopology.h"
00025 #include "Geometry/CaloGeometry/interface/CaloGeometry.h"
00026 #include "Geometry/CaloGeometry/interface/CaloSubdetectorGeometry.h"
00027 #include "Geometry/CaloEventSetup/interface/CaloTopologyRecord.h"
00028 #include "Geometry/CaloTopology/interface/CaloSubdetectorTopology.h"
00029 #include "DataFormats/GeometryVector/interface/GlobalPoint.h"
00030 #include "Geometry/EcalAlgo/interface/EcalBarrelGeometry.h"
00031 #include "Geometry/EcalAlgo/interface/EcalEndcapGeometry.h"
00032
00033 #include "DataFormats/EcalDetId/interface/EBDetId.h"
00034 #include "DataFormats/EcalDetId/interface/EEDetId.h"
00035 #include "DataFormats/EcalDetId/interface/ESDetId.h"
00036 #include "DataFormats/EcalRecHit/interface/EcalRecHitCollections.h"
00037 #include "DataFormats/EgammaReco/interface/SuperCluster.h"
00038 #include "DataFormats/EgammaReco/interface/SuperClusterFwd.h"
00039 #include "DataFormats/EgammaReco/interface/BasicCluster.h"
00040 #include "DataFormats/EgammaReco/interface/BasicClusterFwd.h"
00041 #include "RecoEcal/EgammaCoreTools/interface/EcalClusterTools.h"
00042 #include "RecoEcal/EgammaCoreTools/interface/EcalTools.h"
00043 #include "RecoLocalCalo/EcalRecAlgos/interface/EcalSeverityLevelAlgo.h"
00044 #include "RecoLocalCalo/EcalRecAlgos/interface/EcalCleaningAlgo.h"
00045 #include "RecoEcal/EgammaCoreTools/interface/EcalRecHitLess.h"
00046
00047 #include "DataFormats/TrackReco/interface/Track.h"
00048 #include "DataFormats/TrackReco/interface/TrackFwd.h"
00049
00050 #include "DataFormats/JetReco/interface/CaloJet.h"
00051 #include "DataFormats/JetReco/interface/CaloJetCollection.h"
00052
00053 #include "DQMOffline/Ecal/interface/EBRecoSummary.h"
00054 #include "DataFormats/BeamSpot/interface/BeamSpot.h"
00055 #include "MagneticField/Records/interface/IdealMagneticFieldRecord.h"
00056 #include "MagneticField/Engine/interface/MagneticField.h"
00057
00058 #include "TVector3.h"
00059
00060 #include <iostream>
00061 #include <cmath>
00062 #include <fstream>
00063 #include <string>
00064
00065
00066
00067
00068 EBRecoSummary::EBRecoSummary(const edm::ParameterSet& ps)
00069 {
00070
00071 prefixME_ = ps.getUntrackedParameter<std::string>("prefixME", "");
00072
00073
00074 recHitCollection_EB_ = ps.getParameter<edm::InputTag>("recHitCollection_EB");
00075 redRecHitCollection_EB_ = ps.getParameter<edm::InputTag>("redRecHitCollection_EB");
00076 basicClusterCollection_EB_ = ps.getParameter<edm::InputTag>("basicClusterCollection_EB");
00077 superClusterCollection_EB_ = ps.getParameter<edm::InputTag>("superClusterCollection_EB");
00078
00079 ethrEB_ = ps.getParameter<double>("ethrEB");
00080
00081 scEtThrEB_ = ps.getParameter<double>("scEtThrEB");
00082
00083
00084 dqmStore_= edm::Service<DQMStore>().operator->();
00085
00086
00087 dqmStore_->setCurrentFolder(prefixME_ + "/EBRecoSummary");
00088
00089
00090
00091 h_redRecHits_EB_recoFlag = dqmStore_->book1D("redRecHits_EB_recoFlag","redRecHits_EB_recoFlag",16,-0.5,15.5);
00092
00093
00094
00095 h_recHits_EB_energyMax = dqmStore_->book1D("recHits_EB_energyMax","recHitsEB_energyMax",110,-10,100);
00096 h_recHits_EB_Chi2 = dqmStore_->book1D("recHits_EB_Chi2","recHits_EB_Chi2",200,0,100);
00097 h_recHits_EB_time = dqmStore_->book1D("recHits_EB_time","recHits_EB_time",200,-50,50);
00098 h_recHits_EB_E1oE4 = dqmStore_->book1D("recHits_EB_E1oE4","recHitsEB_E1oE4",150, 0, 1.5);
00099 h_recHits_EB_recoFlag = dqmStore_->book1D("recHits_EB_recoFlag","recHits_EB_recoFlag",16,-0.5,15.5);
00100
00101
00102
00103 h_basicClusters_recHits_EB_recoFlag = dqmStore_->book1D("basicClusters_recHits_EB_recoFlag","basicClusters_recHits_EB_recoFlag",16,-0.5,15.5);
00104
00105
00106
00107 h_superClusters_EB_nBC = dqmStore_->book1D("superClusters_EB_nBC","superClusters_EB_nBC",100,0.,100.);
00108 h_superClusters_EB_E1oE4 = dqmStore_->book1D("superClusters_EB_E1oE4","superClusters_EB_E1oE4",150,0,1.5);
00109
00110 h_superClusters_eta = dqmStore_->book1D("superClusters_eta","superClusters_eta",150,-3.,3.);
00111 h_superClusters_EB_phi = dqmStore_->book1D("superClusters_EB_phi","superClusters_EB_phi",360,-3.1415927,3.1415927);
00112
00113 }
00114
00115
00116
00117 EBRecoSummary::~EBRecoSummary()
00118 {
00119
00120
00121 }
00122
00123
00124
00125
00126
00127
00128
00129 void EBRecoSummary::analyze(const edm::Event& ev, const edm::EventSetup& iSetup)
00130 {
00131
00132
00133 edm::ESHandle<MagneticField> theMagField;
00134 iSetup.get<IdealMagneticFieldRecord>().get(theMagField);
00135
00136
00137 edm::ESHandle<CaloTopology> pTopology;
00138 iSetup.get<CaloTopologyRecord>().get(pTopology);
00139 const CaloTopology *topology = pTopology.product();
00140
00141
00142 edm::Handle<EcalRecHitCollection> redRecHitsEB;
00143 ev.getByLabel( redRecHitCollection_EB_, redRecHitsEB );
00144 const EcalRecHitCollection* theBarrelEcalredRecHits = redRecHitsEB.product () ;
00145 if ( ! redRecHitsEB.isValid() ) {
00146 edm::LogWarning("EBRecoSummary") << "redRecHitsEB not found";
00147 }
00148
00149 for ( EcalRecHitCollection::const_iterator itr = theBarrelEcalredRecHits->begin () ;
00150 itr != theBarrelEcalredRecHits->end () ;++itr)
00151 {
00152
00153 h_redRecHits_EB_recoFlag->Fill( itr -> recoFlag() );
00154
00155 }
00156
00157
00158
00159
00160 edm::Handle<EcalRecHitCollection> recHitsEB;
00161 ev.getByLabel( recHitCollection_EB_, recHitsEB );
00162 const EcalRecHitCollection* theBarrelEcalRecHits = recHitsEB.product () ;
00163 if ( ! recHitsEB.isValid() ) {
00164 edm::LogWarning("EBRecoSummary") << "recHitsEB not found";
00165 }
00166
00167 float maxRecHitEnergyEB = -999.;
00168
00169 EBDetId ebid_MrecHitEB;
00170
00171 for ( EcalRecHitCollection::const_iterator itr = theBarrelEcalRecHits->begin () ;
00172 itr != theBarrelEcalRecHits->end () ;++itr)
00173 {
00174
00175 EBDetId ebid( itr -> id() );
00176
00177 h_recHits_EB_recoFlag -> Fill( itr -> recoFlag() );
00178
00179
00180 if (itr -> energy() > maxRecHitEnergyEB ){
00181 maxRecHitEnergyEB = itr -> energy() ;
00182 }
00183
00184 if ( itr -> energy() > ethrEB_ ){
00185 h_recHits_EB_Chi2 -> Fill( itr -> chi2() );
00186 h_recHits_EB_time -> Fill( itr -> time() );
00187 }
00188
00189 float R4 = EcalTools::swissCross( ebid, *theBarrelEcalRecHits, 0. );
00190
00191 if ( itr -> energy() > 3. && abs(ebid.ieta())!=85 ) h_recHits_EB_E1oE4-> Fill( R4 );
00192
00193 }
00194
00195 h_recHits_EB_energyMax -> Fill( maxRecHitEnergyEB );
00196
00197
00198
00199
00200 edm::Handle<reco::BasicClusterCollection> basicClusters_EB_h;
00201 ev.getByLabel( basicClusterCollection_EB_, basicClusters_EB_h );
00202 const reco::BasicClusterCollection* theBarrelBasicClusters = basicClusters_EB_h.product () ;
00203 if ( ! basicClusters_EB_h.isValid() ) {
00204 edm::LogWarning("EBRecoSummary") << "basicClusters_EB_h not found";
00205 }
00206
00207 for (reco::BasicClusterCollection::const_iterator itBC = theBarrelBasicClusters->begin();
00208 itBC != theBarrelBasicClusters->end(); ++itBC ) {
00209
00210
00211 const std::vector<std::pair<DetId,float> > & hits= itBC->hitsAndFractions();
00212 for (std::vector<std::pair<DetId,float> > ::const_iterator rh = hits.begin(); rh!=hits.end(); ++rh){
00213
00214 EBRecHitCollection::const_iterator itrechit = theBarrelEcalRecHits->find((*rh).first);
00215 if (itrechit==theBarrelEcalRecHits->end()) continue;
00216 h_basicClusters_recHits_EB_recoFlag -> Fill ( itrechit -> recoFlag() );
00217
00218 }
00219
00220 }
00221
00222
00223
00224 edm::Handle<reco::SuperClusterCollection> superClusters_EB_h;
00225 ev.getByLabel( superClusterCollection_EB_, superClusters_EB_h );
00226 const reco::SuperClusterCollection* theBarrelSuperClusters = superClusters_EB_h.product () ;
00227 if ( ! superClusters_EB_h.isValid() ) {
00228 edm::LogWarning("EBRecoSummary") << "superClusters_EB_h not found";
00229 }
00230
00231 for (reco::SuperClusterCollection::const_iterator itSC = theBarrelSuperClusters->begin();
00232 itSC != theBarrelSuperClusters->end(); ++itSC ) {
00233
00234 double scEt = itSC -> energy() * sin(2.*atan( exp(- itSC->position().eta() )));
00235
00236 if (scEt < scEtThrEB_ ) continue;
00237
00238 h_superClusters_EB_nBC -> Fill( itSC -> clustersSize());
00239 h_superClusters_eta -> Fill( itSC -> eta() );
00240 h_superClusters_EB_phi -> Fill( itSC -> phi() );
00241
00242 float E1 = EcalClusterTools::eMax ( *itSC, theBarrelEcalRecHits);
00243 float E4 = EcalClusterTools::eTop ( *itSC, theBarrelEcalRecHits, topology )+
00244 EcalClusterTools::eRight ( *itSC, theBarrelEcalRecHits, topology )+
00245 EcalClusterTools::eBottom( *itSC, theBarrelEcalRecHits, topology )+
00246 EcalClusterTools::eLeft ( *itSC, theBarrelEcalRecHits, topology );
00247
00248 if ( E1 > 3. ) h_superClusters_EB_E1oE4 -> Fill( 1.- E4/E1);
00249
00250 }
00251
00252 }
00253
00254
00255
00256 void
00257 EBRecoSummary::beginJob()
00258 {
00259 }
00260
00261
00262 void
00263 EBRecoSummary::endJob()
00264 {}
00265
00266
00267
00268 void EBRecoSummary::convxtalid(Int_t &nphi,Int_t &neta)
00269 {
00270
00271
00272
00273
00274 if(neta > 0) neta -= 1;
00275 if(nphi > 359) nphi=nphi-360;
00276
00277 }
00278
00279 int EBRecoSummary::diff_neta_s(Int_t neta1, Int_t neta2){
00280 Int_t mdiff;
00281 mdiff=(neta1-neta2);
00282 return mdiff;
00283 }
00284
00285
00286 int EBRecoSummary::diff_nphi_s(Int_t nphi1,Int_t nphi2) {
00287 Int_t mdiff;
00288 if(abs(nphi1-nphi2) < (360-abs(nphi1-nphi2))) {
00289 mdiff=nphi1-nphi2;
00290 }
00291 else {
00292 mdiff=360-abs(nphi1-nphi2);
00293 if(nphi1>nphi2) mdiff=-mdiff;
00294 }
00295 return mdiff;
00296 }
00297
00298
00299 DEFINE_FWK_MODULE(EBRecoSummary);