File indexing completed on 2024-04-06 12:09:20
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020 #include <memory>
0021
0022
0023
0024 #include "DataFormats/Candidate/interface/Candidate.h"
0025 #include "DataFormats/Common/interface/Handle.h"
0026 #include "DataFormats/EgammaCandidates/interface/GsfElectron.h"
0027 #include "DataFormats/EgammaCandidates/interface/GsfElectronFwd.h"
0028 #include "DataFormats/GsfTrackReco/interface/GsfTrack.h"
0029 #include "DataFormats/Math/interface/LorentzVector.h"
0030 #include "FWCore/Framework/interface/Event.h"
0031 #include "FWCore/Framework/interface/EventSetup.h"
0032 #include "FWCore/Framework/interface/MakerMacros.h"
0033 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0034 #include "TLorentzVector.h"
0035 #include "TMath.h"
0036 #include <cmath>
0037 #include <iostream>
0038 #include <string>
0039
0040 #include "DQMServices/Core/interface/DQMEDAnalyzer.h"
0041 #include "DQMServices/Core/interface/DQMStore.h"
0042
0043 class EcalZmassTask : public DQMEDAnalyzer {
0044 public:
0045 explicit EcalZmassTask(const edm::ParameterSet &);
0046 ~EcalZmassTask() override {}
0047
0048 private:
0049 void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override;
0050 void analyze(const edm::Event &, const edm::EventSetup &) override;
0051
0052 const edm::EDGetTokenT<reco::GsfElectronCollection> electronCollectionToken_;
0053 const edm::EDGetTokenT<reco::GsfTrackCollection> trackCollectionToken_;
0054
0055 const std::string prefixME_;
0056
0057 MonitorElement *h_ee_invMass_EB;
0058 MonitorElement *h_ee_invMass_EE;
0059 MonitorElement *h_ee_invMass_BB;
0060 MonitorElement *h_ee_invMass;
0061 MonitorElement *h_e1_et;
0062 MonitorElement *h_e2_et;
0063 MonitorElement *h_e1_eta;
0064 MonitorElement *h_e2_eta;
0065 MonitorElement *h_e1_phi;
0066 MonitorElement *h_e2_phi;
0067 MonitorElement *h_95_ee_invMass_EB;
0068 MonitorElement *h_95_ee_invMass_EE;
0069 MonitorElement *h_95_ee_invMass_BB;
0070 };
0071
0072 EcalZmassTask::EcalZmassTask(const edm::ParameterSet ¶meters)
0073 : electronCollectionToken_(
0074 consumes<reco::GsfElectronCollection>(parameters.getParameter<edm::InputTag>("electronCollection"))),
0075 trackCollectionToken_(
0076 consumes<reco::GsfTrackCollection>(parameters.getParameter<edm::InputTag>("trackCollection"))),
0077 prefixME_(parameters.getUntrackedParameter<std::string>("prefixME", "")) {}
0078
0079
0080 void EcalZmassTask::analyze(const edm::Event &iEvent, const edm::EventSetup &iSetup) {
0081 using namespace edm;
0082 Handle<reco::GsfElectronCollection> electronCollection;
0083 iEvent.getByToken(electronCollectionToken_, electronCollection);
0084 if (!electronCollection.isValid())
0085 return;
0086
0087
0088 Handle<reco::GsfTrackCollection> gsftracks_h;
0089 iEvent.getByToken(trackCollectionToken_, gsftracks_h);
0090
0091 bool isIsolatedBarrel;
0092 bool isIDBarrel;
0093 bool isConvertedBarrel;
0094 bool isIsolatedEndcap;
0095 bool isIDEndcap;
0096 bool isConvertedEndcap;
0097
0098 int elIsAccepted = 0;
0099 int elIsAcceptedEB = 0;
0100 int elIsAcceptedEE = 0;
0101
0102 std::vector<TLorentzVector> LV;
0103
0104 for (reco::GsfElectronCollection::const_iterator recoElectron = electronCollection->begin();
0105 recoElectron != electronCollection->end();
0106 recoElectron++) {
0107 if (recoElectron->et() <= 25)
0108 continue;
0109
0110
0111 double IsoTrk = (recoElectron->dr03TkSumPt() / recoElectron->et());
0112 double IsoEcal = (recoElectron->dr03EcalRecHitSumEt() / recoElectron->et());
0113 double IsoHcal = (recoElectron->dr03HcalTowerSumEt() / recoElectron->et());
0114 double HE = (recoElectron->hcalOverEcal());
0115
0116
0117
0118 float DeltaPhiTkClu = recoElectron->deltaPhiSuperClusterTrackAtVtx();
0119 float DeltaEtaTkClu = recoElectron->deltaEtaSuperClusterTrackAtVtx();
0120 float sigmaIeIe = recoElectron->sigmaIetaIeta();
0121
0122
0123
0124 float Dcot = recoElectron->convDcot();
0125 float Dist = recoElectron->convDist();
0126 int NumberOfExpectedInnerHits =
0127 recoElectron->gsfTrack()->hitPattern().numberOfLostHits(reco::HitPattern::MISSING_INNER_HITS);
0128
0129
0130
0131 isIsolatedBarrel = false;
0132 isIDBarrel = false;
0133 isConvertedBarrel = false;
0134 isIsolatedEndcap = false;
0135 isIDEndcap = false;
0136 isConvertedEndcap = false;
0137
0138
0139
0140 if (fabs(recoElectron->eta()) <= 1.4442) {
0141
0142 if (IsoTrk < 0.09 && IsoEcal < 0.07 && IsoHcal < 0.10) {
0143 isIsolatedBarrel = true;
0144 }
0145
0146
0147 if (fabs(DeltaEtaTkClu) < 0.004 && fabs(DeltaPhiTkClu) < 0.06 && sigmaIeIe < 0.01 && HE < 0.04) {
0148 isIDBarrel = true;
0149 }
0150
0151
0152 if ((fabs(Dist) >= 0.02 || fabs(Dcot) >= 0.02) && NumberOfExpectedInnerHits <= 1.0) {
0153 isConvertedBarrel = true;
0154 }
0155 }
0156
0157 if (isIsolatedBarrel && isIDBarrel && isConvertedBarrel) {
0158 elIsAccepted++;
0159 elIsAcceptedEB++;
0160 TLorentzVector b_e2(
0161 recoElectron->momentum().x(), recoElectron->momentum().y(), recoElectron->momentum().z(), recoElectron->p());
0162 LV.push_back(b_e2);
0163 }
0164
0165
0166
0167 if (fabs(recoElectron->eta()) >= 1.5660 && fabs(recoElectron->eta()) <= 2.5000) {
0168
0169 if (IsoTrk < 0.04 && IsoEcal < 0.05 && IsoHcal < 0.025) {
0170 isIsolatedEndcap = true;
0171 }
0172
0173
0174 if (fabs(DeltaEtaTkClu) < 0.007 && fabs(DeltaPhiTkClu) < 0.03 && sigmaIeIe < 0.031 && HE < 0.15) {
0175 isIDEndcap = true;
0176 }
0177
0178
0179 if ((fabs(Dcot) > 0.02 || fabs(Dist) > 0.02) && NumberOfExpectedInnerHits <= 1.0) {
0180 isConvertedEndcap = true;
0181 }
0182 }
0183 if (isIsolatedEndcap && isIDEndcap && isConvertedEndcap) {
0184 elIsAccepted++;
0185 elIsAcceptedEE++;
0186 TLorentzVector e_e2(
0187 recoElectron->momentum().x(), recoElectron->momentum().y(), recoElectron->momentum().z(), recoElectron->p());
0188 LV.push_back(e_e2);
0189 }
0190 }
0191
0192
0193
0194 if (elIsAccepted > 1) {
0195 double e_ee_invMass = 0;
0196 if (LV.size() == 2) {
0197 TLorentzVector e_pair = LV[0] + LV[1];
0198 e_ee_invMass = e_pair.M();
0199 }
0200
0201 if (elIsAcceptedEB == 2) {
0202 h_ee_invMass_BB->Fill(e_ee_invMass);
0203 }
0204 if (elIsAcceptedEE == 2) {
0205 h_ee_invMass_EE->Fill(e_ee_invMass);
0206 }
0207 if (elIsAcceptedEB == 1 && elIsAcceptedEE == 1) {
0208 h_ee_invMass_EB->Fill(e_ee_invMass);
0209 }
0210
0211 LV.clear();
0212 }
0213 }
0214
0215 void EcalZmassTask::bookHistograms(DQMStore::IBooker &iBooker, edm::Run const &, edm::EventSetup const &) {
0216 std::string logTraceName("EcalZmassTask");
0217
0218 h_ee_invMass_EB = nullptr;
0219 h_ee_invMass_EE = nullptr;
0220 h_ee_invMass_BB = nullptr;
0221
0222 LogTrace(logTraceName) << "Parameters initialization";
0223
0224 iBooker.setCurrentFolder(prefixME_ + "/Zmass");
0225
0226 h_ee_invMass_EB = iBooker.book1D("Z peak - WP80 EB-EE", "Z peak - WP80 EB-EE;InvMass (GeV)", 60, 60.0, 120.0);
0227 h_ee_invMass_EE = iBooker.book1D("Z peak - WP80 EE-EE", "Z peak - WP80 EE-EE;InvMass (Gev)", 60, 60.0, 120.0);
0228 h_ee_invMass_BB = iBooker.book1D("Z peak - WP80 EB-EB", "Z peak - WP80 EB-EB;InvMass (Gev)", 60, 60.0, 120.0);
0229 }
0230
0231
0232 DEFINE_FWK_MODULE(EcalZmassTask);