Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-12-20 03:13:56

0001 /*
0002  * EmulationObserverBase.cc
0003  *
0004  *  Created on: Aug 18, 2021
0005  *      Author: kbunkow
0006  */
0007 
0008 #include <memory>
0009 
0010 #include "L1Trigger/L1TMuonOverlapPhase1/interface/Tools/EmulationObserverBase.h"
0011 
0012 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0013 #include "Math/VectorUtil.h"
0014 
0015 EmulationObserverBase::EmulationObserverBase(const edm::ParameterSet& edmCfg, const OMTFConfiguration* omtfConfig)
0016     : edmCfg(edmCfg), omtfConfig(omtfConfig), simMuon(nullptr) {}
0017 
0018 EmulationObserverBase::~EmulationObserverBase() {
0019   // TODO Auto-generated destructor stub
0020 }
0021 
0022 void EmulationObserverBase::observeProcesorEmulation(unsigned int iProcessor,
0023                                                      l1t::tftype mtfType,
0024                                                      const std::shared_ptr<OMTFinput>& input,
0025                                                      const AlgoMuons& algoCandidates,
0026                                                      const AlgoMuons& gbCandidates,
0027                                                      const std::vector<l1t::RegionalMuonCand>& candMuons) {
0028   unsigned int procIndx = omtfConfig->getProcIndx(iProcessor, mtfType);
0029 
0030   unsigned int i = 0;
0031   for (auto& gbCandidate : gbCandidates) {
0032     if (gbCandidate->getGoldenPatern() != nullptr &&
0033         gbCandidate->getGpResultConstr().getFiredLayerCnt() > omtfCand->getGpResultConstr().getFiredLayerCnt()) {
0034       //LogTrace("l1tOmtfEventPrint") <<__FUNCTION__<<":"<<__LINE__<<" gbCandidate "<<gbCandidate<<" "<<std::endl;
0035       omtfCand = gbCandidate;
0036 
0037       candProcIndx = procIndx;
0038 
0039       //should be good, as the regionalMuonCand is created for every  gbCandidate in OMTFProcessor<GoldenPatternType>::getFinalcandidates
0040       regionalMuonCand = candMuons.at(i);
0041 
0042       //this->algoCandidates = algoCandidates; //TODO uncomment if needed
0043     }
0044     i++;
0045   }
0046 
0047   //////////////////////debug printout/////////////////////////////
0048   /*if(omtfCand->isValid()) { //TODO check this condition
0049     GoldenPatternWithStat* omtfCandGp = static_cast<GoldenPatternWithStat*>(omtfCand.getGoldenPatern());
0050     if( omtfCandGp->key().thePt > 100 && exptCandGp->key().thePt <= 15 ) {
0051       //LogTrace("l1tOmtfEventPrint")  <<iEvent.id()<<std::endl;
0052       LogTrace("l1tOmtfEventPrint") <<" ptSim "<<ptSim<<" chargeSim "<<chargeSim<<std::endl;
0053       LogTrace("l1tOmtfEventPrint") <<"iProcessor "<<iProcessor<<" exptCandGp "<<exptCandGp->key()<<std::endl;
0054       LogTrace("l1tOmtfEventPrint")  <<"iProcessor "<<iProcessor<<" omtfCandGp "<<omtfCandGp->key()<<std::endl;
0055       LogTrace("l1tOmtfEventPrint")  <<"omtfResult "<<std::endl<<omtfResult<<std::endl;
0056       int refHitNum = omtfCand.getRefHitNumber();
0057       LogTrace("l1tOmtfEventPrint")  <<"other gps results"<<endl;
0058       for(auto& gp : goldenPatterns) {
0059         if(omtfResult.getFiredLayerCnt() == gp->getResults()[procIndx][iRefHit].getFiredLayerCnt() )
0060         {
0061           LogTrace("l1tOmtfEventPrint")  <<gp->key()<<std::endl<<gp->getResults()[procIndx][iRefHit]<<std::endl;
0062         }
0063       }
0064       LogTrace("l1tOmtfEventPrint")<<std::endl;
0065     }
0066   }*/
0067 }
0068 
0069 void EmulationObserverBase::observeEventBegin(const edm::Event& iEvent) {
0070   omtfCand = std::make_shared<AlgoMuon>();
0071   candProcIndx = 0xffff;
0072 
0073   simMuon = findSimMuon(iEvent);
0074   //LogTrace("l1tOmtfEventPrint") <<__FUNCTION__<<":"<<__LINE__<<" evevt "<<iEvent.id().event()<<" simMuon pt "<<simMuon->momentum().pt()<<" GeV "<<std::endl;
0075 }
0076 
0077 const SimTrack* EmulationObserverBase::findSimMuon(const edm::Event& event, const SimTrack* previous) {
0078   const SimTrack* result = nullptr;
0079   if (edmCfg.exists("simTracksTag") == false)
0080     return result;
0081 
0082   edm::Handle<edm::SimTrackContainer> simTks;
0083   event.getByLabel(edmCfg.getParameter<edm::InputTag>("simTracksTag"), simTks);
0084 
0085   //LogTrace("l1tOmtfEventPrint")<<__FUNCTION__<<" simTks->size() "<<simTks->size()<<std::endl;
0086   for (std::vector<SimTrack>::const_iterator it = simTks->begin(); it < simTks->end(); it++) {
0087     const SimTrack& aTrack = *it;
0088     if (!(aTrack.type() == 13 || aTrack.type() == -13))
0089       continue;
0090     if (previous && ROOT::Math::VectorUtil::DeltaR(aTrack.momentum(), previous->momentum()) < 0.07)
0091       continue;
0092     if (!result || aTrack.momentum().pt() > result->momentum().pt())
0093       result = &aTrack;
0094   }
0095   return result;
0096 }
0097 
0098 const std::vector<const reco::GenParticle*> EmulationObserverBase::findGenMuon(const edm::Event& event) {
0099   std::vector<const reco::GenParticle*> muons;
0100 
0101   if (edmCfg.exists("genParticleTag") == false)
0102     return muons;
0103 
0104   edm::Handle<reco::GenParticleCollection> genParticles;
0105 
0106   event.getByLabel(edmCfg.getParameter<edm::InputTag>("genParticleTag"), genParticles);
0107 
0108   //todo
0109   float etaCutFrom = 0.8;
0110   float etaCutTo = 1.24;
0111 
0112   for (size_t i = 0; i < genParticles->size(); ++i) {
0113     const reco::GenParticle& genPart = (*genParticles)[i];
0114 
0115     if (std::abs(genPart.pdgId()) == 13) {
0116       if (std::abs(genPart.momentum().eta()) >= etaCutFrom && std::abs(genPart.momentum().eta()) <= etaCutTo) {
0117         genPart.momentum().eta();
0118 
0119         muons.push_back(&genPart);
0120         //int id = p.pdgId();
0121         LogTrace("l1tOmtfEventPrint") << __FUNCTION__ << ":" << __LINE__ << " genParticle: pdgId " << genPart.pdgId()
0122                                       << " px " << genPart.px() << " py " << genPart.py() << " pz " << genPart.pz()
0123                                       << " pt " << genPart.pt() << " vx " << genPart.vx() << " vy " << genPart.vy()
0124                                       << " vz " << genPart.vz() << " eta " << genPart.momentum().eta() << endl;
0125 
0126         LogTrace("l1tOmtfEventPrint") << __FUNCTION__ << ":" << __LINE__ << " genParticle: pdgId " << genPart.pdgId()
0127                                       << " px " << genPart.px() / genPart.p() << " py " << genPart.py() / genPart.p()
0128                                       << " pz " << genPart.pz() / genPart.p() << " pt " << genPart.pt() << endl;
0129         //<<" vx "<<genPart.vx()<<" vy "<<genPart.vy()<<" vz "<<genPart.vz()<<endl;
0130       }
0131     }
0132   }
0133 
0134   return muons;
0135 }