Back to home page

Project CMSSW displayed by LXR

 
 

    


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

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