Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:20:19

0001 #include "FWCore/Framework/interface/Event.h"
0002 #include "FWCore/Framework/interface/MakerMacros.h"
0003 #include "FWCore/Framework/interface/one/EDAnalyzer.h"
0004 
0005 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0006 
0007 #include "FWCore/ServiceRegistry/interface/Service.h"
0008 #include "CommonTools/UtilAlgos/interface/TFileService.h"
0009 
0010 #include "CondFormats/L1TObjects/interface/CaloParams.h"
0011 #include "CondFormats/DataRecord/interface/L1TCaloParamsRcd.h"
0012 
0013 #include "DataFormats/L1TCalorimeter/interface/CaloTower.h"
0014 #include "DataFormats/L1TCalorimeter/interface/CaloCluster.h"
0015 #include "DataFormats/L1Trigger/interface/EGamma.h"
0016 #include "DataFormats/L1Trigger/interface/Tau.h"
0017 #include "DataFormats/L1Trigger/interface/Jet.h"
0018 #include "DataFormats/L1Trigger/interface/EtSum.h"
0019 
0020 #include "TH1F.h"
0021 #include "TH2F.h"
0022 
0023 //
0024 // class declaration
0025 //
0026 
0027 namespace l1t {
0028 
0029   class L1TStage2CaloAnalyzer : public edm::one::EDAnalyzer<edm::one::SharedResources> {
0030   public:
0031     explicit L1TStage2CaloAnalyzer(const edm::ParameterSet&);
0032 
0033     static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0034 
0035   private:
0036     void beginJob() override;
0037     void analyze(const edm::Event&, const edm::EventSetup&) override;
0038 
0039     // ----------member data ---------------------------
0040     edm::EDGetToken m_towerToken;
0041     edm::EDGetToken m_clusterToken;
0042     edm::EDGetToken m_mpEGToken;
0043     edm::EDGetToken m_mpTauToken;
0044     edm::EDGetToken m_mpJetToken;
0045     edm::EDGetToken m_mpSumToken;
0046     edm::EDGetToken m_egToken;
0047     edm::EDGetToken m_tauToken;
0048     edm::EDGetToken m_jetToken;
0049     edm::EDGetToken m_sumToken;
0050 
0051     bool m_doTowers;
0052     bool m_doClusters;
0053     bool m_doMPEGs;
0054     bool m_doMPTaus;
0055     bool m_doMPJets;
0056     bool m_doMPSums;
0057     bool m_doEGs;
0058     bool m_doTaus;
0059     bool m_doJets;
0060     bool m_doSums;
0061 
0062     bool doText_;
0063     bool doHistos_;
0064 
0065     enum ObjectType {
0066       Tower = 1,
0067       Cluster = 2,
0068       EG = 3,
0069       Tau = 4,
0070       Jet = 5,
0071       SumET = 6,
0072       SumETEm = 7,
0073       SumHT = 8,
0074       SumMET = 9,
0075       SumMETHF = 10,
0076       SumMHT = 11,
0077       SumMHTHF = 12,
0078       MPEG = 13,
0079       MPTau = 14,
0080       MPJet = 15,
0081       MPSum = 16,
0082       MPSumET = 17,
0083       MPSumETHF = 18,
0084       MPSumMETx = 19,
0085       MPSumMETxHF = 20,
0086       MPSumMETy = 21,
0087       MPSumMETyHF = 22,
0088       MPSumHT = 23,
0089       MPSumHTHF = 24,
0090       MPSumMHTx = 25,
0091       MPSumMHTxHF = 26,
0092       MPSumMHTy = 27,
0093       MPSumMHTyHF = 28,
0094       MPMinBiasHFP1 = 29,
0095       MPMinBiasHFM1 = 30,
0096       MPMinBiasHFP0 = 31,
0097       MPMinBiasHFM0 = 32,
0098       MinBiasHFP1 = 33,
0099       MinBiasHFM1 = 34,
0100       MinBiasHFP0 = 35,
0101       MinBiasHFM0 = 36,
0102       MPSumETEm = 37,
0103       MPSumHITowCount = 38,
0104       SumHITowCount = 39,
0105       SumCentrality = 40,
0106       SumAsymEt = 41,
0107       SumAsymEtHF = 42,
0108       SumAsymHt = 43,
0109       SumAsymHtHF = 44
0110     };
0111 
0112     std::vector<ObjectType> types_;
0113     std::vector<std::string> typeStr_;
0114 
0115     std::map<ObjectType, TFileDirectory> dirs_;
0116     std::map<ObjectType, TH1F*> het_;
0117     std::map<ObjectType, TH1F*> heta_;
0118     std::map<ObjectType, TH1F*> hphi_;
0119     std::map<ObjectType, TH1F*> hbx_;
0120     std::map<ObjectType, TH1F*> hem_;
0121     std::map<ObjectType, TH1F*> hhad_;
0122     std::map<ObjectType, TH1F*> hratio_;
0123     std::map<ObjectType, TH2F*> hetaphi_;
0124     std::map<ObjectType, TH1F*> hiso_;
0125 
0126     TFileDirectory evtDispDir_;
0127 
0128     TH1F *hsortMP_, *hsort_;
0129 
0130     int m_mpBx = 0;
0131     int m_dmxBx = 0;
0132     bool m_allBx = false;
0133     bool m_doEvtDisp = false;
0134   };
0135 
0136   //
0137   // constants, enums and typedefs
0138   //
0139 
0140   //
0141   // static data member definitions
0142   //
0143 
0144   //
0145   // constructors and destructor
0146   //
0147   L1TStage2CaloAnalyzer::L1TStage2CaloAnalyzer(const edm::ParameterSet& iConfig)
0148       : doText_(iConfig.getUntrackedParameter<bool>("doText", true)),
0149         doHistos_(iConfig.getUntrackedParameter<bool>("doHistos", true)) {
0150     usesResource(TFileService::kSharedResource);
0151     //now do what ever initialization is needed
0152 
0153     m_mpBx = iConfig.getParameter<int>("mpBx");
0154     m_dmxBx = iConfig.getParameter<int>("dmxBx");
0155     m_allBx = iConfig.getParameter<bool>("allBx");
0156     m_doEvtDisp = iConfig.getParameter<bool>("doEvtDisp");
0157 
0158     // register what you consume and keep token for later access:
0159     edm::InputTag nullTag("None");
0160 
0161     edm::InputTag towerTag = iConfig.getParameter<edm::InputTag>("towerToken");
0162     m_towerToken = consumes<l1t::CaloTowerBxCollection>(towerTag);
0163     m_doTowers = !(towerTag == nullTag);
0164 
0165     edm::InputTag clusterTag = iConfig.getParameter<edm::InputTag>("clusterToken");
0166     m_clusterToken = consumes<l1t::CaloClusterBxCollection>(clusterTag);
0167     m_doClusters = !(clusterTag == nullTag);
0168 
0169     edm::InputTag mpEGTag = iConfig.getParameter<edm::InputTag>("mpEGToken");
0170     m_mpEGToken = consumes<l1t::EGammaBxCollection>(mpEGTag);
0171     m_doMPEGs = !(mpEGTag == nullTag);
0172 
0173     edm::InputTag mpTauTag = iConfig.getParameter<edm::InputTag>("mpTauToken");
0174     m_mpTauToken = consumes<l1t::TauBxCollection>(mpTauTag);
0175     m_doMPTaus = !(mpTauTag == nullTag);
0176 
0177     edm::InputTag mpJetTag = iConfig.getParameter<edm::InputTag>("mpJetToken");
0178     m_mpJetToken = consumes<l1t::JetBxCollection>(mpJetTag);
0179     m_doMPJets = !(mpJetTag == nullTag);
0180 
0181     edm::InputTag mpSumTag = iConfig.getParameter<edm::InputTag>("mpEtSumToken");
0182     m_mpSumToken = consumes<l1t::EtSumBxCollection>(mpSumTag);
0183     m_doMPSums = !(mpSumTag == nullTag);
0184 
0185     edm::InputTag egTag = iConfig.getParameter<edm::InputTag>("egToken");
0186     m_egToken = consumes<l1t::EGammaBxCollection>(egTag);
0187     m_doEGs = !(egTag == nullTag);
0188 
0189     edm::InputTag tauTag = iConfig.getParameter<edm::InputTag>("tauToken");
0190     m_tauToken = consumes<l1t::TauBxCollection>(tauTag);
0191     m_doTaus = !(tauTag == nullTag);
0192 
0193     edm::InputTag jetTag = iConfig.getParameter<edm::InputTag>("jetToken");
0194     m_jetToken = consumes<l1t::JetBxCollection>(jetTag);
0195     m_doJets = !(jetTag == nullTag);
0196 
0197     edm::InputTag sumTag = iConfig.getParameter<edm::InputTag>("etSumToken");
0198     m_sumToken = consumes<l1t::EtSumBxCollection>(sumTag);
0199     m_doSums = !(sumTag == nullTag);
0200 
0201     std::cout << "Processing " << sumTag.label() << std::endl;
0202 
0203     types_.push_back(Tower);
0204     types_.push_back(Cluster);
0205     types_.push_back(MPEG);
0206     types_.push_back(MPTau);
0207     types_.push_back(MPJet);
0208     types_.push_back(MPSumET);
0209     types_.push_back(MPSumETHF);
0210     types_.push_back(MPSumETEm);
0211     types_.push_back(MPSumMETx);
0212     types_.push_back(MPSumMETxHF);
0213     types_.push_back(MPSumMETy);
0214     types_.push_back(MPSumMETyHF);
0215     types_.push_back(MPSumHT);
0216     types_.push_back(MPSumHTHF);
0217     types_.push_back(MPSumMHTx);
0218     types_.push_back(MPSumMHTxHF);
0219     types_.push_back(MPSumMHTy);
0220     types_.push_back(MPSumMHTyHF);
0221     types_.push_back(EG);
0222     types_.push_back(Tau);
0223     types_.push_back(Jet);
0224     types_.push_back(SumET);
0225     types_.push_back(SumETEm);
0226     types_.push_back(SumHT);
0227     types_.push_back(SumMET);
0228     types_.push_back(SumMETHF);
0229     types_.push_back(SumMHT);
0230     types_.push_back(SumMHTHF);
0231     types_.push_back(MPMinBiasHFP0);
0232     types_.push_back(MPMinBiasHFM0);
0233     types_.push_back(MPMinBiasHFP1);
0234     types_.push_back(MPMinBiasHFM1);
0235     types_.push_back(MinBiasHFP0);
0236     types_.push_back(MinBiasHFM0);
0237     types_.push_back(MinBiasHFP1);
0238     types_.push_back(MinBiasHFM1);
0239     types_.push_back(MPSumHITowCount);
0240     types_.push_back(SumHITowCount);
0241     types_.push_back(SumCentrality);
0242     types_.push_back(SumAsymEt);
0243     types_.push_back(SumAsymEtHF);
0244     types_.push_back(SumAsymHt);
0245     types_.push_back(SumAsymHtHF);
0246 
0247     typeStr_.push_back("tower");
0248     typeStr_.push_back("cluster");
0249     typeStr_.push_back("mpeg");
0250     typeStr_.push_back("mptau");
0251     typeStr_.push_back("mpjet");
0252     typeStr_.push_back("mpsumet");
0253     typeStr_.push_back("mpsumethf");
0254     typeStr_.push_back("mpsumetem");
0255     typeStr_.push_back("mpsummetx");
0256     typeStr_.push_back("mpsummetxhf");
0257     typeStr_.push_back("mpsummety");
0258     typeStr_.push_back("mpsummetyhf");
0259     typeStr_.push_back("mpsumht");
0260     typeStr_.push_back("mpsumhthf");
0261     typeStr_.push_back("mpsummhtx");
0262     typeStr_.push_back("mpsummhtxhf");
0263     typeStr_.push_back("mpsummhty");
0264     typeStr_.push_back("mpsummhtyhf");
0265     typeStr_.push_back("eg");
0266     typeStr_.push_back("tau");
0267     typeStr_.push_back("jet");
0268     typeStr_.push_back("sumet");
0269     typeStr_.push_back("sumetem");
0270     typeStr_.push_back("sumht");
0271     typeStr_.push_back("summet");
0272     typeStr_.push_back("summethf");
0273     typeStr_.push_back("summht");
0274     typeStr_.push_back("summhthf");
0275     typeStr_.push_back("mpminbiashfp0");
0276     typeStr_.push_back("mpminbiashfm0");
0277     typeStr_.push_back("mpminbiashfp1");
0278     typeStr_.push_back("mpminbiashfm1");
0279     typeStr_.push_back("minbiashfp0");
0280     typeStr_.push_back("minbiashfm0");
0281     typeStr_.push_back("minbiashfp1");
0282     typeStr_.push_back("minbiashfm1");
0283     typeStr_.push_back("mpsumhitowercount");
0284     typeStr_.push_back("sumhitowercount");
0285     typeStr_.push_back("sumcentrality");
0286     typeStr_.push_back("sumasymet");
0287     typeStr_.push_back("sumasymethf");
0288     typeStr_.push_back("sumasymht");
0289     typeStr_.push_back("sumasymhthf");
0290   }
0291 
0292   //
0293   // member functions
0294   //
0295 
0296   // ------------ method called for each event  ------------
0297   void L1TStage2CaloAnalyzer::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0298     using namespace edm;
0299 
0300     std::stringstream text;
0301 
0302     TH2F* hEvtTow = new TH2F();
0303     TH2F* hEvtMPEG = new TH2F();
0304     TH2F* hEvtMPTau = new TH2F();
0305     TH2F* hEvtMPJet = new TH2F();
0306     TH2F* hEvtDemuxEG = new TH2F();
0307     TH2F* hEvtDemuxTau = new TH2F();
0308     TH2F* hEvtDemuxJet = new TH2F();
0309 
0310     if (m_doEvtDisp) {
0311       std::stringstream ss;
0312       ss << iEvent.run() << "-" << iEvent.id().event();
0313       TFileDirectory dir = evtDispDir_.mkdir(ss.str());
0314       hEvtTow = dir.make<TH2F>("Tower", "", 83, -41.5, 41.5, 72, .5, 72.5);
0315       hEvtMPEG = dir.make<TH2F>("MPEG", "", 83, -41.5, 41.5, 72, .5, 72.5);
0316       hEvtMPTau = dir.make<TH2F>("MPTau", "", 83, -41.5, 41.5, 72, .5, 72.5);
0317       hEvtMPJet = dir.make<TH2F>("MPJet", "", 83, -41.5, 41.5, 72, .5, 72.5);
0318       hEvtDemuxEG = dir.make<TH2F>("DemuxEG", "", 227, -113.5, 113.5, 144, -0.5, 143.5);
0319       hEvtDemuxTau = dir.make<TH2F>("DemuxTau", "", 227, -113.5, 113.5, 144, -0.5, 143.5);
0320       hEvtDemuxJet = dir.make<TH2F>("DemuxJet", "", 227, -113.5, 113.5, 144, -0.5, 143.5);
0321     }
0322 
0323     // get TPs ?
0324     // get regions ?
0325     // get RCT clusters ?
0326 
0327     //check mpbx and dmxbx
0328     if (m_mpBx < -2 || m_mpBx > 2 || m_dmxBx < -2 || m_dmxBx > 2)
0329       edm::LogError("L1T")
0330           << "Selected MP Bx or Demux Bx to fill histograms is outside of range -2,2. Histos will be empty!";
0331 
0332     // get towers
0333     if (m_doTowers) {
0334       Handle<BXVector<l1t::CaloTower> > towers;
0335       iEvent.getByToken(m_towerToken, towers);
0336 
0337       for (int ibx = towers->getFirstBX(); ibx <= towers->getLastBX(); ++ibx) {
0338         if (!m_allBx && ibx != m_mpBx)
0339           continue;
0340 
0341         for (auto itr = towers->begin(ibx); itr != towers->end(ibx); ++itr) {
0342           if (itr->hwPt() <= 0)
0343             continue;
0344 
0345           hbx_.at(Tower)->Fill(ibx);
0346           het_.at(Tower)->Fill(itr->hwPt());
0347           heta_.at(Tower)->Fill(itr->hwEta());
0348           hphi_.at(Tower)->Fill(itr->hwPhi());
0349           hem_.at(Tower)->Fill(itr->hwEtEm());
0350           hhad_.at(Tower)->Fill(itr->hwEtHad());
0351           hratio_.at(Tower)->Fill(itr->hwEtRatio());
0352           hetaphi_.at(Tower)->Fill(itr->hwEta(), itr->hwPhi(), itr->hwPt());
0353 
0354           text << "Tower : "
0355                << " BX=" << ibx << " ipt=" << itr->hwPt() << " ieta=" << itr->hwEta() << " iphi=" << itr->hwPhi()
0356                << " iem=" << itr->hwEtEm() << " ihad=" << itr->hwEtHad() << " iratio=" << itr->hwEtRatio() << std::endl;
0357 
0358           if (m_doEvtDisp)
0359             hEvtTow->Fill(itr->hwEta(), itr->hwPhi(), itr->hwPt());
0360         }
0361       }
0362     }
0363     // get cluster
0364     if (m_doClusters) {
0365       Handle<BXVector<l1t::CaloCluster> > clusters;
0366       iEvent.getByToken(m_clusterToken, clusters);
0367 
0368       for (int ibx = clusters->getFirstBX(); ibx <= clusters->getLastBX(); ++ibx) {
0369         if (!m_allBx && ibx != m_mpBx)
0370           continue;
0371 
0372         for (auto itr = clusters->begin(ibx); itr != clusters->end(ibx); ++itr) {
0373           hbx_.at(Cluster)->Fill(ibx);
0374           het_.at(Cluster)->Fill(itr->hwPt());
0375           heta_.at(Cluster)->Fill(itr->hwEta());
0376           hphi_.at(Cluster)->Fill(itr->hwPhi());
0377           hetaphi_.at(Cluster)->Fill(itr->hwEta(), itr->hwPhi(), itr->hwPt());
0378           text << "Cluster : "
0379                << " BX=" << ibx << " ipt=" << itr->hwPt() << " ieta=" << itr->hwEta() << " iphi=" << itr->hwPhi()
0380                << std::endl;
0381         }
0382       }
0383     }
0384 
0385     // get EG
0386     if (m_doMPEGs) {
0387       Handle<BXVector<l1t::EGamma> > mpegs;
0388       iEvent.getByToken(m_mpEGToken, mpegs);
0389 
0390       for (int ibx = mpegs->getFirstBX(); ibx <= mpegs->getLastBX(); ++ibx) {
0391         if (!m_allBx && ibx != m_mpBx)
0392           continue;
0393 
0394         for (auto itr = mpegs->begin(ibx); itr != mpegs->end(ibx); ++itr) {
0395           hbx_.at(MPEG)->Fill(ibx);
0396           het_.at(MPEG)->Fill(itr->hwPt());
0397           heta_.at(MPEG)->Fill(itr->hwEta());
0398           hphi_.at(MPEG)->Fill(itr->hwPhi());
0399           hetaphi_.at(MPEG)->Fill(itr->hwEta(), itr->hwPhi(), itr->hwPt());
0400           hiso_.at(MPEG)->Fill(itr->hwIso());
0401 
0402           text << "MP EG : "
0403                << " BX=" << ibx << " ipt=" << itr->hwPt() << " ieta=" << itr->hwEta() << " iphi=" << itr->hwPhi()
0404                << std::endl;
0405 
0406           if (m_doEvtDisp)
0407             hEvtMPEG->Fill(itr->hwEta(), itr->hwPhi(), itr->hwPt());
0408         }
0409       }
0410     }
0411 
0412     // get tau
0413     if (m_doMPTaus) {
0414       Handle<BXVector<l1t::Tau> > mptaus;
0415       iEvent.getByToken(m_mpTauToken, mptaus);
0416 
0417       for (int ibx = mptaus->getFirstBX(); ibx <= mptaus->getLastBX(); ++ibx) {
0418         if (!m_allBx && ibx != m_mpBx)
0419           continue;
0420 
0421         for (auto itr = mptaus->begin(ibx); itr != mptaus->end(ibx); ++itr) {
0422           hbx_.at(MPTau)->Fill(ibx);
0423           het_.at(MPTau)->Fill(itr->hwPt());
0424           heta_.at(MPTau)->Fill(itr->hwEta());
0425           hphi_.at(MPTau)->Fill(itr->hwPhi());
0426           hetaphi_.at(MPTau)->Fill(itr->hwEta(), itr->hwPhi(), itr->hwPt());
0427 
0428           text << "MP Tau : "
0429                << " BX=" << ibx << " ipt=" << itr->hwPt() << " ieta=" << itr->hwEta() << " iphi=" << itr->hwPhi()
0430                << std::endl;
0431 
0432           if (m_doEvtDisp)
0433             hEvtMPTau->Fill(itr->hwEta(), itr->hwPhi(), itr->hwPt());
0434         }
0435       }
0436     }
0437 
0438     // get jet
0439     std::vector<l1t::Jet> thejets_poseta;
0440     std::vector<l1t::Jet> thejets_negeta;
0441 
0442     if (m_doMPJets) {
0443       Handle<BXVector<l1t::Jet> > mpjets;
0444       iEvent.getByToken(m_mpJetToken, mpjets);
0445 
0446       //Handle<BXVector<l1t::Jet>> jets;
0447       //iEvent.getByToken(m_jetToken,jets);
0448       //if (mpjets->size(0) == jets->size(0)) std::cout<<"******notequal"<<std::endl;
0449       for (int ibx = mpjets->getFirstBX(); ibx <= mpjets->getLastBX(); ++ibx) {
0450         if (!m_allBx && ibx != m_mpBx)
0451           continue;
0452 
0453         for (auto itr = mpjets->begin(ibx); itr != mpjets->end(ibx); ++itr) {
0454           hbx_.at(MPJet)->Fill(ibx);
0455           het_.at(MPJet)->Fill(itr->hwPt());
0456           heta_.at(MPJet)->Fill(itr->hwEta());
0457           hphi_.at(MPJet)->Fill(itr->hwPhi());
0458           hetaphi_.at(MPJet)->Fill(itr->hwEta(), itr->hwPhi(), itr->hwPt());
0459 
0460           text << "MP Jet : "
0461                << " BX=" << ibx << " ipt=" << itr->hwPt() << " ieta=" << itr->hwEta() << " iphi=" << itr->hwPhi()
0462                << std::endl;
0463 
0464           if (m_doEvtDisp)
0465             hEvtMPJet->Fill(itr->hwEta(), itr->hwPhi(), itr->hwPt());
0466 
0467           itr->hwEta() > 0 ? thejets_poseta.push_back(*itr) : thejets_negeta.push_back(*itr);
0468         }
0469       }
0470     }
0471 
0472     if (!thejets_poseta.empty()) {
0473       for (unsigned int i = 0; i < thejets_poseta.size() - 1; i++) {
0474         for (unsigned int j = i + 1; j < thejets_poseta.size(); j++) {
0475           hsortMP_->Fill(thejets_poseta.at(i).hwPt() - thejets_poseta.at(j).hwPt());
0476         }
0477       }
0478     }
0479 
0480     if (!thejets_negeta.empty()) {
0481       for (unsigned int i = 0; i < thejets_negeta.size() - 1; i++) {
0482         for (unsigned int j = i + 1; j < thejets_negeta.size(); j++) {
0483           hsortMP_->Fill(thejets_negeta.at(i).hwPt() - thejets_negeta.at(j).hwPt());
0484         }
0485       }
0486     }
0487 
0488     // get sums
0489     if (m_doMPSums) {
0490       Handle<BXVector<l1t::EtSum> > mpsums;
0491       iEvent.getByToken(m_mpSumToken, mpsums);
0492 
0493       for (int ibx = mpsums->getFirstBX(); ibx <= mpsums->getLastBX(); ++ibx) {
0494         if (!m_allBx && ibx != m_mpBx)
0495           continue;
0496 
0497         for (auto itr = mpsums->begin(ibx); itr != mpsums->end(ibx); ++itr) {
0498           switch (itr->getType()) {
0499             case l1t::EtSum::EtSumType::kTotalEt:
0500               het_.at(MPSumET)->Fill(itr->hwPt());
0501               break;
0502             case l1t::EtSum::EtSumType::kTotalEtHF:
0503               het_.at(MPSumETHF)->Fill(itr->hwPt());
0504               break;
0505             case l1t::EtSum::EtSumType::kTotalEtEm:
0506               het_.at(MPSumETEm)->Fill(itr->hwPt());
0507               break;
0508             case l1t::EtSum::EtSumType::kTotalEtx:
0509               het_.at(MPSumMETx)->Fill(itr->hwPt());
0510               break;
0511             case l1t::EtSum::EtSumType::kTotalEtxHF:
0512               het_.at(MPSumMETxHF)->Fill(itr->hwPt());
0513               break;
0514             case l1t::EtSum::EtSumType::kTotalEty:
0515               het_.at(MPSumMETy)->Fill(itr->hwPt());
0516               break;
0517             case l1t::EtSum::EtSumType::kTotalEtyHF:
0518               het_.at(MPSumMETyHF)->Fill(itr->hwPt());
0519               break;
0520             case l1t::EtSum::EtSumType::kTotalHt:
0521               het_.at(MPSumHT)->Fill(itr->hwPt());
0522               break;
0523             case l1t::EtSum::EtSumType::kTotalHtHF:
0524               het_.at(MPSumHTHF)->Fill(itr->hwPt());
0525               break;
0526             case l1t::EtSum::EtSumType::kTotalHtx:
0527               het_.at(MPSumMHTx)->Fill(itr->hwPt());
0528               break;
0529             case l1t::EtSum::EtSumType::kTotalHtxHF:
0530               het_.at(MPSumMHTxHF)->Fill(itr->hwPt());
0531               break;
0532             case l1t::EtSum::EtSumType::kTotalHty:
0533               het_.at(MPSumMHTy)->Fill(itr->hwPt());
0534               break;
0535             case l1t::EtSum::EtSumType::kTotalHtyHF:
0536               het_.at(MPSumMHTyHF)->Fill(itr->hwPt());
0537               break;
0538             case l1t::EtSum::EtSumType::kMinBiasHFP0:
0539               het_.at(MPMinBiasHFP0)->Fill(itr->hwPt());
0540               break;
0541             case l1t::EtSum::EtSumType::kMinBiasHFM0:
0542               het_.at(MPMinBiasHFM0)->Fill(itr->hwPt());
0543               break;
0544             case l1t::EtSum::EtSumType::kMinBiasHFP1:
0545               het_.at(MPMinBiasHFP1)->Fill(itr->hwPt());
0546               break;
0547             case l1t::EtSum::EtSumType::kMinBiasHFM1:
0548               het_.at(MPMinBiasHFM1)->Fill(itr->hwPt());
0549               break;
0550             case l1t::EtSum::EtSumType::kTowerCount:
0551               het_.at(MPSumHITowCount)->Fill(itr->hwPt());
0552               break;
0553 
0554             default:
0555               std::cout << "wrong type of MP sum" << std::endl;
0556           }
0557           text << "MP Sum : "
0558                << " type=" << itr->getType() << " BX=" << ibx << " ipt=" << itr->hwPt() << " ieta=" << itr->hwEta()
0559                << " iphi=" << itr->hwPhi() << std::endl;
0560         }
0561       }
0562     }
0563 
0564     // get EG
0565     if (m_doEGs) {
0566       Handle<BXVector<l1t::EGamma> > egs;
0567       iEvent.getByToken(m_egToken, egs);
0568 
0569       for (int ibx = egs->getFirstBX(); ibx <= egs->getLastBX(); ++ibx) {
0570         if (!m_allBx && ibx != m_dmxBx)
0571           continue;
0572 
0573         for (auto itr = egs->begin(ibx); itr != egs->end(ibx); ++itr) {
0574           hbx_.at(EG)->Fill(ibx);
0575           het_.at(EG)->Fill(itr->hwPt());
0576           heta_.at(EG)->Fill(itr->hwEta());
0577           hphi_.at(EG)->Fill(itr->hwPhi());
0578           hetaphi_.at(EG)->Fill(itr->hwEta(), itr->hwPhi(), itr->hwPt());
0579           hiso_.at(EG)->Fill(itr->hwIso());
0580 
0581           text << "EG : "
0582                << " BX=" << ibx << " ipt=" << itr->hwPt() << " ieta=" << itr->hwEta() << " iphi=" << itr->hwPhi()
0583                << std::endl;
0584 
0585           if (m_doEvtDisp)
0586             hEvtDemuxEG->Fill(itr->hwEta(), itr->hwPhi(), itr->hwPt());
0587         }
0588       }
0589     }
0590 
0591     // get tau
0592     if (m_doTaus) {
0593       Handle<BXVector<l1t::Tau> > taus;
0594       iEvent.getByToken(m_tauToken, taus);
0595 
0596       for (int ibx = taus->getFirstBX(); ibx <= taus->getLastBX(); ++ibx) {
0597         if (!m_allBx && ibx != m_dmxBx)
0598           continue;
0599 
0600         for (auto itr = taus->begin(ibx); itr != taus->end(ibx); ++itr) {
0601           hbx_.at(Tau)->Fill(ibx);
0602           het_.at(Tau)->Fill(itr->hwPt());
0603           heta_.at(Tau)->Fill(itr->hwEta());
0604           hphi_.at(Tau)->Fill(itr->hwPhi());
0605           hetaphi_.at(Tau)->Fill(itr->hwEta(), itr->hwPhi(), itr->hwPt());
0606 
0607           text << "Tau : "
0608                << " BX=" << ibx << " ipt=" << itr->hwPt() << " ieta=" << itr->hwEta() << " iphi=" << itr->hwPhi()
0609                << std::endl;
0610 
0611           if (m_doEvtDisp)
0612             hEvtDemuxTau->Fill(itr->hwEta(), itr->hwPhi(), itr->hwPt());
0613         }
0614       }
0615     }
0616 
0617     // get jet
0618     std::vector<l1t::Jet> thejets;
0619 
0620     if (m_doJets) {
0621       Handle<BXVector<l1t::Jet> > jets;
0622       iEvent.getByToken(m_jetToken, jets);
0623 
0624       for (int ibx = jets->getFirstBX(); ibx <= jets->getLastBX(); ++ibx) {
0625         if (!m_allBx && ibx != m_dmxBx)
0626           continue;
0627 
0628         for (auto itr = jets->begin(ibx); itr != jets->end(ibx); ++itr) {
0629           hbx_.at(Jet)->Fill(ibx);
0630           het_.at(Jet)->Fill(itr->hwPt());
0631           heta_.at(Jet)->Fill(itr->hwEta());
0632           hphi_.at(Jet)->Fill(itr->hwPhi());
0633           hetaphi_.at(Jet)->Fill(itr->hwEta(), itr->hwPhi(), itr->hwPt());
0634 
0635           text << "Jet : "
0636                << " BX=" << ibx << " ipt=" << itr->hwPt() << " ieta=" << itr->hwEta() << " iphi=" << itr->hwPhi()
0637                << std::endl;
0638 
0639           if (m_doEvtDisp)
0640             hEvtDemuxJet->Fill(itr->hwEta(), itr->hwPhi(), itr->hwPt());
0641 
0642           thejets.push_back(*itr);
0643         }
0644       }
0645     }
0646 
0647     if (!thejets.empty()) {
0648       for (unsigned int i = 0; i < thejets.size() - 1; i++) {
0649         for (unsigned int j = i + 1; j < thejets.size(); j++) {
0650           hsort_->Fill(thejets.at(i).hwPt() - thejets.at(j).hwPt());
0651         }
0652       }
0653     }
0654 
0655     // get sums
0656     if (m_doSums) {
0657       Handle<BXVector<l1t::EtSum> > sums;
0658       iEvent.getByToken(m_sumToken, sums);
0659 
0660       for (int ibx = sums->getFirstBX(); ibx <= sums->getLastBX(); ++ibx) {
0661         if (!m_allBx && ibx != m_dmxBx)
0662           continue;
0663 
0664         for (auto itr = sums->begin(ibx); itr != sums->end(ibx); ++itr) {
0665           switch (itr->getType()) {
0666             case l1t::EtSum::EtSumType::kTotalEt:
0667               het_.at(SumET)->Fill(itr->hwPt());
0668               break;
0669             case l1t::EtSum::EtSumType::kTotalEtHF:
0670               break;
0671             case l1t::EtSum::EtSumType::kTotalEtEm:
0672               het_.at(SumETEm)->Fill(itr->hwPt());
0673               break;
0674             case l1t::EtSum::EtSumType::kTotalHt:
0675               het_.at(SumHT)->Fill(itr->hwPt());
0676               break;
0677             case l1t::EtSum::EtSumType::kTotalHtHF:
0678               break;
0679             case l1t::EtSum::EtSumType::kMissingEt:
0680               het_.at(SumMET)->Fill(itr->hwPt());
0681               hphi_.at(SumMET)->Fill(itr->hwPhi());
0682               break;
0683             case l1t::EtSum::EtSumType::kMissingEtHF:
0684               het_.at(SumMETHF)->Fill(itr->hwPt());
0685               hphi_.at(SumMETHF)->Fill(itr->hwPhi());
0686               break;
0687             case l1t::EtSum::EtSumType::kMissingHt:
0688               het_.at(SumMHT)->Fill(itr->hwPt());
0689               hphi_.at(SumMHT)->Fill(itr->hwPhi());
0690               break;
0691             case l1t::EtSum::EtSumType::kMissingHtHF:
0692               het_.at(SumMHTHF)->Fill(itr->hwPt());
0693               hphi_.at(SumMHTHF)->Fill(itr->hwPhi());
0694               break;
0695             case l1t::EtSum::EtSumType::kMinBiasHFP0:
0696               het_.at(MinBiasHFP0)->Fill(itr->hwPt());
0697               break;
0698             case l1t::EtSum::EtSumType::kMinBiasHFM0:
0699               het_.at(MinBiasHFM0)->Fill(itr->hwPt());
0700               break;
0701             case l1t::EtSum::EtSumType::kMinBiasHFP1:
0702               het_.at(MinBiasHFP1)->Fill(itr->hwPt());
0703               break;
0704             case l1t::EtSum::EtSumType::kMinBiasHFM1:
0705               het_.at(MinBiasHFM1)->Fill(itr->hwPt());
0706               break;
0707             case l1t::EtSum::EtSumType::kTowerCount:
0708               het_.at(SumHITowCount)->Fill(itr->hwPt());
0709               break;
0710             case l1t::EtSum::EtSumType::kCentrality:
0711               het_.at(SumCentrality)->Fill(itr->hwPt());
0712               break;
0713             case l1t::EtSum::EtSumType::kAsymEt:
0714               het_.at(SumAsymEt)->Fill(itr->hwPt());
0715               break;
0716             case l1t::EtSum::EtSumType::kAsymHt:
0717               het_.at(SumAsymHt)->Fill(itr->hwPt());
0718               break;
0719             case l1t::EtSum::EtSumType::kAsymEtHF:
0720               het_.at(SumAsymEtHF)->Fill(itr->hwPt());
0721               break;
0722             case l1t::EtSum::EtSumType::kAsymHtHF:
0723               het_.at(SumAsymHtHF)->Fill(itr->hwPt());
0724               break;
0725 
0726             default:
0727               std::cout << "wrong type of demux sum: " << itr->getType() << std::endl;
0728           }
0729 
0730           text << "Sum : "
0731                << " type=" << itr->getType() << " BX=" << ibx << " ipt=" << itr->hwPt() << " ieta=" << itr->hwEta()
0732                << " iphi=" << itr->hwPhi() << std::endl;
0733         }
0734       }
0735     }
0736 
0737     if (doText_)
0738       edm::LogVerbatim("L1TCaloEvents") << text.str();
0739   }
0740 
0741   // ------------ method called once each job just before starting event loop  ------------
0742   void L1TStage2CaloAnalyzer::beginJob() {
0743     edm::Service<TFileService> fs;
0744 
0745     auto itr = types_.cbegin();
0746     auto str = typeStr_.cbegin();
0747 
0748     for (; itr != types_.end(); ++itr, ++str) {
0749       dirs_.insert(std::pair<ObjectType, TFileDirectory>(*itr, fs->mkdir(*str)));
0750 
0751       if (*itr == MPSumMETx || *itr == MPSumMETxHF || *itr == MPSumMETy || *itr == MPSumMETyHF) {
0752         het_.insert(
0753             std::pair<ObjectType, TH1F*>(*itr, dirs_.at(*itr).make<TH1F>("et", "", 2200000, -2199999500, 2200000500)));
0754       } else if (*itr == MPSumMHTx || *itr == MPSumMHTxHF || *itr == MPSumMHTy || *itr == MPSumMHTyHF) {
0755         het_.insert(
0756             std::pair<ObjectType, TH1F*>(*itr, dirs_.at(*itr).make<TH1F>("et", "", 2200000, -2199999500, 2200000500)));
0757       } else if (*itr == SumET || *itr == SumETEm || *itr == MPSumETEm || *itr == MPSumET || *itr == MPSumETHF ||
0758                  *itr == SumHT || *itr == MPSumHT || *itr == MPSumHTHF) {
0759         het_.insert(std::pair<ObjectType, TH1F*>(*itr, dirs_.at(*itr).make<TH1F>("et", "", 100000, -0.5, 99999.5)));
0760       } else if (*itr == MPMinBiasHFP0 || *itr == MPMinBiasHFM0 || *itr == MPMinBiasHFP1 || *itr == MPMinBiasHFM1 ||
0761                  *itr == MinBiasHFP0 || *itr == MinBiasHFM0 || *itr == MinBiasHFP1 || *itr == MinBiasHFM1) {
0762         het_.insert(std::pair<ObjectType, TH1F*>(*itr, dirs_.at(*itr).make<TH1F>("et", "", 16, -0.5, 15.5)));
0763       } else if (*itr == MPSumHITowCount || *itr == SumHITowCount) {
0764         het_.insert(std::pair<ObjectType, TH1F*>(*itr, dirs_.at(*itr).make<TH1F>("et", "", 5904, -0.5, 5903.5)));
0765       } else if (*itr == EG || *itr == Tau || *itr == MPEG || *itr == MPTau) {
0766         het_.insert(std::pair<ObjectType, TH1F*>(*itr, dirs_.at(*itr).make<TH1F>("et", "", 1500, -0.5, 1499.5)));
0767       } else {
0768         het_.insert(std::pair<ObjectType, TH1F*>(*itr, dirs_.at(*itr).make<TH1F>("et", "", 100000, -0.5, 99999.5)));
0769       }
0770 
0771       hbx_.insert(std::pair<ObjectType, TH1F*>(*itr, dirs_.at(*itr).make<TH1F>("bx", "", 11, -5.5, 5.5)));
0772 
0773       if (*itr == EG || *itr == Jet || *itr == Tau) {
0774         heta_.insert(std::pair<ObjectType, TH1F*>(*itr, dirs_.at(*itr).make<TH1F>("eta", "", 227, -113.5, 113.5)));
0775         hphi_.insert(std::pair<ObjectType, TH1F*>(*itr, dirs_.at(*itr).make<TH1F>("phi", "", 144, -0.5, 143.5)));
0776         hetaphi_.insert(std::pair<ObjectType, TH2F*>(
0777             *itr, dirs_.at(*itr).make<TH2F>("etaphi", "", 227, -113.5, 113.5, 144, -0.5, 143.5)));
0778         if (*itr == EG)
0779           hiso_.insert(std::pair<ObjectType, TH1F*>(*itr, dirs_.at(*itr).make<TH1F>("iso", "", 4, -0.5, 3.5)));
0780       } else if (*itr == Tower || *itr == Cluster || *itr == MPEG || *itr == MPJet || *itr == MPTau) {
0781         heta_.insert(std::pair<ObjectType, TH1F*>(*itr, dirs_.at(*itr).make<TH1F>("eta", "", 83, -41.5, 41.5)));
0782         hphi_.insert(std::pair<ObjectType, TH1F*>(*itr, dirs_.at(*itr).make<TH1F>("phi", "", 72, 0.5, 72.5)));
0783         hetaphi_.insert(
0784             std::pair<ObjectType, TH2F*>(*itr, dirs_.at(*itr).make<TH2F>("etaphi", "", 83, -41.5, 41.5, 72, .5, 72.5)));
0785         if (*itr == MPEG)
0786           hiso_.insert(std::pair<ObjectType, TH1F*>(*itr, dirs_.at(*itr).make<TH1F>("iso", "", 4, -0.5, 3.5)));
0787       } else if (*itr == SumMET || *itr == SumMETHF || *itr == SumMHT || *itr == SumMHTHF) {
0788         hphi_.insert(std::pair<ObjectType, TH1F*>(*itr, dirs_.at(*itr).make<TH1F>("phi", "", 1008, -0.5, 1007.5)));
0789       }
0790 
0791       if (*itr == Tower) {
0792         hem_.insert(std::pair<ObjectType, TH1F*>(*itr, dirs_.at(*itr).make<TH1F>("em", "", 101, -0.5, 100.5)));
0793         hhad_.insert(std::pair<ObjectType, TH1F*>(*itr, dirs_.at(*itr).make<TH1F>("had", "", 101, -0.5, 100.5)));
0794         hratio_.insert(std::pair<ObjectType, TH1F*>(*itr, dirs_.at(*itr).make<TH1F>("ratio", "", 11, -0.5, 10.5)));
0795       }
0796     }
0797 
0798     if (m_doEvtDisp) {
0799       evtDispDir_ = fs->mkdir("Events");
0800     }
0801 
0802     hsort_ = fs->make<TH1F>("sort", "", 201, -100.5, 100.5);
0803     hsortMP_ = fs->make<TH1F>("sortMP", "", 201, -100.5, 100.5);
0804   }
0805 
0806   // ------------ method fills 'descriptions' with the allowed parameters for the module  ------------
0807   void L1TStage2CaloAnalyzer::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0808     //The following says we do not know what parameters are allowed so do no validation
0809     // Please change this to state exactly what you do use, even if it is no parameters
0810     edm::ParameterSetDescription desc;
0811     desc.setUnknown();
0812     descriptions.addDefault(desc);
0813   }
0814 
0815 }  // namespace l1t
0816 
0817 using namespace l1t;
0818 
0819 //define this as a plug-in
0820 DEFINE_FWK_MODULE(L1TStage2CaloAnalyzer);