CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_2_9_HLT1_bphpatch4/src/L1Trigger/DTTrigger/src/DTTrigProd.cc

Go to the documentation of this file.
00001 //-------------------------------------------------
00002 //
00013 //
00014 //--------------------------------------------------
00015 
00016 // This class's header
00017 #include "L1Trigger/DTTrigger/interface/DTTrigProd.h"
00018 
00019 // Framework related classes
00020 #include "FWCore/Framework/interface/EventSetup.h"
00021 #include "FWCore/Framework/interface/ESHandle.h"
00022 #include "FWCore/Framework/interface/Run.h"
00023 
00024 #include "L1TriggerConfig/DTTPGConfig/interface/DTConfigManager.h"
00025 #include "L1TriggerConfig/DTTPGConfig/interface/DTConfigManagerRcd.h"
00026 
00027 // Data Formats classes
00028 #include "L1Trigger/DTSectorCollector/interface/DTSectCollPhSegm.h"
00029 #include "L1Trigger/DTSectorCollector/interface/DTSectCollThSegm.h"
00030 #include "DataFormats/L1DTTrackFinder/interface/L1MuDTChambPhContainer.h"
00031 #include "DataFormats/L1DTTrackFinder/interface/L1MuDTChambPhDigi.h"
00032 #include "DataFormats/L1DTTrackFinder/interface/L1MuDTChambThContainer.h"
00033 #include "DataFormats/L1DTTrackFinder/interface/L1MuDTChambThDigi.h"
00034 
00035 // Collaborating classes
00036 #include <iostream>
00037 
00038 using namespace edm;
00039 using namespace std;
00040 
00041 // DataFormats interface
00042 typedef vector<DTSectCollPhSegm>  SectCollPhiColl;
00043 typedef SectCollPhiColl::const_iterator SectCollPhiColl_iterator;
00044 typedef vector<DTSectCollThSegm>  SectCollThetaColl;
00045 typedef SectCollThetaColl::const_iterator SectCollThetaColl_iterator;
00046 
00047 DTTrigProd::DTTrigProd(const ParameterSet& pset) : my_trig(0) {
00048   
00049   produces<L1MuDTChambPhContainer>();
00050   produces<L1MuDTChambThContainer>();
00051 
00052   my_debug = pset.getUntrackedParameter<bool>("debug");
00053   my_DTTFnum = pset.getParameter<bool>("DTTFSectorNumbering");
00054   my_params = pset;
00055 
00056   my_lut_dump_flag = pset.getUntrackedParameter<bool>("lutDumpFlag");
00057   my_lut_btic = pset.getUntrackedParameter<int>("lutBtic");
00058 
00059 }
00060 
00061 DTTrigProd::~DTTrigProd(){
00062 
00063   if (my_trig) delete my_trig;
00064 
00065 }
00066 
00067 void DTTrigProd::beginRun(edm::Run& iRun, const edm::EventSetup& iEventSetup) {
00068 
00069   if(my_debug)
00070     cout << "DTTrigProd::beginRun  " << iRun.id().run() << endl;
00071   
00072   ESHandle< DTConfigManager > dtConfig ;
00073   iEventSetup.get< DTConfigManagerRcd >().get( dtConfig ) ;
00074 
00075   my_CCBValid = dtConfig->CCBConfigValidity();
00076 
00077   if (!my_trig) {
00078     my_trig = new DTTrig(my_params);
00079     my_trig->createTUs(iEventSetup);
00080     if (my_debug)
00081       cout << "[DTTrigProd] TU's Created" << endl;
00082     
00083     if(my_lut_dump_flag) {
00084       cout << "Dumping luts...." << endl;
00085       my_trig->dumpLuts(my_lut_btic, dtConfig.product());
00086     }   
00087   }
00088 
00089 
00090   
00091 }
00092 
00093 
00094 void DTTrigProd::produce(Event & iEvent, const EventSetup& iEventSetup){
00095 
00096   // SV check if CCB configuration is valid, otherwise just produce empty collections
00097   if(!my_CCBValid)  {
00098     if (my_debug)
00099       cout << "[DTTrigProd] CCB configuration is not valid for this run, emulator will not run ! " << endl;
00100   }  else  {
00101     my_trig->triggerReco(iEvent,iEventSetup);
00102     my_BXoffset = my_trig->getBXOffset();
00103   
00104     if (my_debug)
00105       cout << "[DTTrigProd] Trigger algorithm run for " <<iEvent.id() << endl;
00106   
00107     // Convert Phi Segments
00108     SectCollPhiColl myPhiSegments;
00109     myPhiSegments = my_trig->SCPhTrigs();
00110     vector<L1MuDTChambPhDigi> outPhi;
00111 
00112     SectCollPhiColl_iterator SCPCend = myPhiSegments.end();
00113     for (SectCollPhiColl_iterator it=myPhiSegments.begin();it!=SCPCend;++it){
00114       int step = (*it).step() - my_BXoffset; // Shift correct BX to 0 (needed for DTTF data processing)
00115       int sc_sector = (*it).SCId().sector();
00116       if (my_DTTFnum == true) sc_sector--; // Modified for DTTF numbering [0-11]
00117         outPhi.push_back(L1MuDTChambPhDigi(step,
00118                                        (*it).ChamberId().wheel(),
00119                                        sc_sector,
00120                                        (*it).ChamberId().station(),
00121                                        (*it).phi(),
00122                                        (*it).phiB(),
00123                                        (*it).code(),
00124                                        !(*it).isFirst(),
00125                                        0
00126                                        ));
00127     }
00128 
00129     // Convert Theta Segments
00130     SectCollThetaColl myThetaSegments;
00131     myThetaSegments = my_trig->SCThTrigs();
00132     vector<L1MuDTChambThDigi> outTheta;
00133   
00134     SectCollThetaColl_iterator SCTCend = myThetaSegments.end();
00135     for (SectCollThetaColl_iterator it=myThetaSegments.begin();it!=SCTCend;++it){
00136       int pos[7], qual[7];
00137       for (int i=0; i<7; i++){
00138         pos[i] =(*it).position(i);
00139         qual[i]=(*it).quality(i);
00140       }
00141       int step =(*it).step() - my_BXoffset; // Shift correct BX to 0 (needed for DTTF data processing)
00142       int sc_sector =  (*it).SCId().sector();
00143       if (my_DTTFnum == true) sc_sector--; // Modified for DTTF numbering [0-11]
00144       outTheta.push_back(L1MuDTChambThDigi( step,
00145                                          (*it).ChamberId().wheel(),
00146                                          sc_sector,
00147                                          (*it).ChamberId().station(),
00148                                          pos,
00149                                          qual
00150                                          ));
00151     }
00152 
00153    // Write everything into the event
00154    std::auto_ptr<L1MuDTChambPhContainer> resultPhi (new L1MuDTChambPhContainer);
00155    resultPhi->setContainer(outPhi);
00156    iEvent.put(resultPhi);
00157    std::auto_ptr<L1MuDTChambThContainer> resultTheta (new L1MuDTChambThContainer);
00158    resultTheta->setContainer(outTheta);
00159    iEvent.put(resultTheta);
00160   }
00161 }
00162