Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:09:27

0001 #include "DQMOffline/JetMET/interface/BeamHaloAnalyzer.h"
0002 //author : Ronny Remington, University of Florida
0003 //date : 11/11/09
0004 
0005 using namespace edm;
0006 using namespace reco;
0007 
0008 int Phi_To_iPhi(float phi) {
0009   phi = phi < 0 ? phi + 2. * TMath::Pi() : phi;
0010   float phi_degrees = phi * (360.) / (2. * TMath::Pi());
0011   int iPhi = (int)((phi_degrees / 5.) + 1.);
0012 
0013   return iPhi < 73 ? iPhi : 73;
0014 }
0015 
0016 BeamHaloAnalyzer::BeamHaloAnalyzer(const edm::ParameterSet& iConfig) {
0017   OutputFileName = iConfig.getParameter<std::string>("OutputFile");
0018   TextFileName = iConfig.getParameter<std::string>("TextFile");
0019 
0020   if (!TextFileName.empty())
0021     out = new std::ofstream(TextFileName.c_str());
0022 
0023   if (iConfig.exists(
0024           "StandardDQM"))  // If StandardDQM == true , coarse binning is used on selected (important) histograms
0025     StandardDQM = iConfig.getParameter<bool>("StandardDQM");
0026   else
0027     StandardDQM = false;
0028 
0029   //Get Input Tags
0030   //Digi Level
0031   IT_L1MuGMTReadout = iConfig.getParameter<edm::InputTag>("L1MuGMTReadoutLabel");
0032 
0033   //RecHit Level
0034   IT_CSCRecHit = consumes<CSCRecHit2DCollection>(iConfig.getParameter<edm::InputTag>("CSCRecHitLabel"));
0035   IT_EBRecHit = consumes<EBRecHitCollection>(iConfig.getParameter<edm::InputTag>("EBRecHitLabel"));
0036   IT_EERecHit = consumes<EERecHitCollection>(iConfig.getParameter<edm::InputTag>("EERecHitLabel"));
0037   IT_ESRecHit = consumes<ESRecHitCollection>(iConfig.getParameter<edm::InputTag>("ESRecHitLabel"));
0038   IT_HBHERecHit = consumes<HBHERecHitCollection>(iConfig.getParameter<edm::InputTag>("HBHERecHitLabel"));
0039   IT_HFRecHit = consumes<HFRecHitCollection>(iConfig.getParameter<edm::InputTag>("HFRecHitLabel"));
0040   IT_HORecHit = consumes<HORecHitCollection>(iConfig.getParameter<edm::InputTag>("HORecHitLabel"));
0041 
0042   //Higher Level Reco
0043   IT_CSCSegment = consumes<CSCSegmentCollection>(iConfig.getParameter<edm::InputTag>("CSCSegmentLabel"));
0044   IT_CosmicStandAloneMuon =
0045       consumes<reco::MuonCollection>(iConfig.getParameter<edm::InputTag>("CosmicStandAloneMuonLabel"));
0046   IT_BeamHaloMuon = consumes<reco::MuonCollection>(iConfig.getParameter<edm::InputTag>("BeamHaloMuonLabel"));
0047   IT_CollisionMuon = consumes<reco::MuonCollection>(iConfig.getParameter<edm::InputTag>("CollisionMuonLabel"));
0048   IT_CollisionStandAloneMuon =
0049       consumes<reco::MuonCollection>(iConfig.getParameter<edm::InputTag>("CollisionStandAloneMuonLabel"));
0050   IT_met = consumes<reco::CaloMETCollection>(iConfig.getParameter<edm::InputTag>("metLabel"));
0051   IT_CaloTower = consumes<edm::View<reco::Candidate> >(iConfig.getParameter<edm::InputTag>("CaloTowerLabel"));
0052   IT_SuperCluster = consumes<SuperClusterCollection>(iConfig.getParameter<edm::InputTag>("SuperClusterLabel"));
0053   IT_Photon = consumes<reco::PhotonCollection>(iConfig.getParameter<edm::InputTag>("PhotonLabel"));
0054 
0055   //Halo Data
0056   IT_CSCHaloData = consumes<reco::CSCHaloData>(iConfig.getParameter<edm::InputTag>("CSCHaloDataLabel"));
0057   IT_EcalHaloData = consumes<reco::EcalHaloData>(iConfig.getParameter<edm::InputTag>("EcalHaloDataLabel"));
0058   IT_HcalHaloData = consumes<reco::HcalHaloData>(iConfig.getParameter<edm::InputTag>("HcalHaloDataLabel"));
0059   IT_GlobalHaloData = consumes<reco::GlobalHaloData>(iConfig.getParameter<edm::InputTag>("GlobalHaloDataLabel"));
0060   IT_BeamHaloSummary = consumes<BeamHaloSummary>(iConfig.getParameter<edm::InputTag>("BeamHaloSummaryLabel"));
0061 
0062   cscGeomToken_ = esConsumes();
0063 
0064   edm::InputTag CosmicSAMuonLabel = iConfig.getParameter<edm::InputTag>("CosmicStandAloneMuonLabel");
0065   IT_CSCTimeMapToken = consumes<reco::MuonTimeExtraMap>(edm::InputTag(CosmicSAMuonLabel.label(), std::string("csc")));
0066   FolderName = iConfig.getParameter<std::string>("folderName");
0067   DumpMET = iConfig.getParameter<double>("DumpMET");
0068 
0069   //Muon to Segment Matching
0070   edm::ParameterSet matchParameters = iConfig.getParameter<edm::ParameterSet>("MatchParameters");
0071   edm::ConsumesCollector iC = consumesCollector();
0072   TheMatcher = new MuonSegmentMatcher(matchParameters, iC);
0073 }
0074 
0075 void BeamHaloAnalyzer::bookHistograms(DQMStore::IBooker& ibooker, edm::Run const& iRun, edm::EventSetup const&) {
0076   // EcalHaloData
0077   ibooker.setCurrentFolder(FolderName + "/EcalHaloData");
0078   if (StandardDQM) {
0079     hEcalHaloData_PhiWedgeMultiplicity = ibooker.book1D("EcalHaloData_PhiWedgeMultiplicity", "", 20, -0.5, 19.5);
0080     hEcalHaloData_PhiWedgeConstituents = ibooker.book1D("EcalHaloData_PhiWedgeConstituents", "", 20, -0.5, 19.5);
0081     //  hEcalHaloData_PhiWedgeiPhi         = ibooker.book1D("EcalHaloData_PhiWedgeiPhi","", 360, 0.5, 360.5) ;
0082     hEcalHaloData_PhiWedgeZDirectionConfidence =
0083         ibooker.book1D("EcalHaloData_ZDirectionConfidence", "", 120, -1.2, 1.2);
0084     hEcalHaloData_SuperClusterShowerShapes =
0085         ibooker.book2D("EcalHaloData_SuperClusterShowerShapes", "", 30, 0, 3.2, 25, 0.0, 2.0);
0086     hEcalHaloData_SuperClusterEnergy = ibooker.book1D("EcalHaloData_SuperClusterEnergy", "", 50, -0.5, 99.5);
0087     hEcalHaloData_SuperClusterNHits = ibooker.book1D("EcalHaloData_SuperClusterNHits", "", 20, -0.5, 19.5);
0088   } else {
0089     hEcalHaloData_PhiWedgeMultiplicity = ibooker.book1D("EcalHaloData_PhiWedgeMultiplicity", "", 20, -0.5, 19.5);
0090     hEcalHaloData_PhiWedgeEnergy = ibooker.book1D("EcalHaloData_PhiWedgeEnergy", "", 50, -0.5, 199.5);
0091     hEcalHaloData_PhiWedgeConstituents = ibooker.book1D("EcalHaloData_PhiWedgeConstituents", "", 20, -0.5, 19.5);
0092     hEcalHaloData_PhiWedgeMinTime = ibooker.book1D("EcalHaloData_PhiWedgeMinTime", "", 100, -225.0, 225.0);
0093     hEcalHaloData_PhiWedgeMaxTime = ibooker.book1D("EcalHaloData_PhiWedgeMaxTime", "", 100, -225.0, 225.0);
0094     hEcalHaloData_PhiWedgeiPhi = ibooker.book1D("EcalHaloData_PhiWedgeiPhi", "", 360, 0.5, 360.5);
0095     hEcalHaloData_PhiWedgePlusZDirectionConfidence =
0096         ibooker.book1D("EcalHaloData_PlusZDirectionConfidence", "", 50, 0., 1.0);
0097     hEcalHaloData_PhiWedgeZDirectionConfidence =
0098         ibooker.book1D("EcalHaloData_ZDirectionConfidence", "", 120, -1.2, 1.2);
0099     hEcalHaloData_PhiWedgeMinVsMaxTime =
0100         ibooker.book2D("EcalHaloData_PhiWedgeMinVsMaxTime", "", 50, -100.0, 100.0, 50, -100.0, 100.0);
0101     hEcalHaloData_SuperClusterShowerShapes =
0102         ibooker.book2D("EcalHaloData_SuperClusterShowerShapes", "", 30, 0, 3.2, 25, 0.0, 2.0);
0103     hEcalHaloData_SuperClusterEnergy = ibooker.book1D("EcalHaloData_SuperClusterEnergy", "", 100, -0.5, 99.5);
0104     hEcalHaloData_SuperClusterNHits = ibooker.book1D("EcalHaloData_SuperClusterNHits", "", 20, -0.5, 19.5);
0105     hEcalHaloData_SuperClusterPhiVsEta =
0106         ibooker.book2D("EcalHaloData_SuperClusterPhiVsEta", "", 60, -3.0, 3.0, 60, -3.2, 3.2);
0107   }
0108 
0109   // HcalHaloData
0110   ibooker.setCurrentFolder(FolderName + "/HcalHaloData");
0111   if (StandardDQM) {
0112     hHcalHaloData_PhiWedgeMultiplicity = ibooker.book1D("HcalHaloData_PhiWedgeMultiplicity", "", 20, -0.5, 19.5);
0113     hHcalHaloData_PhiWedgeConstituents = ibooker.book1D("HcalHaloData_PhiWedgeConstituents", "", 20, -0.5, 19.5);
0114     //hHcalHaloData_PhiWedgeiPhi         = ibooker.book1D("HcalHaloData_PhiWedgeiPhi","", 72, 0.5,72.5);
0115     hHcalHaloData_PhiWedgeZDirectionConfidence =
0116         ibooker.book1D("HcalHaloData_ZDirectionConfidence", "", 120, -1.2, 1.2);
0117   } else {
0118     hHcalHaloData_PhiWedgeMultiplicity = ibooker.book1D("HcalHaloData_PhiWedgeMultiplicity", "", 20, -0.5, 19.5);
0119     hHcalHaloData_PhiWedgeEnergy = ibooker.book1D("HcalHaloData_PhiWedgeEnergy", "", 50, -0.5, 199.5);
0120     hHcalHaloData_PhiWedgeConstituents = ibooker.book1D("HcalHaloData_PhiWedgeConstituents", "", 20, -0.5, 19.5);
0121     hHcalHaloData_PhiWedgeiPhi = ibooker.book1D("HcalHaloData_PhiWedgeiPhi", "", 72, 0.5, 72.5);
0122     hHcalHaloData_PhiWedgeMinTime = ibooker.book1D("HcalHaloData_PhiWedgeMinTime", "", 50, -100.0, 100.0);
0123     hHcalHaloData_PhiWedgeMaxTime = ibooker.book1D("HcalHaloData_PhiWedgeMaxTime", "", 50, -100.0, 100.0);
0124     hHcalHaloData_PhiWedgePlusZDirectionConfidence =
0125         ibooker.book1D("HcalHaloData_PlusZDirectionConfidence", "", 50, 0., 1.0);
0126     hHcalHaloData_PhiWedgeZDirectionConfidence =
0127         ibooker.book1D("HcalHaloData_ZDirectionConfidence", "", 120, -1.2, 1.2);
0128     hHcalHaloData_PhiWedgeMinVsMaxTime =
0129         ibooker.book2D("HcalHaloData_PhiWedgeMinVsMaxTime", "", 50, -100.0, 100.0, 50, -100.0, 100.0);
0130   }
0131 
0132   // CSCHaloData
0133   ibooker.setCurrentFolder(FolderName + "/CSCHaloData");
0134   if (StandardDQM) {
0135     hCSCHaloData_TrackMultiplicity = ibooker.book1D("CSCHaloData_TrackMultiplicity", "", 15, -0.5, 14.5);
0136     hCSCHaloData_TrackMultiplicityMEPlus = ibooker.book1D("CSCHaloData_TrackMultiplicityMEPlus", "", 15, -0.5, 14.5);
0137     hCSCHaloData_TrackMultiplicityMEMinus = ibooker.book1D("CSCHaloData_TrackMultiplicityMEMinus", "", 15, -0.5, 14.5);
0138     hCSCHaloData_InnerMostTrackHitR = ibooker.book1D("CSCHaloData_InnerMostTrackHitR", "", 70, 99.5, 799.5);
0139     hCSCHaloData_InnerMostTrackHitPhi = ibooker.book1D("CSCHaloData_InnerMostTrackHitPhi", "", 60, -3.2, 3.2);
0140     hCSCHaloData_L1HaloTriggersMEPlus = ibooker.book1D("CSCHaloData_L1HaloTriggersMEPlus", "", 10, -0.5, 9.5);
0141     hCSCHaloData_L1HaloTriggersMEMinus = ibooker.book1D("CSCHaloData_L1HaloTriggersMEMinus", "", 10, -0.5, 9.5);
0142     hCSCHaloData_L1HaloTriggers = ibooker.book1D("CSCHaloData_L1HaloTriggers", "", 10, -0.5, 9.5);
0143     hCSCHaloData_HLHaloTriggers = ibooker.book1D("CSCHaloData_HLHaloTriggers", "", 2, -0.5, 1.5);
0144     hCSCHaloData_NOutOfTimeTriggersvsL1HaloExists =
0145         ibooker.book2D("CSCHaloData_NOutOfTimeTriggersvsL1HaloExists", "", 20, -0.5, 19.5, 2, -0.5, 1.5);
0146     hCSCHaloData_NOutOfTimeTriggersMEPlus = ibooker.book1D("CSCHaloData_NOutOfTimeTriggersMEPlus", "", 20, -0.5, 19.5);
0147     hCSCHaloData_NOutOfTimeTriggersMEMinus =
0148         ibooker.book1D("CSCHaloData_NOutOfTimeTriggersMEMinus", "", 20, -0.5, 19.5);
0149     hCSCHaloData_NOutOfTimeTriggers = ibooker.book1D("CSCHaloData_NOutOfTimeTriggers", "", 20, -0.5, 19.5);
0150     hCSCHaloData_NOutOfTimeHits = ibooker.book1D("CSCHaloData_NOutOfTimeHits", "", 60, -0.5, 59.5);
0151     hCSCHaloData_NTracksSmalldT = ibooker.book1D("CSCHaloData_NTracksSmalldT", "", 15, -0.5, 14.5);
0152     hCSCHaloData_NTracksSmallBeta = ibooker.book1D("CSCHaloData_NTracksSmallBeta", "", 15, -0.5, 14.5);
0153     hCSCHaloData_NTracksSmallBetaAndSmalldT =
0154         ibooker.book1D("CSCHaloData_NTracksSmallBetaAndSmalldT", "", 15, -0.5, 14.5);
0155     hCSCHaloData_NTracksSmalldTvsNHaloTracks =
0156         ibooker.book2D("CSCHaloData_NTracksSmalldTvsNHaloTracks", "", 15, -0.5, 14.5, 15, -0.5, 14.5);
0157     hCSCHaloData_SegmentdT = ibooker.book1D("CSCHaloData_SegmentdT", "", 100, -100, 100);
0158     hCSCHaloData_FreeInverseBeta = ibooker.book1D("CSCHaloData_FreeInverseBeta", "", 80, -4, 4);
0159     hCSCHaloData_FreeInverseBetaVsSegmentdT =
0160         ibooker.book2D("CSCHaloData_FreeInverseBetaVsSegmentdT", "", 100, -100, 100, 80, -4, 4);
0161     // MLR
0162     hCSCHaloData_NFlatHaloSegments = ibooker.book1D("CSCHaloData_NFlatHaloSegments", "", 20, 0, 20);
0163     hCSCHaloData_SegmentsInBothEndcaps = ibooker.book1D("CSCHaloData_SegmentsInBothEndcaps", "", 2, 0, 2);
0164     hCSCHaloData_NFlatSegmentsInBothEndcaps = ibooker.book1D("CSCHaloData_NFlatSegmentsInBothEndcaps", "", 20, 0, 20);
0165     // End MLR
0166   } else {
0167     hCSCHaloData_TrackMultiplicity = ibooker.book1D("CSCHaloData_TrackMultiplicity", "", 15, -0.5, 14.5);
0168     hCSCHaloData_TrackMultiplicityMEPlus = ibooker.book1D("CSCHaloData_TrackMultiplicityMEPlus", "", 15, -0.5, 14.5);
0169     hCSCHaloData_TrackMultiplicityMEMinus = ibooker.book1D("CSCHaloData_TrackMultiplicityMEMinus", "", 15, -0.5, 14.5);
0170     hCSCHaloData_InnerMostTrackHitXY =
0171         ibooker.book2D("CSCHaloData_InnerMostTrackHitXY", "", 100, -700, 700, 100, -700, 700);
0172     hCSCHaloData_InnerMostTrackHitR = ibooker.book1D("CSCHaloData_InnerMostTrackHitR", "", 400, -0.5, 799.5);
0173     hCSCHaloData_InnerMostTrackHitRPlusZ =
0174         ibooker.book2D("CSCHaloData_InnerMostTrackHitRPlusZ", "", 400, 400, 1200, 400, -0.5, 799.5);
0175     hCSCHaloData_InnerMostTrackHitRMinusZ =
0176         ibooker.book2D("CSCHaloData_InnerMostTrackHitRMinusZ", "", 400, -1200, -400, 400, -0.5, 799.5);
0177     hCSCHaloData_InnerMostTrackHitiPhi = ibooker.book1D("CSCHaloData_InnerMostTrackHitiPhi", "", 72, 0.5, 72.5);
0178     hCSCHaloData_InnerMostTrackHitPhi = ibooker.book1D("CSCHaloData_InnerMostTrackHitPhi", "", 60, -3.2, 3.2);
0179     hCSCHaloData_L1HaloTriggersMEPlus = ibooker.book1D("CSCHaloData_L1HaloTriggersMEPlus", "", 10, -0.5, 9.5);
0180     hCSCHaloData_L1HaloTriggersMEMinus = ibooker.book1D("CSCHaloData_L1HaloTriggersMEMinus", "", 10, -0.5, 9.5);
0181     hCSCHaloData_L1HaloTriggers = ibooker.book1D("CSCHaloData_L1HaloTriggers", "", 10, -0.5, 9.5);
0182     hCSCHaloData_HLHaloTriggers = ibooker.book1D("CSCHaloData_HLHaloTriggers", "", 2, -0.5, 1.5);
0183     hCSCHaloData_NOutOfTimeTriggersvsL1HaloExists =
0184         ibooker.book2D("CSCHaloData_NOutOfTimeTriggersvsL1HaloExists", "", 20, -0.5, 19.5, 2, -0.5, 1.5);
0185     hCSCHaloData_NOutOfTimeTriggers = ibooker.book1D("CSCHaloData_NOutOfTimeTriggers", "", 20, -0.5, 19.5);
0186     hCSCHaloData_NOutOfTimeHits = ibooker.book1D("CSCHaloData_NOutOfTimeHits", "", 60, -0.5, 59.5);
0187     hCSCHaloData_NTracksSmalldT = ibooker.book1D("CSCHaloData_NTracksSmalldT", "", 15, -0.5, 14.5);
0188     hCSCHaloData_NTracksSmallBeta = ibooker.book1D("CSCHaloData_NTracksSmallBeta", "", 15, -0.5, 14.5);
0189     hCSCHaloData_NTracksSmallBetaAndSmalldT =
0190         ibooker.book1D("CSCHaloData_NTracksSmallBetaAndSmalldT", "", 15, -0.5, 14.5);
0191     hCSCHaloData_NTracksSmalldTvsNHaloTracks =
0192         ibooker.book2D("CSCHaloData_NTracksSmalldTvsNHaloTracks", "", 15, -0.5, 14.5, 15, -0.5, 14.5);
0193     hCSCHaloData_SegmentdT = ibooker.book1D("CSCHaloData_SegmentdT", "", 100, -100, 100);
0194     hCSCHaloData_FreeInverseBeta = ibooker.book1D("CSCHaloData_FreeInverseBeta", "", 80, -4, 4);
0195     hCSCHaloData_FreeInverseBetaVsSegmentdT =
0196         ibooker.book2D("CSCHaloData_FreeInverseBetaVsSegmentdT", "", 100, -100, 100, 80, -4, 4);
0197     // MLR
0198     hCSCHaloData_NFlatHaloSegments = ibooker.book1D("CSCHaloData_NFlatHaloSegments", "", 20, 0, 20);
0199     hCSCHaloData_SegmentsInBothEndcaps = ibooker.book1D("CSCHaloData_SegmentsInBothEndcaps", "", 2, 0, 2);
0200     hCSCHaloData_NFlatSegmentsInBothEndcaps = ibooker.book1D("CSCHaloData_NFlatSegmentsInBothEndcaps", "", 20, 0, 20);
0201     // End MLR
0202   }
0203 
0204   // GlobalHaloData
0205   ibooker.setCurrentFolder(FolderName + "/GlobalHaloData");
0206   if (!StandardDQM) {
0207     hGlobalHaloData_MExCorrection = ibooker.book1D("GlobalHaloData_MExCorrection", "", 200, -200., 200.);
0208     hGlobalHaloData_MEyCorrection = ibooker.book1D("GlobalHaloData_MEyCorrection", "", 200, -200., 200.);
0209     hGlobalHaloData_SumEtCorrection = ibooker.book1D("GlobalHaloData_SumEtCorrection", "", 200, -0.5, 399.5);
0210     hGlobalHaloData_HaloCorrectedMET = ibooker.book1D("GlobalHaloData_HaloCorrectedMET", "", 500, -0.5, 1999.5);
0211     hGlobalHaloData_RawMETMinusHaloCorrectedMET =
0212         ibooker.book1D("GlobalHaloData_RawMETMinusHaloCorrectedMET", "", 250, -500., 500.);
0213     hGlobalHaloData_RawMETOverSumEt = ibooker.book1D("GlobalHaloData_RawMETOverSumEt", "", 100, 0.0, 1.0);
0214     hGlobalHaloData_MatchedHcalPhiWedgeMultiplicity =
0215         ibooker.book1D("GlobalHaloData_MatchedHcalPhiWedgeMultiplicity", "", 15, -0.5, 14.5);
0216     hGlobalHaloData_MatchedHcalPhiWedgeEnergy =
0217         ibooker.book1D("GlobalHaloData_MatchedHcalPhiWedgeEnergy", "", 50, -0.5, 199.5);
0218     hGlobalHaloData_MatchedHcalPhiWedgeConstituents =
0219         ibooker.book1D("GlobalHaloData_MatchedHcalPhiWedgeConstituents", "", 20, -0.5, 19.5);
0220     hGlobalHaloData_MatchedHcalPhiWedgeiPhi =
0221         ibooker.book1D("GlobalHaloData_MatchedHcalPhiWedgeiPhi", "", 1, 0.5, 72.5);
0222     hGlobalHaloData_MatchedHcalPhiWedgeMinTime =
0223         ibooker.book1D("GlobalHaloData_MatchedHcalPhiWedgeMinTime", "", 50, -100.0, 100.0);
0224     hGlobalHaloData_MatchedHcalPhiWedgeMaxTime =
0225         ibooker.book1D("GlobalHaloData_MatchedHcalPhiWedgeMaxTime", "", 50, -100.0, 100.0);
0226     hGlobalHaloData_MatchedHcalPhiWedgeZDirectionConfidence =
0227         ibooker.book1D("GlobalHaloData_MatchedHcalPhiWedgeZDirectionConfidence", "", 120, -1.2, 1.2);
0228     hGlobalHaloData_MatchedEcalPhiWedgeMultiplicity =
0229         ibooker.book1D("GlobalHaloData_MatchedEcalPhiWedgeMultiplicity", "", 15, -0.5, 14.5);
0230     hGlobalHaloData_MatchedEcalPhiWedgeEnergy =
0231         ibooker.book1D("GlobalHaloData_MatchedEcalPhiWedgeEnergy", "", 50, -0.5, 199.5);
0232     hGlobalHaloData_MatchedEcalPhiWedgeConstituents =
0233         ibooker.book1D("GlobalHaloData_MatchedEcalPhiWedgeConstituents", "", 20, -0.5, 19.5);
0234     hGlobalHaloData_MatchedEcalPhiWedgeiPhi =
0235         ibooker.book1D("GlobalHaloData_MatchedEcalPhiWedgeiPhi", "", 360, 0.5, 360.5);
0236     hGlobalHaloData_MatchedEcalPhiWedgeMinTime =
0237         ibooker.book1D("GlobalHaloData_MatchedEcalPhiWedgeMinTime", "", 50, -100.0, 100.0);
0238     hGlobalHaloData_MatchedEcalPhiWedgeMaxTime =
0239         ibooker.book1D("GlobalHaloData_MatchedEcalPhiWedgeMaxTime", "", 50, -100.0, 100.0);
0240     hGlobalHaloData_MatchedEcalPhiWedgeZDirectionConfidence =
0241         ibooker.book1D("GlobalHaloData_MatchedEcalPhiWedgeZDirectionConfidence", "", 120, 1.2, 1.2);
0242   }
0243   // BeamHaloSummary
0244   ibooker.setCurrentFolder(FolderName + "/BeamHaloSummary");
0245 
0246   hBeamHaloSummary_Id = ibooker.book1D("BeamHaloSumamry_Id", "", 11, 0.5, 11.5);
0247   hBeamHaloSummary_Id->setBinLabel(1, "CSC Loose");
0248   hBeamHaloSummary_Id->setBinLabel(2, "CSC Tight");
0249   hBeamHaloSummary_Id->setBinLabel(3, "Ecal Loose");
0250   hBeamHaloSummary_Id->setBinLabel(4, "Ecal Tight");
0251   hBeamHaloSummary_Id->setBinLabel(5, "Hcal Loose");
0252   hBeamHaloSummary_Id->setBinLabel(6, "Hcal Tight");
0253   hBeamHaloSummary_Id->setBinLabel(7, "Global Loose");
0254   hBeamHaloSummary_Id->setBinLabel(8, "Global Tight");
0255   hBeamHaloSummary_Id->setBinLabel(9, "Event Loose");
0256   hBeamHaloSummary_Id->setBinLabel(10, "Event Tight");
0257   hBeamHaloSummary_Id->setBinLabel(11, "Nothing");
0258   if (!StandardDQM) {
0259     hBeamHaloSummary_BXN = ibooker.book2D("BeamHaloSummary_BXN", "", 11, 0.5, 11.5, 4000, -0.5, 3999.5);
0260     hBeamHaloSummary_BXN->setBinLabel(1, "CSC Loose");
0261     hBeamHaloSummary_BXN->setBinLabel(2, "CSC Tight");
0262     hBeamHaloSummary_BXN->setBinLabel(3, "Ecal Loose");
0263     hBeamHaloSummary_BXN->setBinLabel(4, "Ecal Tight");
0264     hBeamHaloSummary_BXN->setBinLabel(5, "Hcal Loose");
0265     hBeamHaloSummary_BXN->setBinLabel(6, "Hcal Tight");
0266     hBeamHaloSummary_BXN->setBinLabel(7, "Global Loose");
0267     hBeamHaloSummary_BXN->setBinLabel(8, "Global Tight");
0268     hBeamHaloSummary_BXN->setBinLabel(9, "Event Loose");
0269     hBeamHaloSummary_BXN->setBinLabel(10, "Event Tight");
0270     hBeamHaloSummary_BXN->setBinLabel(11, "Nothing");
0271   }
0272   // Extra
0273   ibooker.setCurrentFolder(FolderName + "/ExtraHaloData");
0274   if (StandardDQM) {
0275     hExtra_CSCTrackInnerOuterDPhi = ibooker.book1D("Extra_CSCTrackInnerOuterDPhi", "", 30, 0, 3.2);
0276     hExtra_CSCTrackInnerOuterDEta = ibooker.book1D("Extra_CSCTrackInnerOuterDEta", "", 100, 0, 3.0);
0277     hExtra_CSCTrackChi2Ndof = ibooker.book1D("Extra_CSCTrackChi2Ndof", "", 25, 0, 10);
0278     hExtra_CSCTrackNHits = ibooker.book1D("Extra_CSCTrackNHits", "", 75, 0, 75);
0279     hExtra_CSCActivityWithMET = ibooker.book2D("Extra_CSCActivityWithMET", "", 4, 0.5, 4.5, 4, 0.5, 4.5);
0280     hExtra_CSCActivityWithMET->setBinLabel(1, "Track", 1);
0281     hExtra_CSCActivityWithMET->setBinLabel(1, "Track", 2);
0282     hExtra_CSCActivityWithMET->setBinLabel(2, "Segments", 1);
0283     hExtra_CSCActivityWithMET->setBinLabel(2, "Segments", 2);
0284     hExtra_CSCActivityWithMET->setBinLabel(3, "RecHits", 1);
0285     hExtra_CSCActivityWithMET->setBinLabel(3, "RecHits", 2);
0286     hExtra_CSCActivityWithMET->setBinLabel(4, "Nothing", 1);
0287     hExtra_CSCActivityWithMET->setBinLabel(4, "Nothing", 2);
0288     hExtra_InnerMostTrackHitR = ibooker.book1D("Extra_InnerMostTrackHitR", "", 70, 99.5, 799.5);
0289     hExtra_InnerMostTrackHitPhi = ibooker.book1D("Extra_InnerMostTrackHitPhi", "", 60, -3.2, 3.2);
0290   } else {
0291     hExtra_CSCActivityWithMET = ibooker.book2D("Extra_CSCActivityWithMET", "", 4, 0.5, 4.5, 4, 0.5, 4.5);
0292     hExtra_CSCActivityWithMET->setBinLabel(1, "Track", 1);
0293     hExtra_CSCActivityWithMET->setBinLabel(1, "Track", 2);
0294     hExtra_CSCActivityWithMET->setBinLabel(2, "Segments", 1);
0295     hExtra_CSCActivityWithMET->setBinLabel(2, "Segments", 2);
0296     hExtra_CSCActivityWithMET->setBinLabel(3, "RecHits", 1);
0297     hExtra_CSCActivityWithMET->setBinLabel(3, "RecHits", 2);
0298     hExtra_CSCActivityWithMET->setBinLabel(4, "Nothing", 1);
0299     hExtra_CSCActivityWithMET->setBinLabel(4, "Nothing", 2);
0300     hExtra_HcalToF = ibooker.book2D("Extra_HcalToF", "", 83, -41.5, 41.5, 1000, -125., 125.);
0301     hExtra_HcalToF_HaloId = ibooker.book2D("Extra_HcalToF_HaloId", "", 83, -41.5, 41.5, 1000, -125., 125.);
0302     hExtra_EcalToF = ibooker.book2D("Extra_EcalToF", "", 171, -85.5, 85.5, 2000, -225., 225.);
0303     hExtra_EcalToF_HaloId = ibooker.book2D("Extra_EcalToF_HaloId", "", 171, -85.5, 85.5, 2000, -225., 225.);
0304     hExtra_CSCTrackInnerOuterDPhi = ibooker.book1D("Extra_CSCTrackInnerOuterDPhi", "", 30, 0, 3.2);
0305     hExtra_CSCTrackInnerOuterDEta = ibooker.book1D("Extra_CSCTrackInnerOuterDEta", "", 30, 0, 3.2);
0306     hExtra_CSCTrackChi2Ndof = ibooker.book1D("Extra_CSCTrackChi2Ndof", "", 100, 0, 10);
0307     hExtra_CSCTrackNHits = ibooker.book1D("Extra_CSCTrackNHits", "", 75, 0, 75);
0308     hExtra_InnerMostTrackHitXY = ibooker.book2D("Extra_InnerMostTrackHitXY", "", 100, -700, 700, 100, -700, 700);
0309     hExtra_InnerMostTrackHitR = ibooker.book1D("Extra_InnerMostTrackHitR", "", 400, -0.5, 799.5);
0310     hExtra_InnerMostTrackHitRPlusZ =
0311         ibooker.book2D("Extra_InnerMostTrackHitRPlusZ", "", 400, 400, 1200, 400, -0.5, 799.5);
0312     hExtra_InnerMostTrackHitRMinusZ =
0313         ibooker.book2D("Extra_InnerMostTrackHitRMinusZ", "", 400, -1200, -400, 400, -0.5, 799.5);
0314     hExtra_InnerMostTrackHitiPhi = ibooker.book1D("Extra_InnerMostTrackHitiPhi", "", 72, 0.5, 72.5);
0315     hExtra_InnerMostTrackHitPhi = ibooker.book1D("Extra_InnerMostTrackHitPhi", "", 60, -3.2, 3.2);
0316     hExtra_BXN = ibooker.book1D("Extra_BXN", "BXN Occupancy", 4000, 0.5, 4000.5);
0317   }
0318 }
0319 
0320 void BeamHaloAnalyzer::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0321   EventID TheEvent = iEvent.id();
0322   int BXN = iEvent.bunchCrossing();
0323   bool Dump = !TextFileName.empty();
0324   edm::EventNumber_t TheEventNumber = TheEvent.event();
0325   edm::LuminosityBlockNumber_t Lumi = iEvent.luminosityBlock();
0326   edm::RunNumber_t Run = iEvent.run();
0327 
0328   //Get CSC Geometry
0329   const auto& TheCSCGeometry = iSetup.getHandle(cscGeomToken_);
0330   //Note - removed getting calogeometry since it was unused
0331   //Get Stand-alone Muons from Cosmic Muon Reconstruction
0332   edm::Handle<reco::MuonCollection> TheCosmics;
0333   iEvent.getByToken(IT_CosmicStandAloneMuon, TheCosmics);
0334   edm::Handle<reco::MuonTimeExtraMap> TheCSCTimeMap;
0335   iEvent.getByToken(IT_CSCTimeMapToken, TheCSCTimeMap);
0336   bool CSCTrackPlus = false;
0337   bool CSCTrackMinus = false;
0338   int imucount = 0;
0339   if (TheCosmics.isValid()) {
0340     for (reco::MuonCollection::const_iterator iMuon = TheCosmics->begin(); iMuon != TheCosmics->end();
0341          iMuon++, imucount++) {
0342       reco::TrackRef Track = iMuon->outerTrack();
0343       if (!Track)
0344         continue;
0345 
0346       if (!CSCTrackPlus || !CSCTrackMinus) {
0347         if (Track->eta() > 0 || Track->outerPosition().z() > 0 || Track->innerPosition().z() > 0)
0348           CSCTrackPlus = true;
0349         else if (Track->eta() < 0 || Track->outerPosition().z() < 0 || Track->innerPosition().z() < 0)
0350           CSCTrackMinus = true;
0351       }
0352 
0353       float innermost_phi = 0.;
0354       float outermost_phi = 0.;
0355       float innermost_z = 99999.;
0356       float outermost_z = 0.;
0357       float innermost_eta = 0.;
0358       float outermost_eta = 0.;
0359       float innermost_x = 0.;
0360       float innermost_y = 0.;
0361       float innermost_r = 0.;
0362       for (unsigned int j = 0; j < Track->extra()->recHitsSize(); j++) {
0363         auto hit = Track->extra()->recHitRef(j);
0364         DetId TheDetUnitId(hit->geographicalId());
0365         if (TheDetUnitId.det() != DetId::Muon)
0366           continue;
0367         if (TheDetUnitId.subdetId() != MuonSubdetId::CSC)
0368           continue;
0369 
0370         const GeomDetUnit* TheUnit = TheCSCGeometry->idToDetUnit(TheDetUnitId);
0371         LocalPoint TheLocalPosition = hit->localPosition();
0372         const BoundPlane& TheSurface = TheUnit->surface();
0373         const GlobalPoint TheGlobalPosition = TheSurface.toGlobal(TheLocalPosition);
0374 
0375         float z = TheGlobalPosition.z();
0376         if (TMath::Abs(z) < innermost_z) {
0377           innermost_phi = TheGlobalPosition.phi();
0378           innermost_eta = TheGlobalPosition.eta();
0379           innermost_z = TheGlobalPosition.z();
0380           innermost_x = TheGlobalPosition.x();
0381           innermost_y = TheGlobalPosition.y();
0382           innermost_r = TMath::Sqrt(innermost_x * innermost_x + innermost_y * innermost_y);
0383         }
0384         if (TMath::Abs(z) > outermost_z) {
0385           outermost_phi = TheGlobalPosition.phi();
0386           outermost_eta = TheGlobalPosition.eta();
0387           outermost_z = TheGlobalPosition.z();
0388         }
0389       }
0390       float dphi = TMath::Abs(outermost_phi - innermost_phi);
0391       float deta = TMath::Abs(outermost_eta - innermost_eta);
0392       hExtra_CSCTrackInnerOuterDPhi->Fill(dphi);
0393       hExtra_CSCTrackInnerOuterDEta->Fill(deta);
0394       hExtra_CSCTrackChi2Ndof->Fill(Track->normalizedChi2());
0395       hExtra_CSCTrackNHits->Fill(Track->numberOfValidHits());
0396       hExtra_InnerMostTrackHitR->Fill(innermost_r);
0397       hExtra_InnerMostTrackHitPhi->Fill(innermost_phi);
0398       if (!StandardDQM) {
0399         hExtra_InnerMostTrackHitXY->Fill(innermost_x, innermost_y);
0400         hExtra_InnerMostTrackHitiPhi->Fill(Phi_To_iPhi(innermost_phi));
0401         if (innermost_z > 0)
0402           hExtra_InnerMostTrackHitRPlusZ->Fill(innermost_z, innermost_r);
0403         else
0404           hExtra_InnerMostTrackHitRMinusZ->Fill(innermost_z, innermost_r);
0405       }
0406 
0407       std::vector<const CSCSegment*> MatchedSegments = TheMatcher->matchCSC(*Track, iEvent);
0408       // Find the inner and outer segments separately in case they don't agree completely with recHits
0409       // Plan for the possibility segments in both endcaps
0410       float InnerSegmentTime[2] = {0, 0};
0411       float OuterSegmentTime[2] = {0, 0};
0412       float innermost_seg_z[2] = {1500, 1500};
0413       float outermost_seg_z[2] = {0, 0};
0414       for (std::vector<const CSCSegment*>::const_iterator segment = MatchedSegments.begin();
0415            segment != MatchedSegments.end();
0416            ++segment) {
0417         CSCDetId TheCSCDetId((*segment)->cscDetId());
0418         const CSCChamber* TheCSCChamber = TheCSCGeometry->chamber(TheCSCDetId);
0419         LocalPoint TheLocalPosition = (*segment)->localPosition();
0420         const GlobalPoint TheGlobalPosition = TheCSCChamber->toGlobal(TheLocalPosition);
0421         float z = TheGlobalPosition.z();
0422         int TheEndcap = TheCSCDetId.endcap();
0423         if (TMath::Abs(z) < innermost_seg_z[TheEndcap - 1]) {
0424           innermost_seg_z[TheEndcap - 1] = TMath::Abs(z);
0425           InnerSegmentTime[TheEndcap - 1] = (*segment)->time();
0426         }
0427         if (TMath::Abs(z) > outermost_seg_z[TheEndcap - 1]) {
0428           outermost_seg_z[TheEndcap - 1] = TMath::Abs(z);
0429           OuterSegmentTime[TheEndcap - 1] = (*segment)->time();
0430         }
0431       }
0432 
0433       float dT_Segment = 0;                         // default safe value, looks like collision muon
0434       if (innermost_seg_z[0] < outermost_seg_z[0])  // two segments in ME+
0435         dT_Segment = OuterSegmentTime[0] - InnerSegmentTime[0];
0436       if (innermost_seg_z[1] < outermost_seg_z[1])  // two segments in ME-
0437       {
0438         // replace the measurement if there weren't segments in ME+ or
0439         // if the track in ME- has timing more consistent with an incoming particle
0440         if (dT_Segment == 0.0 || OuterSegmentTime[1] - InnerSegmentTime[1] < dT_Segment)
0441           dT_Segment = OuterSegmentTime[1] - InnerSegmentTime[1];
0442       }
0443       hCSCHaloData_SegmentdT->Fill(dT_Segment);
0444 
0445       // Analyze the MuonTimeExtra information
0446       reco::MuonRef muonR(TheCosmics, imucount);
0447       if (TheCSCTimeMap.isValid()) {
0448         const reco::MuonTimeExtraMap& timeMapCSC = *TheCSCTimeMap;
0449         reco::MuonTimeExtra timecsc = timeMapCSC[muonR];
0450         float freeInverseBeta = timecsc.freeInverseBeta();
0451         hCSCHaloData_FreeInverseBeta->Fill(freeInverseBeta);
0452         hCSCHaloData_FreeInverseBetaVsSegmentdT->Fill(dT_Segment, freeInverseBeta);
0453       }
0454     }
0455   }
0456 
0457   //Get CSC Segments
0458   edm::Handle<CSCSegmentCollection> TheCSCSegments;
0459   iEvent.getByToken(IT_CSCSegment, TheCSCSegments);
0460 
0461   // Group segments according to endcaps
0462   std::vector<CSCSegment> vCSCSegments_Plus;
0463   std::vector<CSCSegment> vCSCSegments_Minus;
0464 
0465   bool CSCSegmentPlus = false;
0466   bool CSCSegmentMinus = false;
0467   if (TheCSCSegments.isValid()) {
0468     for (CSCSegmentCollection::const_iterator iSegment = TheCSCSegments->begin(); iSegment != TheCSCSegments->end();
0469          iSegment++) {
0470       const std::vector<CSCRecHit2D> vCSCRecHits = iSegment->specificRecHits();
0471       CSCDetId iDetId = (CSCDetId)(*iSegment).cscDetId();
0472 
0473       if (iDetId.endcap() == 1)
0474         vCSCSegments_Plus.push_back(*iSegment);
0475       else
0476         vCSCSegments_Minus.push_back(*iSegment);
0477     }
0478   }
0479 
0480   // Are there segments on the plus/minus side?
0481   if (!vCSCSegments_Plus.empty())
0482     CSCSegmentPlus = true;
0483   if (!vCSCSegments_Minus.empty())
0484     CSCSegmentMinus = true;
0485 
0486   //Get CSC RecHits
0487   Handle<CSCRecHit2DCollection> TheCSCRecHits;
0488   iEvent.getByToken(IT_CSCRecHit, TheCSCRecHits);
0489   bool CSCRecHitPlus = false;
0490   bool CSCRecHitMinus = false;
0491   if (TheCSCRecHits.isValid()) {
0492     for (CSCRecHit2DCollection::const_iterator iCSCRecHit = TheCSCRecHits->begin(); iCSCRecHit != TheCSCRecHits->end();
0493          iCSCRecHit++) {
0494       DetId TheDetUnitId(iCSCRecHit->geographicalId());
0495       const GeomDetUnit* TheUnit = (*TheCSCGeometry).idToDetUnit(TheDetUnitId);
0496       LocalPoint TheLocalPosition = iCSCRecHit->localPosition();
0497       const BoundPlane& TheSurface = TheUnit->surface();
0498       GlobalPoint TheGlobalPosition = TheSurface.toGlobal(TheLocalPosition);
0499 
0500       //Are there hits on the plus/minus side?
0501       if (TheGlobalPosition.z() > 0)
0502         CSCRecHitPlus = true;
0503       else
0504         CSCRecHitMinus = true;
0505     }
0506   }
0507 
0508   //Get  EB RecHits
0509   edm::Handle<EBRecHitCollection> TheEBRecHits;
0510   iEvent.getByToken(IT_EBRecHit, TheEBRecHits);
0511   if (TheEBRecHits.isValid()) {
0512     for (EBRecHitCollection::const_iterator iEBRecHit = TheEBRecHits->begin(); iEBRecHit != TheEBRecHits->end();
0513          iEBRecHit++) {
0514       if (iEBRecHit->energy() < 0.5)
0515         continue;
0516       DetId id = DetId(iEBRecHit->id());
0517       EBDetId EcalId(id.rawId());
0518       int ieta = EcalId.ieta();
0519       if (!StandardDQM)
0520         hExtra_EcalToF->Fill(ieta, iEBRecHit->time());
0521     }
0522   }
0523 
0524   //Get HB/HE RecHits
0525   edm::Handle<HBHERecHitCollection> TheHBHERecHits;
0526   iEvent.getByToken(IT_HBHERecHit, TheHBHERecHits);
0527   if (TheHBHERecHits.isValid()) {
0528     for (HBHERecHitCollection::const_iterator iHBHERecHit = TheHBHERecHits->begin();
0529          iHBHERecHit != TheHBHERecHits->end();
0530          iHBHERecHit++) {
0531       if (iHBHERecHit->energy() < 1.)
0532         continue;
0533       HcalDetId id = HcalDetId(iHBHERecHit->id());
0534       if (!StandardDQM)
0535         hExtra_HcalToF->Fill(id.ieta(), iHBHERecHit->time());
0536     }
0537   }
0538 
0539   //Get MET
0540   edm::Handle<reco::CaloMETCollection> TheCaloMET;
0541   iEvent.getByToken(IT_met, TheCaloMET);
0542 
0543   //Get CSCHaloData
0544   edm::Handle<reco::CSCHaloData> TheCSCDataHandle;
0545   iEvent.getByToken(IT_CSCHaloData, TheCSCDataHandle);
0546   int TheHaloOrigin = 0;
0547   if (TheCSCDataHandle.isValid()) {
0548     const CSCHaloData CSCData = (*TheCSCDataHandle.product());
0549     if (CSCData.NumberOfOutOfTimeTriggers(HaloData::plus) && !CSCData.NumberOfOutOfTimeTriggers(HaloData::minus))
0550       TheHaloOrigin = 1;
0551     else if (CSCData.NumberOfOutOfTimeTriggers(HaloData::minus) && !CSCData.NumberOfOutOfTimeTriggers(HaloData::plus))
0552       TheHaloOrigin = -1;
0553 
0554     for (std::vector<GlobalPoint>::const_iterator i = CSCData.GetCSCTrackImpactPositions().begin();
0555          i != CSCData.GetCSCTrackImpactPositions().end();
0556          i++) {
0557       float r = TMath::Sqrt(i->x() * i->x() + i->y() * i->y());
0558       if (!StandardDQM) {
0559         hCSCHaloData_InnerMostTrackHitXY->Fill(i->x(), i->y());
0560         hCSCHaloData_InnerMostTrackHitiPhi->Fill(Phi_To_iPhi(i->phi()));
0561         if (i->z() > 0)
0562           hCSCHaloData_InnerMostTrackHitRPlusZ->Fill(i->z(), r);
0563         else
0564           hCSCHaloData_InnerMostTrackHitRMinusZ->Fill(i->z(), r);
0565       }
0566       hCSCHaloData_InnerMostTrackHitR->Fill(r);
0567       hCSCHaloData_InnerMostTrackHitPhi->Fill(i->phi());
0568     }
0569     hCSCHaloData_L1HaloTriggersMEPlus->Fill(CSCData.NumberOfHaloTriggers(HaloData::plus));
0570     hCSCHaloData_L1HaloTriggersMEMinus->Fill(CSCData.NumberOfHaloTriggers(HaloData::minus));
0571     hCSCHaloData_L1HaloTriggers->Fill(CSCData.NumberOfHaloTriggers(HaloData::both));
0572     hCSCHaloData_HLHaloTriggers->Fill(CSCData.CSCHaloHLTAccept());
0573     hCSCHaloData_TrackMultiplicityMEPlus->Fill(CSCData.NumberOfHaloTracks(HaloData::plus));
0574     hCSCHaloData_TrackMultiplicityMEMinus->Fill(CSCData.NumberOfHaloTracks(HaloData::minus));
0575     hCSCHaloData_TrackMultiplicity->Fill(CSCData.GetTracks().size());
0576     hCSCHaloData_NOutOfTimeTriggersMEPlus->Fill(CSCData.NOutOfTimeTriggers(HaloData::plus));
0577     hCSCHaloData_NOutOfTimeTriggersMEMinus->Fill(CSCData.NOutOfTimeTriggers(HaloData::minus));
0578     hCSCHaloData_NOutOfTimeTriggers->Fill(CSCData.NOutOfTimeTriggers(HaloData::both));
0579     hCSCHaloData_NOutOfTimeHits->Fill(CSCData.NOutOfTimeHits());
0580     hCSCHaloData_NOutOfTimeTriggersvsL1HaloExists->Fill(CSCData.NOutOfTimeTriggers(HaloData::both),
0581                                                         CSCData.NumberOfHaloTriggers(HaloData::both) > 0);
0582     hCSCHaloData_NTracksSmalldT->Fill(CSCData.NTracksSmalldT());
0583     hCSCHaloData_NTracksSmallBeta->Fill(CSCData.NTracksSmallBeta());
0584     hCSCHaloData_NTracksSmallBetaAndSmalldT->Fill(CSCData.NTracksSmallBetaAndSmalldT());
0585     hCSCHaloData_NTracksSmalldTvsNHaloTracks->Fill(CSCData.GetTracks().size(), CSCData.NTracksSmalldT());
0586     // MLR
0587     hCSCHaloData_NFlatHaloSegments->Fill(CSCData.NFlatHaloSegments());
0588     hCSCHaloData_SegmentsInBothEndcaps->Fill(CSCData.GetSegmentsInBothEndcaps());
0589     if (CSCData.GetSegmentsInBothEndcaps())
0590       hCSCHaloData_NFlatSegmentsInBothEndcaps->Fill(CSCData.NFlatHaloSegments());
0591     // End MLR
0592   }
0593 
0594   //Get EcalHaloData
0595   edm::Handle<reco::EcalHaloData> TheEcalHaloData;
0596   iEvent.getByToken(IT_EcalHaloData, TheEcalHaloData);
0597   if (TheEcalHaloData.isValid()) {
0598     const EcalHaloData EcalData = (*TheEcalHaloData.product());
0599     std::vector<PhiWedge> EcalWedges = EcalData.GetPhiWedges();
0600     for (std::vector<PhiWedge>::const_iterator iWedge = EcalWedges.begin(); iWedge != EcalWedges.end(); iWedge++) {
0601       if (!StandardDQM) {
0602         hEcalHaloData_PhiWedgeEnergy->Fill(iWedge->Energy());
0603         hEcalHaloData_PhiWedgeMinTime->Fill(iWedge->MinTime());
0604         hEcalHaloData_PhiWedgeMaxTime->Fill(iWedge->MaxTime());
0605         hEcalHaloData_PhiWedgeMinVsMaxTime->Fill(iWedge->MinTime(), iWedge->MaxTime());
0606         hEcalHaloData_PhiWedgePlusZDirectionConfidence->Fill(iWedge->PlusZDirectionConfidence());
0607         hEcalHaloData_PhiWedgeiPhi->Fill(iWedge->iPhi());
0608       }
0609       hEcalHaloData_PhiWedgeZDirectionConfidence->Fill(iWedge->ZDirectionConfidence());
0610       hEcalHaloData_PhiWedgeConstituents->Fill(iWedge->NumberOfConstituents());
0611     }
0612 
0613     hEcalHaloData_PhiWedgeMultiplicity->Fill(EcalWedges.size());
0614 
0615     edm::ValueMap<float> vm_Angle = EcalData.GetShowerShapesAngle();
0616     edm::ValueMap<float> vm_Roundness = EcalData.GetShowerShapesRoundness();
0617     //Access selected SuperClusters
0618     for (unsigned int n = 0; n < EcalData.GetSuperClusters().size(); n++) {
0619       edm::Ref<SuperClusterCollection> cluster(EcalData.GetSuperClusters()[n]);
0620       float angle = vm_Angle[cluster];
0621       float roundness = vm_Roundness[cluster];
0622       hEcalHaloData_SuperClusterShowerShapes->Fill(angle, roundness);
0623       hEcalHaloData_SuperClusterNHits->Fill(cluster->size());
0624       hEcalHaloData_SuperClusterEnergy->Fill(cluster->energy());
0625 
0626       if (!StandardDQM) {
0627         hEcalHaloData_SuperClusterPhiVsEta->Fill(cluster->eta(), cluster->phi());
0628       }
0629     }
0630   }
0631 
0632   //Get HcalHaloData
0633   edm::Handle<reco::HcalHaloData> TheHcalHaloData;
0634   iEvent.getByToken(IT_HcalHaloData, TheHcalHaloData);
0635   if (TheHcalHaloData.isValid()) {
0636     const HcalHaloData HcalData = (*TheHcalHaloData.product());
0637     std::vector<PhiWedge> HcalWedges = HcalData.GetPhiWedges();
0638     hHcalHaloData_PhiWedgeMultiplicity->Fill(HcalWedges.size());
0639     for (std::vector<PhiWedge>::const_iterator iWedge = HcalWedges.begin(); iWedge != HcalWedges.end(); iWedge++) {
0640       if (!StandardDQM) {
0641         hHcalHaloData_PhiWedgeEnergy->Fill(iWedge->Energy());
0642         hHcalHaloData_PhiWedgeMinTime->Fill(iWedge->MinTime());
0643         hHcalHaloData_PhiWedgeMaxTime->Fill(iWedge->MaxTime());
0644         hHcalHaloData_PhiWedgePlusZDirectionConfidence->Fill(iWedge->PlusZDirectionConfidence());
0645         hHcalHaloData_PhiWedgeMinVsMaxTime->Fill(iWedge->MinTime(), iWedge->MaxTime());
0646         hHcalHaloData_PhiWedgeiPhi->Fill(iWedge->iPhi());
0647       }
0648 
0649       hHcalHaloData_PhiWedgeConstituents->Fill(iWedge->NumberOfConstituents());
0650       hHcalHaloData_PhiWedgeZDirectionConfidence->Fill(iWedge->ZDirectionConfidence());
0651     }
0652   }
0653 
0654   if (!StandardDQM) {
0655     //Get GlobalHaloData
0656     edm::Handle<reco::GlobalHaloData> TheGlobalHaloData;
0657     iEvent.getByToken(IT_GlobalHaloData, TheGlobalHaloData);
0658     if (TheGlobalHaloData.isValid()) {
0659       const GlobalHaloData GlobalData = (*TheGlobalHaloData.product());
0660       if (TheCaloMET.isValid()) {
0661         // Get Raw Uncorrected CaloMET
0662         const CaloMETCollection* calometcol = TheCaloMET.product();
0663         const CaloMET* RawMET = &(calometcol->front());
0664 
0665         // Get BeamHalo Corrected CaloMET
0666         const CaloMET CorrectedMET = GlobalData.GetCorrectedCaloMET(*RawMET);
0667         hGlobalHaloData_MExCorrection->Fill(GlobalData.DeltaMEx());
0668         hGlobalHaloData_MEyCorrection->Fill(GlobalData.DeltaMEy());
0669         hGlobalHaloData_HaloCorrectedMET->Fill(CorrectedMET.pt());
0670         hGlobalHaloData_RawMETMinusHaloCorrectedMET->Fill(RawMET->pt() - CorrectedMET.pt());
0671         if (RawMET->sumEt())
0672           hGlobalHaloData_RawMETOverSumEt->Fill(RawMET->pt() / RawMET->sumEt());
0673       }
0674 
0675       // Get Matched Hcal Phi Wedges
0676       std::vector<PhiWedge> HcalWedges = GlobalData.GetMatchedHcalPhiWedges();
0677       hGlobalHaloData_MatchedHcalPhiWedgeMultiplicity->Fill(HcalWedges.size());
0678       // Loop over Matched Hcal Phi Wedges
0679       for (std::vector<PhiWedge>::const_iterator iWedge = HcalWedges.begin(); iWedge != HcalWedges.end(); iWedge++) {
0680         hGlobalHaloData_MatchedHcalPhiWedgeEnergy->Fill(iWedge->Energy());
0681         hGlobalHaloData_MatchedHcalPhiWedgeConstituents->Fill(iWedge->NumberOfConstituents());
0682         hGlobalHaloData_MatchedHcalPhiWedgeiPhi->Fill(iWedge->iPhi());
0683         hGlobalHaloData_MatchedHcalPhiWedgeMinTime->Fill(iWedge->MinTime());
0684         hGlobalHaloData_MatchedHcalPhiWedgeMaxTime->Fill(iWedge->MaxTime());
0685         hGlobalHaloData_MatchedHcalPhiWedgeZDirectionConfidence->Fill(iWedge->ZDirectionConfidence());
0686         if (TheHBHERecHits.isValid()) {
0687           for (HBHERecHitCollection::const_iterator iHBHERecHit = TheHBHERecHits->begin();
0688                iHBHERecHit != TheHBHERecHits->end();
0689                iHBHERecHit++) {
0690             HcalDetId id = HcalDetId(iHBHERecHit->id());
0691             int iphi = id.iphi();
0692             if (iphi != iWedge->iPhi())
0693               continue;
0694             if (iHBHERecHit->energy() < 1.0)
0695               continue;  // Otherwise there are thousands of hits per event (even with negative energies)
0696 
0697             float time = iHBHERecHit->time();
0698             int ieta = id.ieta();
0699             hExtra_HcalToF_HaloId->Fill(ieta, time);
0700           }
0701         }
0702       }
0703 
0704       // Get Matched Hcal Phi Wedges
0705       std::vector<PhiWedge> EcalWedges = GlobalData.GetMatchedEcalPhiWedges();
0706       hGlobalHaloData_MatchedEcalPhiWedgeMultiplicity->Fill(EcalWedges.size());
0707       for (std::vector<PhiWedge>::const_iterator iWedge = EcalWedges.begin(); iWedge != EcalWedges.end(); iWedge++) {
0708         hGlobalHaloData_MatchedEcalPhiWedgeEnergy->Fill(iWedge->Energy());
0709         hGlobalHaloData_MatchedEcalPhiWedgeConstituents->Fill(iWedge->NumberOfConstituents());
0710         hGlobalHaloData_MatchedEcalPhiWedgeiPhi->Fill(iWedge->iPhi());
0711         hGlobalHaloData_MatchedEcalPhiWedgeMinTime->Fill(iWedge->MinTime());
0712         hGlobalHaloData_MatchedEcalPhiWedgeMaxTime->Fill(iWedge->MaxTime());
0713         hGlobalHaloData_MatchedEcalPhiWedgeZDirectionConfidence->Fill(iWedge->ZDirectionConfidence());
0714         if (TheEBRecHits.isValid()) {
0715           for (EBRecHitCollection::const_iterator iEBRecHit = TheEBRecHits->begin(); iEBRecHit != TheEBRecHits->end();
0716                iEBRecHit++) {
0717             if (iEBRecHit->energy() < 0.5)
0718               continue;
0719             DetId id = DetId(iEBRecHit->id());
0720             EBDetId EcalId(id.rawId());
0721             int iPhi = EcalId.iphi();
0722             iPhi = (iPhi - 1) / 5 + 1;
0723             if (iPhi != iWedge->iPhi())
0724               continue;
0725             hExtra_EcalToF_HaloId->Fill(EcalId.ieta(), iEBRecHit->time());
0726           }
0727         }
0728       }
0729     }
0730   }
0731 
0732   // Get BeamHaloSummary
0733   edm::Handle<BeamHaloSummary> TheBeamHaloSummary;
0734   iEvent.getByToken(IT_BeamHaloSummary, TheBeamHaloSummary);
0735   if (TheBeamHaloSummary.isValid()) {
0736     const BeamHaloSummary TheSummary = (*TheBeamHaloSummary.product());
0737     if (TheSummary.CSCLooseHaloId()) {
0738       hBeamHaloSummary_Id->Fill(1);
0739       if (!StandardDQM)
0740         hBeamHaloSummary_BXN->Fill(1, BXN);
0741       if (Dump)
0742         *out << std::setw(15) << "CSCLoose" << std::setw(15) << Run << std::setw(15) << Lumi << std::setw(15)
0743              << TheEventNumber << std::endl;
0744     }
0745     if (TheSummary.CSCTightHaloId()) {
0746       hBeamHaloSummary_Id->Fill(2);
0747       if (!StandardDQM)
0748         hBeamHaloSummary_BXN->Fill(2, BXN);
0749     }
0750     if (TheSummary.EcalLooseHaloId()) {
0751       hBeamHaloSummary_Id->Fill(3);
0752       if (!StandardDQM)
0753         hBeamHaloSummary_BXN->Fill(3, BXN);
0754       if (Dump)
0755         *out << std::setw(15) << "EcalLoose" << std::setw(15) << Run << std::setw(15) << Lumi << std::setw(15)
0756              << TheEventNumber << std::endl;
0757     }
0758     if (TheSummary.EcalTightHaloId()) {
0759       hBeamHaloSummary_Id->Fill(4);
0760       if (!StandardDQM)
0761         hBeamHaloSummary_BXN->Fill(4, BXN);
0762     }
0763     if (TheSummary.HcalLooseHaloId()) {
0764       hBeamHaloSummary_Id->Fill(5);
0765       if (!StandardDQM)
0766         hBeamHaloSummary_BXN->Fill(5, BXN);
0767       if (Dump)
0768         *out << std::setw(15) << "HcalLoose" << std::setw(15) << Run << std::setw(15) << Lumi << std::setw(15)
0769              << TheEventNumber << std::endl;
0770     }
0771     if (TheSummary.HcalTightHaloId()) {
0772       hBeamHaloSummary_Id->Fill(6);
0773       if (!StandardDQM)
0774         hBeamHaloSummary_BXN->Fill(6, BXN);
0775     }
0776     if (TheSummary.GlobalLooseHaloId()) {
0777       hBeamHaloSummary_Id->Fill(7);
0778       if (!StandardDQM)
0779         hBeamHaloSummary_BXN->Fill(7, BXN);
0780       if (Dump)
0781         *out << std::setw(15) << "GlobalLoose" << std::setw(15) << Run << std::setw(15) << Lumi << std::setw(15)
0782              << TheEventNumber << std::endl;
0783     }
0784     if (TheSummary.GlobalTightHaloId()) {
0785       hBeamHaloSummary_Id->Fill(8);
0786       if (!StandardDQM)
0787         hBeamHaloSummary_BXN->Fill(8, BXN);
0788     }
0789     if (TheSummary.LooseId()) {
0790       hBeamHaloSummary_Id->Fill(9);
0791       if (!StandardDQM)
0792         hBeamHaloSummary_BXN->Fill(9, BXN);
0793     }
0794     if (TheSummary.TightId()) {
0795       hBeamHaloSummary_Id->Fill(10);
0796       if (!StandardDQM)
0797         hBeamHaloSummary_BXN->Fill(10, BXN);
0798     }
0799     if (!TheSummary.EcalLooseHaloId() && !TheSummary.HcalLooseHaloId() && !TheSummary.CSCLooseHaloId() &&
0800         !TheSummary.GlobalLooseHaloId()) {
0801       hBeamHaloSummary_Id->Fill(11);
0802       if (!StandardDQM)
0803         hBeamHaloSummary_BXN->Fill(11, BXN);
0804     }
0805   }
0806 
0807   if (TheCaloMET.isValid()) {
0808     const CaloMETCollection* calometcol = TheCaloMET.product();
0809     const CaloMET* calomet = &(calometcol->front());
0810 
0811     if (calomet->pt() > DumpMET)
0812       if (Dump)
0813         *out << std::setw(15) << "HighMET" << std::setw(15) << Run << std::setw(15) << Lumi << std::setw(15)
0814              << TheEventNumber << std::endl;
0815 
0816     //Fill CSC Activity Plot
0817     if (calomet->pt() > 15.0) {
0818       if (TheHaloOrigin > 0) {
0819         if (CSCTrackPlus && CSCTrackMinus)
0820           hExtra_CSCActivityWithMET->Fill(1, 1);
0821         else if (CSCTrackPlus && CSCSegmentMinus)
0822           hExtra_CSCActivityWithMET->Fill(1, 2);
0823         else if (CSCTrackPlus && CSCRecHitMinus)
0824           hExtra_CSCActivityWithMET->Fill(1, 3);
0825         else if (CSCTrackPlus)
0826           hExtra_CSCActivityWithMET->Fill(1, 4);
0827         else if (CSCSegmentPlus && CSCTrackMinus)
0828           hExtra_CSCActivityWithMET->Fill(2, 1);
0829         else if (CSCSegmentPlus && CSCSegmentMinus)
0830           hExtra_CSCActivityWithMET->Fill(2, 2);
0831         else if (CSCSegmentPlus && CSCRecHitMinus)
0832           hExtra_CSCActivityWithMET->Fill(2, 3);
0833         else if (CSCSegmentPlus)
0834           hExtra_CSCActivityWithMET->Fill(2, 4);
0835         else if (CSCRecHitPlus && CSCTrackMinus)
0836           hExtra_CSCActivityWithMET->Fill(3, 1);
0837         else if (CSCRecHitPlus && CSCSegmentMinus)
0838           hExtra_CSCActivityWithMET->Fill(3, 2);
0839         else if (CSCRecHitPlus && CSCRecHitMinus)
0840           hExtra_CSCActivityWithMET->Fill(3, 3);
0841         else if (CSCRecHitPlus)
0842           hExtra_CSCActivityWithMET->Fill(3, 4);
0843         else
0844           hExtra_CSCActivityWithMET->Fill(4, 4);
0845       } else if (TheHaloOrigin < 0) {
0846         if (CSCTrackMinus && CSCTrackPlus)
0847           hExtra_CSCActivityWithMET->Fill(1, 1);
0848         else if (CSCTrackMinus && CSCSegmentPlus)
0849           hExtra_CSCActivityWithMET->Fill(1, 2);
0850         else if (CSCTrackMinus && CSCRecHitPlus)
0851           hExtra_CSCActivityWithMET->Fill(1, 3);
0852         else if (CSCTrackMinus)
0853           hExtra_CSCActivityWithMET->Fill(1, 4);
0854         else if (CSCSegmentMinus && CSCTrackPlus)
0855           hExtra_CSCActivityWithMET->Fill(2, 1);
0856         else if (CSCSegmentMinus && CSCSegmentPlus)
0857           hExtra_CSCActivityWithMET->Fill(2, 2);
0858         else if (CSCSegmentMinus && CSCRecHitPlus)
0859           hExtra_CSCActivityWithMET->Fill(2, 3);
0860         else if (CSCSegmentMinus)
0861           hExtra_CSCActivityWithMET->Fill(2, 4);
0862         else if (CSCRecHitMinus && CSCTrackPlus)
0863           hExtra_CSCActivityWithMET->Fill(3, 1);
0864         else if (CSCRecHitMinus && CSCSegmentPlus)
0865           hExtra_CSCActivityWithMET->Fill(3, 2);
0866         else if (CSCRecHitMinus && CSCRecHitPlus)
0867           hExtra_CSCActivityWithMET->Fill(3, 3);
0868         else if (CSCRecHitMinus)
0869           hExtra_CSCActivityWithMET->Fill(3, 4);
0870         else
0871           hExtra_CSCActivityWithMET->Fill(4, 4);
0872       }
0873     }
0874   }
0875 }
0876 
0877 BeamHaloAnalyzer::~BeamHaloAnalyzer() {}
0878 
0879 //DEFINE_FWK_MODULE(CMSEventAnalyzer);