CMS 3D CMS Logo

FastL1MuonProducer.cc

Go to the documentation of this file.
00001 //
00002 // Package:    FastL1MuonProducer
00003 // Class:      FastL1MuonProducer
00004 // 
00012 //
00013 // Original Author:  Andrea Perrotta
00014 // Modifications: Patrick Janot.
00015 //         Created:  Mon Oct 30 14:37:24 CET 2006
00016 //
00017 //
00018 
00019 // CMSSW headers 
00020 #include "FWCore/ServiceRegistry/interface/Service.h"
00021 #include "FWCore/Utilities/interface/RandomNumberGenerator.h"
00022 #include "FWCore/Framework/interface/Event.h"
00023 #include "FWCore/Framework/interface/ESHandle.h"
00024 
00025 // Regional eta scales...
00026 #include "CondFormats/L1TObjects/interface/L1MuTriggerScales.h"
00027 #include "CondFormats/DataRecord/interface/L1MuTriggerScalesRcd.h"
00028 #include "CondFormats/L1TObjects/interface/L1MuTriggerPtScale.h"
00029 #include "CondFormats/DataRecord/interface/L1MuTriggerPtScaleRcd.h"
00030 
00031 // Fast Simulation headers
00032 #include "FastSimulation/Utilities/interface/RandomEngine.h"
00033 #include "FastSimulation/Muons/plugins/FastL1MuonProducer.h"
00034 
00035 // SimTrack
00036 #include "SimDataFormats/Track/interface/SimTrack.h"
00037 #include "SimDataFormats/TrackingHit/interface/PSimHitContainer.h"
00038 #include "SimDataFormats/TrackingHit/interface/PSimHit.h"
00039 
00040 // L1
00041 #include "FastSimDataFormats/L1GlobalMuonTrigger/interface/SimpleL1MuGMTCand.h"
00042 #include "FastSimulation/Muons/interface/FML1EfficiencyHandler.h"
00043 #include "FastSimulation/Muons/interface/FML1PtSmearer.h"
00044 
00045 // STL headers 
00046 #include <iostream>
00047 
00048 // Data Formats
00049 #include "DataFormats/L1GlobalMuonTrigger/interface/L1MuGMTReadoutCollection.h"
00050 #include "DataFormats/L1Trigger/interface/L1MuonParticle.h"
00051 #include "DataFormats/GeometryVector/interface/LocalPoint.h"
00052 #include "DataFormats/GeometryVector/interface/GlobalPoint.h"
00053 
00054 // Geometry
00055 #include "Geometry/Records/interface/MuonGeometryRecord.h"
00056 
00057 #include <map>
00058 
00059 // constants, enums and typedefs
00060 
00061 //
00062 // static data member definitions
00063 //
00064 
00065 double FastL1MuonProducer::muonMassGeV_ = 0.105658369 ; // PDG06
00066 
00067 //
00068 // constructors and destructor
00069 //
00070 FastL1MuonProducer::FastL1MuonProducer(const edm::ParameterSet& iConfig)
00071 {
00072 
00073   readParameters(iConfig.getParameter<edm::ParameterSet>("MUONS"));
00074 
00075   //register your products
00076   produces<L1MuonCollection> ();
00077   produces<L1ExtraCollection> ();
00078   produces<L1MuGMTReadoutCollection>();
00079 
00080   // Initialize the random number generator service
00081   edm::Service<edm::RandomNumberGenerator> rng;
00082   if ( ! rng.isAvailable() ) {
00083     throw cms::Exception("Configuration") <<
00084       "ParamMuonProducer requires the RandomGeneratorService \n"
00085       "which is not present in the configuration file. \n"
00086       "You must add the service in the configuration file\n"
00087       "or remove the module that requires it.";
00088   }
00089 
00090   random = new RandomEngine(&(*rng));
00091 
00092 }
00093 
00094 
00095 FastL1MuonProducer::~FastL1MuonProducer()
00096 {
00097  
00098   // do anything here that needs to be done at destruction time
00099   // (e.g. close files, deallocate resources etc.)
00100   
00101   if ( random ) { 
00102     delete random;
00103   }
00104 }
00105 
00106 
00107 //
00108 // member functions
00109 //
00110 
00111 // ------------ method called to produce the data  ------------
00112 
00113 void 
00114 FastL1MuonProducer::produce(edm::Event& iEvent, const edm::EventSetup& iSetup) 
00115 {
00116   using namespace edm;
00117 
00118   Handle<std::vector<SimTrack> > simMuons;
00119   iEvent.getByLabel(theSimModule,simMuons);
00120   unsigned nmuons = simMuons->size();
00121   //  Handle<std::vector<SimVertex> > simVertices;
00122   //  iEvent.getByLabel(theSimModuleLabel_,simVertices);
00123 
00124   Handle<PSimHitContainer> muonDTHits;
00125   iEvent.getByLabel(theDTHits,muonDTHits);
00126   Handle<PSimHitContainer> muonCSCHits;
00127   iEvent.getByLabel(theCSCHits,muonCSCHits);
00128   Handle<PSimHitContainer> muonRPCHits;
00129   iEvent.getByLabel(theRPCHits,muonRPCHits);
00130 
00131 
00132 //
00133 // Loop over generated muons and reconstruct L1, L3 and Global muons
00134 //
00135   
00136   int nMu = 0;
00137   mySimpleL1MuonCands.clear();
00138   mySimpleL1MuonExtraCands.clear();
00139 
00140   std::multimap<float,SimpleL1MuGMTCand*>  mySimpleL1MuonCandsTemp;
00141 
00142   for( unsigned fsimi=0; fsimi < nmuons; ++fsimi) {
00143     // The sim track can be a muon or a decaying hadron
00144     const SimTrack& mySimTrack = (*simMuons)[fsimi];
00145     // Keep only the muons at L1 (either primary or secondary)
00146     int pid = mySimTrack.type();        
00147     if ( fabs(pid) != 13 ) continue;
00148 
00149     // Check whether there are hits in DT/CSC/RPC,
00150     // and keep for the L1 mu the position of first such hit:
00151     bool hasPSimHits = false;
00152     GlobalPoint glbPosition;
00153     PSimHitContainer::const_iterator simDTHit=muonDTHits->begin();
00154     PSimHitContainer::const_iterator endDTHit=muonDTHits->end();
00155     for ( ; simDTHit!=endDTHit; ++simDTHit) {
00156       if ( simDTHit->trackId() == mySimTrack.trackId() ) {
00157         glbPosition = dtGeometry->idToDet(simDTHit->detUnitId())->surface().toGlobal(simDTHit->localPosition());
00158         hasPSimHits = true;
00159         break;
00160       }
00161     }
00162 
00163     if (! hasPSimHits) {
00164       PSimHitContainer::const_iterator simCSCHit=muonCSCHits->begin();
00165       PSimHitContainer::const_iterator endCSCHit=muonCSCHits->end();
00166       for ( ; simCSCHit!=endCSCHit; ++simCSCHit) {
00167         if ( simCSCHit->trackId() == mySimTrack.trackId() ) {
00168           glbPosition = cscGeometry->idToDet(simCSCHit->detUnitId())->surface().toGlobal(simCSCHit->localPosition());
00169           hasPSimHits = true;
00170           break;
00171         }
00172       }
00173     }
00174     if (! hasPSimHits) {
00175       PSimHitContainer::const_iterator simRPCHit=muonRPCHits->begin();
00176       PSimHitContainer::const_iterator endRPCHit=muonRPCHits->end();
00177       for ( ; simRPCHit!=endRPCHit; ++simRPCHit) {
00178         if ( simRPCHit->trackId() == mySimTrack.trackId() ) {
00179           glbPosition = rpcGeometry->idToDet(simRPCHit->detUnitId())->surface().toGlobal(simRPCHit->localPosition());
00180           hasPSimHits = true;
00181           break;
00182         }
00183       }
00184     }
00185 
00186     // *** Reconstruct parameterized muons starting from undecayed simulated muons
00187     if (hasPSimHits) {
00188 
00189       nMu++;
00190 
00191       //
00192       // Now L1 parametrization
00193       //
00194       double pT = mySimTrack.momentum().pt();
00195       double eta = glbPosition.eta();
00196       // Avoid L1MuScales complains if |eta|>2.4:
00197       if (eta > 2.4) eta = 2.4-1e-6; else if (eta < -2.4) eta = -2.4+1e-6;
00198       double phi = glbPosition.phi();
00199       if ( phi < 0. ) phi = 2* M_PI + phi;
00200 
00201       unsigned etaIndex = theMuScales->getGMTEtaScale()->getPacked(eta);
00202       unsigned phiIndex = theMuScales->getPhiScale()->getPacked(phi);
00203       unsigned pTIndex = theMuPtScale->getPtScale()->getPacked(pT);
00204       float etaValue = theMuScales->getGMTEtaScale()->getLowEdge(etaIndex);
00205       float phiValue = theMuScales->getPhiScale()->getLowEdge(phiIndex);
00206       float pTValue = theMuPtScale->getPtScale()->getLowEdge(pTIndex) + 1e-6;
00207       float etaValue2 = theMuScales->getGMTEtaScale()->getLowEdge(etaIndex+1);
00208       float phiValue2 = theMuScales->getPhiScale()->getLowEdge(phiIndex+1);
00209       float pTValue2 = theMuPtScale->getPtScale()->getLowEdge(pTIndex+1) + 1e-6;
00210      
00211       // Choose the closest index. (Not sure it is what is to be done)
00212       if ( fabs(etaValue2 - eta) < fabs(etaValue-eta) ) { 
00213         etaValue = etaValue2;
00214         ++etaIndex;
00215       }
00216       if ( fabs(phiValue2-phi) < fabs(phiValue-phi) ) { 
00217         phiValue = phiValue2;
00218         ++phiIndex;
00219       }
00220       if ( fabs(pTValue2-pT) < fabs(pTValue-pT) ) { 
00221         pTValue = pTValue2;
00222         ++pTIndex;
00223       }
00224 
00225       SimpleL1MuGMTCand * thisL1MuonCand = 
00226         new SimpleL1MuGMTCand(&mySimTrack,
00227                               etaIndex, phiIndex, pTIndex,
00228                               etaValue,phiValue,pTValue);
00229       bool hasL1 = myL1EfficiencyHandler->kill(thisL1MuonCand);
00230       if (hasL1) {
00231         bool status2 = myL1PtSmearer->smear(thisL1MuonCand);
00232         if (!status2) { std::cout << "Pt smearing of L1 muon went wrong!!" << std::endl; }
00233         if (status2) {
00234           mySimpleL1MuonCandsTemp.insert(
00235                std::pair<float,SimpleL1MuGMTCand*>(thisL1MuonCand->ptValue(),thisL1MuonCand));
00236         }
00237         else {
00238           delete thisL1MuonCand;
00239         }
00240       }
00241 
00242     }
00243 
00244   }
00245   
00246   // kill low ranked L1 muons, and fill L1extra muons -->
00247   std::multimap<float,SimpleL1MuGMTCand*>::const_reverse_iterator L1mu = mySimpleL1MuonCandsTemp.rbegin();
00248   std::multimap<float,SimpleL1MuGMTCand*>::const_reverse_iterator lastL1mu = mySimpleL1MuonCandsTemp.rend();
00249   unsigned rank=0;
00250   unsigned rankb=0;
00251   unsigned rankf=0;
00252   for ( ; L1mu!=lastL1mu; ++L1mu ) {
00253     SimpleL1MuGMTCand* theMuon = L1mu->second;
00254     theMuon->setRPCBit(0);
00255     ++rank;
00256     bool addMu = false;
00257     if  (theMuon->isFwd() ) { 
00258       if ( rankf < 4 ) addMu = true;
00259       theMuon->setRPCIndex(rankf);
00260       theMuon->setDTCSCIndex(rankf);
00261       rankf++;
00262     } else {  
00263       if ( rankb < 4 ) addMu = true;
00264       theMuon->setRPCIndex(rankb);
00265       theMuon->setDTCSCIndex(rankb);
00266       rankb++;
00267     } 
00268 
00269     if ( addMu ) {
00270       theMuon->setRank(rank);
00271       mySimpleL1MuonCands.push_back(theMuon);
00272       double pt = theMuon->ptValue() + 1.e-6 ;
00273       double eta = theMuon->etaValue();
00274       double phi = theMuon->phiValue();
00275       math::PtEtaPhiMLorentzVector PtEtaPhiMP4(pt,eta,phi,muonMassGeV_);
00276       math::XYZTLorentzVector myL1P4(PtEtaPhiMP4);
00277       // math::PtEtaPhiMLorentzVector myL1P4(pt,eta,phi,muonMassGeV_);
00278       mySimpleL1MuonExtraCands.push_back( l1extra::L1MuonParticle( theMuon->charge(), myL1P4, *theMuon ) );
00279     }
00280     else {
00281       theMuon->setRank(0);
00282     }
00283 
00284   }
00285 
00286 // end killing of low ranked L1 muons -->
00287 
00288 
00289   int nL1 =  mySimpleL1MuonCands.size();
00290   nMuonTot   += nMu;
00291   nL1MuonTot += nL1;
00292 
00293   std::auto_ptr<L1MuonCollection> l1Out(new L1MuonCollection);
00294   std::auto_ptr<L1ExtraCollection> l1ExtraOut(new L1ExtraCollection);
00295   std::auto_ptr<L1MuGMTReadoutCollection> l1ReadOut(new L1MuGMTReadoutCollection(1));
00296   
00297   loadL1Muons(*l1Out,*l1ExtraOut,*l1ReadOut);
00298   iEvent.put(l1Out);
00299   iEvent.put(l1ExtraOut);
00300   iEvent.put(l1ReadOut);
00301 
00302 }
00303 
00304 
00305 void FastL1MuonProducer::loadL1Muons(L1MuonCollection & c , 
00306                                      L1ExtraCollection & d,
00307                                      L1MuGMTReadoutCollection & e) const
00308 {
00309 
00310   FML1Muons::const_iterator l1mu;
00311   L1ExtraCollection::const_iterator l1ex;
00312   L1MuGMTReadoutRecord rc = L1MuGMTReadoutRecord(0);
00313 
00314   // Add L1 muons:
00315   // int nmuons = mySimpleL1MuonCands.size();
00316   for (l1mu=mySimpleL1MuonCands.begin();l1mu!=mySimpleL1MuonCands.end();++l1mu) {
00317     c.push_back(*(*l1mu));
00318   }
00319   
00320 // Add extra particles.
00321   int nr = 0; 
00322   int nrb = 0; 
00323   int nrf = 0;
00324   for (l1ex=mySimpleL1MuonExtraCands.begin();l1ex!=mySimpleL1MuonExtraCands.end();++l1ex) {
00325 
00326     d.push_back(*l1ex);
00327     
00328     L1MuGMTExtendedCand aMuon(l1ex->gmtMuonCand());
00329 
00330     rc.setGMTCand(nr,aMuon);
00331     // Find the regional eta index
00332     double etaPilePoil = mySimpleL1MuonCands[nr]->getMomentum().Eta();
00333 
00334     nr++;
00335     unsigned typeRPC=0;
00336     unsigned typeDTCSC=0;
00337     unsigned RPCIndex=0;
00338     unsigned DTCSCIndex=0;
00339     unsigned RPCRegionalEtaIndex=0;
00340     unsigned DTCSCRegionalEtaIndex=0;
00341     float etaRPCValue=-10.;
00342     float etaDTCSCValue=-10.;
00343     // Forward muons
00344     if ( aMuon.isFwd() ) { 
00345 
00346       rc.setGMTFwdCand(nrf,aMuon);
00347 
00348       // CSC
00349       typeDTCSC = 2;
00350       DTCSCIndex = 8+nrf;
00351       DTCSCRegionalEtaIndex = theMuScales->getRegionalEtaScale(2)->getPacked(etaPilePoil); 
00352       etaDTCSCValue = theMuScales->getRegionalEtaScale(2)->getLowEdge(DTCSCRegionalEtaIndex);
00353       float etaDTCSCValue2 = theMuScales->getRegionalEtaScale(2)->getLowEdge(DTCSCRegionalEtaIndex+1);
00354       if ( fabs(etaDTCSCValue2-etaPilePoil) < fabs(etaDTCSCValue-etaPilePoil) ) { 
00355         etaDTCSCValue = etaDTCSCValue2;
00356         ++DTCSCRegionalEtaIndex;
00357       }
00358       // RPC (limited to the RPC acceptance)
00359       if ( fabs(etaPilePoil) < 2.1 ) { 
00360         RPCIndex = 12+nrf;
00361         typeRPC = 3;
00362         RPCRegionalEtaIndex = theMuScales->getRegionalEtaScale(3)->getPacked(etaPilePoil); 
00363         etaRPCValue = theMuScales->getRegionalEtaScale(3)->getLowEdge(RPCRegionalEtaIndex);
00364         float etaRPCValue2 = theMuScales->getRegionalEtaScale(3)->getLowEdge(RPCRegionalEtaIndex+1);
00365         if ( fabs(etaRPCValue2-etaPilePoil) < fabs(etaRPCValue-etaPilePoil) ) { 
00366           etaRPCValue = etaRPCValue2;
00367           ++RPCRegionalEtaIndex;
00368         }
00369       }
00370       // Next muon
00371       nrf++;
00372 
00373     // Barrel muons
00374     } else { 
00375 
00376       rc.setGMTBrlCand(nrb,aMuon);
00377 
00378       // DT
00379       typeDTCSC = 0;
00380       DTCSCIndex = 0+nrb;
00381       DTCSCRegionalEtaIndex = theMuScales->getRegionalEtaScale(0)->getPacked(etaPilePoil); 
00382       etaDTCSCValue = theMuScales->getRegionalEtaScale(0)->getLowEdge(DTCSCRegionalEtaIndex);
00383       float etaDTCSCValue2 = theMuScales->getRegionalEtaScale(0)->getLowEdge(DTCSCRegionalEtaIndex+1);
00384       if ( fabs(etaDTCSCValue2-etaPilePoil) < fabs(etaDTCSCValue-etaPilePoil) ) { 
00385         etaDTCSCValue = etaDTCSCValue2;
00386         ++DTCSCRegionalEtaIndex;
00387       }
00388 
00389       // RPC
00390       typeRPC = 1;
00391       RPCIndex = 4+nrb;
00392       RPCRegionalEtaIndex = theMuScales->getRegionalEtaScale(1)->getPacked(etaPilePoil); 
00393       etaRPCValue = theMuScales->getRegionalEtaScale(1)->getLowEdge(RPCRegionalEtaIndex);
00394       float etaRPCValue2 = theMuScales->getRegionalEtaScale(1)->getLowEdge(RPCRegionalEtaIndex+1);
00395       if ( fabs(etaRPCValue2-etaPilePoil) < fabs(etaRPCValue-etaPilePoil) ) { 
00396         etaRPCValue = etaRPCValue2;
00397         ++RPCRegionalEtaIndex;
00398       }
00399 
00400       // Next muon
00401       nrb++;
00402 
00403     }
00404 
00405     // Add a muon regional candidate - first DT/CSC
00406     L1MuRegionalCand regionalMuonDTCSC = 
00407       L1MuRegionalCand(typeDTCSC,
00408                        aMuon.phiIndex(),
00409                        DTCSCRegionalEtaIndex,
00410                        aMuon.ptIndex(),
00411                        (1-aMuon.charge())/2,
00412                        aMuon.charge_valid(),
00413                        1,   // FineHalo
00414                        aMuon.quality(),
00415                        aMuon.bx());
00416     regionalMuonDTCSC.setPhiValue(aMuon.phiValue());
00417     regionalMuonDTCSC.setEtaValue(etaDTCSCValue);
00418     regionalMuonDTCSC.setPtValue(aMuon.ptValue());    
00419     
00420     rc.setInputCand(DTCSCIndex,regionalMuonDTCSC);
00421  
00422     // Then RPC (if in RPC acceptance)
00423     if ( fabs(etaPilePoil) < 2.1 ) { 
00424       L1MuRegionalCand regionalMuonRPC = 
00425         L1MuRegionalCand(typeRPC,
00426                          aMuon.phiIndex(),
00427                          RPCRegionalEtaIndex,
00428                          aMuon.ptIndex(),
00429                          (1-aMuon.charge())/2,
00430                          aMuon.charge_valid(),
00431                          0,   // FineHalo
00432                          aMuon.quality(),
00433                          aMuon.bx());
00434       regionalMuonRPC.setPhiValue(aMuon.phiValue());
00435       regionalMuonRPC.setEtaValue(etaRPCValue);
00436       regionalMuonRPC.setPtValue(aMuon.ptValue());
00437 
00438       rc.setInputCand(RPCIndex,regionalMuonRPC);
00439 
00440     }
00441 
00442   }
00443 
00444   // Update the event
00445   e.addRecord(rc);
00446 
00447 }
00448 
00449 // ------------ method called once each job just before starting event loop  ------------
00450 void 
00451 FastL1MuonProducer::beginJob(const edm::EventSetup& es)
00452 {
00453 
00454   // Initialize
00455   nMuonTot = 0;
00456 
00457   nL1MuonTot = 0;
00458   mySimpleL1MuonCands.clear();
00459   mySimpleL1MuonExtraCands.clear();
00460   myL1EfficiencyHandler = new FML1EfficiencyHandler(random);
00461   myL1PtSmearer = new FML1PtSmearer(random);
00462 
00463 // Get the DT Geometry
00464   es.get<MuonGeometryRecord>().get(dtGeometry);
00465 // Get the CSC Geometry
00466   es.get<MuonGeometryRecord>().get(cscGeometry);
00467 // Get the RPC Geometry
00468   es.get<MuonGeometryRecord>().get(rpcGeometry);
00469 
00470 }
00471 
00472 void 
00473 FastL1MuonProducer::beginRun(edm::Run & run, 
00474                              const edm::EventSetup & es) {
00475 
00476   // Read trigger scales
00477   edm::ESHandle< L1MuTriggerScales > muScales ;
00478   es.get< L1MuTriggerScalesRcd >().get( muScales ) ;
00479   theMuScales = &(*muScales);
00480 
00481   edm::ESHandle< L1MuTriggerPtScale > muPtScale ;
00482   es.get< L1MuTriggerPtScaleRcd >().get( muPtScale ) ;
00483   theMuPtScale = &(*muPtScale);
00484 
00485 }
00486 
00487 // ------------ method called once each job just after ending the event loop  ------------
00488 void 
00489 FastL1MuonProducer::endJob() {
00490 
00491   std::cout << " ===> FastL1MuonProducer , final report." << std::endl;
00492   std::cout << " ===> Number of total -> L1 muons in the whole run : "
00493             <<   nMuonTot << " -> " << nL1MuonTot << std::endl;
00494 }
00495 
00496 
00497 void FastL1MuonProducer::readParameters(const edm::ParameterSet& fastMuons) {
00498   // Muons
00499   theSimModule = fastMuons.getParameter<edm::InputTag>("simModule");
00500   theDTHits = fastMuons.getParameter<edm::InputTag>("dtSimHits");
00501   theCSCHits = fastMuons.getParameter<edm::InputTag>("cscSimHits");
00502   theRPCHits = fastMuons.getParameter<edm::InputTag>("rpcSimHits");
00503 
00504 }
00505 
00506 

Generated on Tue Jun 9 17:35:09 2009 for CMSSW by  doxygen 1.5.4