CMS 3D CMS Logo

/data/doxygen/doxygen-1.7.3/gen/CMSSW_4_2_8/src/L1Trigger/RegionalCaloTrigger/plugins/L1RCTProducer.cc

Go to the documentation of this file.
00001 #include "L1Trigger/RegionalCaloTrigger/interface/L1RCTProducer.h" 
00002 
00003 
00004 #include <vector>
00005 using std::vector;
00006 #include <iostream>
00007 
00008 
00009 
00010 using std::cout;
00011 using std::endl;
00012 const int L1RCTProducer::crateFED[18][5]=
00013       {{613, 614, 603, 702, 718},
00014     {611, 612, 602, 700, 718},
00015     {627, 610, 601,716,   722},
00016     {625, 626, 609, 714, 722},
00017     {623, 624, 608, 712, 722},
00018     {621, 622, 607, 710, 720},
00019     {619, 620, 606, 708, 720},
00020     {617, 618, 605, 706, 720},
00021     {615, 616, 604, 704, 718},
00022     {631, 632, 648, 703, 719},
00023     {629, 630, 647, 701, 719},
00024     {645, 628, 646, 717, 723},
00025     {643, 644, 654, 715, 723},
00026     {641, 642, 653, 713, 723},
00027     {639, 640, 652, 711, 721},
00028     {637, 638, 651, 709, 721},
00029     {635, 636, 650, 707, 721},
00030     {633, 634, 649, 705, 719}};
00031 
00032 
00033 
00034 L1RCTProducer::L1RCTProducer(const edm::ParameterSet& conf) : 
00035   rctLookupTables(new L1RCTLookupTables),
00036   rct(new L1RCT(rctLookupTables)),
00037   useEcal(conf.getParameter<bool>("useEcal")),
00038   useHcal(conf.getParameter<bool>("useHcal")),
00039   ecalDigis(conf.getParameter<std::vector<edm::InputTag> >("ecalDigis")),
00040   hcalDigis(conf.getParameter<std::vector<edm::InputTag> >("hcalDigis")),
00041   bunchCrossings(conf.getParameter<std::vector<int> >("BunchCrossings")),
00042   fedUpdatedMask(0)
00043 {
00044   produces<L1CaloEmCollection>();
00045   produces<L1CaloRegionCollection>();
00046 
00047  
00048 
00049 }
00050 
00051 L1RCTProducer::~L1RCTProducer()
00052 {
00053   if(rct != 0) delete rct;
00054   if(rctLookupTables != 0) delete rctLookupTables;
00055   if(fedUpdatedMask != 0) delete fedUpdatedMask;
00056 }
00057 
00058 
00059 void L1RCTProducer::beginRun(edm::Run& run, const edm::EventSetup& eventSetup)
00060 {
00061 
00062 }
00063 
00064 void L1RCTProducer::beginLuminosityBlock(edm::LuminosityBlock& lumiSeg,const edm::EventSetup& context)
00065 {
00066   updateConfiguration(context);
00067 } 
00068 
00069 
00070 
00071 void L1RCTProducer::updateConfiguration(const edm::EventSetup& eventSetup)
00072 {
00073   // Refresh configuration information every event
00074   // Hopefully, this does not take too much time
00075   // There should be a call back function in future to
00076   // handle changes in configuration
00077   // parameters to configure RCT (thresholds, etc)
00078   edm::ESHandle<L1RCTParameters> rctParameters;
00079   eventSetup.get<L1RCTParametersRcd>().get(rctParameters);
00080   const L1RCTParameters* r = rctParameters.product();
00081 
00082   // list of RCT channels to mask
00083   edm::ESHandle<L1RCTChannelMask> channelMask;
00084   eventSetup.get<L1RCTChannelMaskRcd>().get(channelMask);
00085   const L1RCTChannelMask* cEs = channelMask.product();
00086 
00087 
00088   // list of Noisy RCT channels to mask
00089   edm::ESHandle<L1RCTNoisyChannelMask> hotChannelMask;
00090   eventSetup.get<L1RCTNoisyChannelMaskRcd>().get(hotChannelMask);
00091   const L1RCTNoisyChannelMask* cEsNoise = hotChannelMask.product();
00092   rctLookupTables->setNoisyChannelMask(cEsNoise);
00093 
00094 
00095   
00096   //Update the channel mask according to the FED VECTOR
00097   //This is the beginning of run. We delete the old
00098   //create the new and set it in the LUTs
00099 
00100   if(fedUpdatedMask!=0) delete fedUpdatedMask;
00101 
00102   fedUpdatedMask = new L1RCTChannelMask();
00103   // copy a constant object
00104    for (int i = 0; i < 18; i++)
00105      {
00106        for (int j = 0; j < 2; j++)
00107          {
00108            for (int k = 0; k < 28; k++)
00109              {
00110                fedUpdatedMask->ecalMask[i][j][k] = cEs->ecalMask[i][j][k];
00111                fedUpdatedMask->hcalMask[i][j][k] =cEs->hcalMask[i][j][k] ;
00112              }
00113            for (int k = 0; k < 4; k++)
00114              {
00115                fedUpdatedMask->hfMask[i][j][k] = cEs->hfMask[i][j][k];
00116              }
00117          }
00118      }
00119 
00120 
00121   // adding fed mask into channel mask
00122   edm::ESHandle<RunInfo> sum;
00123   eventSetup.get<RunInfoRcd>().get(sum);
00124   const RunInfo* summary=sum.product();
00125 
00126 
00127   std::vector<int> caloFeds;  // pare down the feds to the intresting ones
00128   // is this unneccesary?
00129   // Mike B : This will decrease the find speed so better do it
00130   const std::vector<int> Feds = summary->m_fed_in;
00131   for(std::vector<int>::const_iterator cf = Feds.begin(); cf != Feds.end(); ++cf){
00132     int fedNum = *cf;
00133     if(fedNum > 600 && fedNum <724) 
00134       caloFeds.push_back(fedNum);
00135   }
00136 
00137   for(int  cr = 0; cr < 18; ++cr){
00138 
00139     for(crateSection cs = c_min; cs <= c_max; cs = crateSection(cs +1)) {
00140       bool fedFound = false;
00141 
00142 
00143       //Try to find the FED
00144       std::vector<int>::iterator fv = std::find(caloFeds.begin(),caloFeds.end(),crateFED[cr][cs]);
00145       if(fv!=caloFeds.end())
00146         fedFound = true;
00147 
00148       if(!fedFound) {
00149         int eta_min=0;
00150         int eta_max=0;
00151         bool phi_even[2] = {false};//, phi_odd = false;
00152         bool ecal=false;
00153         
00154         switch (cs) {
00155         case ebEvenFed :
00156           eta_min = minBarrel;
00157           eta_max = maxBarrel;
00158           phi_even[0] = true;
00159           ecal = true;  
00160           break;
00161           
00162         case ebOddFed:
00163           eta_min = minBarrel;
00164           eta_max = maxBarrel;
00165           phi_even[1] = true;
00166           ecal = true;  
00167           break;
00168           
00169         case eeFed:
00170           eta_min = minEndcap;
00171           eta_max = maxEndcap;
00172           phi_even[0] = true;
00173           phi_even[1] = true;
00174           ecal = true;  
00175           break;        
00176           
00177         case hbheFed:
00178           eta_min = minBarrel;
00179           eta_max = maxEndcap;
00180           phi_even[0] = true;
00181           phi_even[1] = true;
00182           ecal = false;
00183           break;
00184 
00185         case hfFed:     
00186           eta_min = minHF;
00187           eta_max = maxHF;
00188 
00189           phi_even[0] = true;
00190           phi_even[1] = true;
00191           ecal = false;
00192           break;
00193         default:
00194           break;
00195 
00196         }
00197         for(int ieta = eta_min; ieta <= eta_max; ++ieta){
00198           if(ieta<=28) // barrel and endcap
00199             for(int even = 0; even<=1 ; even++){         
00200               if(phi_even[even]){
00201                 if(ecal)
00202                   fedUpdatedMask->ecalMask[cr][even][ieta-1] = true;
00203                 else
00204                   fedUpdatedMask->hcalMask[cr][even][ieta-1] = true;
00205               }
00206             }
00207           else
00208             for(int even = 0; even<=1 ; even++)
00209               if(phi_even[even])
00210                   fedUpdatedMask->hfMask[cr][even][ieta-29] = true;
00211         
00212         }
00213       }
00214     }
00215   }
00216 
00217 
00218   //SCALES
00219 
00220   // energy scale to convert eGamma output
00221   edm::ESHandle<L1CaloEtScale> emScale;
00222   eventSetup.get<L1EmEtScaleRcd>().get(emScale);
00223   const L1CaloEtScale* s = emScale.product();
00224 
00225 
00226  // get energy scale to convert input from ECAL
00227   edm::ESHandle<L1CaloEcalScale> ecalScale;
00228   eventSetup.get<L1CaloEcalScaleRcd>().get(ecalScale);
00229   const L1CaloEcalScale* e = ecalScale.product();
00230   
00231   // get energy scale to convert input from HCAL
00232   edm::ESHandle<L1CaloHcalScale> hcalScale;
00233   eventSetup.get<L1CaloHcalScaleRcd>().get(hcalScale);
00234   const L1CaloHcalScale* h = hcalScale.product();
00235 
00236   // set scales
00237   rctLookupTables->setEcalScale(e);
00238   rctLookupTables->setHcalScale(h);
00239 
00240 
00241   rctLookupTables->setRCTParameters(r);
00242   rctLookupTables->setChannelMask(fedUpdatedMask);
00243   rctLookupTables->setL1CaloEtScale(s);
00244 }
00245 
00246 
00247 void L1RCTProducer::produce(edm::Event& event, const edm::EventSetup& eventSetup)
00248 {
00249 
00250 
00251   std::auto_ptr<L1CaloEmCollection> rctEmCands (new L1CaloEmCollection);
00252   std::auto_ptr<L1CaloRegionCollection> rctRegions (new L1CaloRegionCollection);
00253 
00254 
00255   if(!(ecalDigis.size()==hcalDigis.size()&&hcalDigis.size()==bunchCrossings.size()))
00256       throw cms::Exception("BadInput")
00257         << "From what I see the number of your your ECAL input digi collections.\n"
00258         <<"is different from the size of your HCAL digi input collections\n"
00259         <<"or the size of your BX factor collection" 
00260         <<"They must be the same to correspond to the same Bxs\n"
00261         << "It does not matter if one of them is empty\n"; 
00262 
00263 
00264 
00265 
00266   // loop through and process each bx
00267     for (unsigned short sample = 0; sample < bunchCrossings.size(); sample++)
00268       {
00269         edm::Handle<EcalTrigPrimDigiCollection> ecal;
00270         edm::Handle<HcalTrigPrimDigiCollection> hcal;
00271 
00272         EcalTrigPrimDigiCollection ecalIn;
00273         HcalTrigPrimDigiCollection hcalIn;
00274 
00275 
00276         if(useHcal&&event.getByLabel(hcalDigis[sample], hcal))
00277           hcalIn = *hcal;
00278 
00279         if(useEcal&&event.getByLabel(ecalDigis[sample],ecal))
00280           ecalIn = *ecal;
00281 
00282         rct->digiInput(ecalIn,hcalIn);
00283         rct->processEvent();
00284 
00285       // Stuff to create
00286         for (int j = 0; j<18; j++)
00287           {
00288             L1CaloEmCollection isolatedEGObjects = rct->getIsolatedEGObjects(j);
00289             L1CaloEmCollection nonisolatedEGObjects = rct->getNonisolatedEGObjects(j);
00290             for (int i = 0; i<4; i++) 
00291               {
00292                 isolatedEGObjects.at(i).setBx(bunchCrossings[sample]);
00293                 nonisolatedEGObjects.at(i).setBx(bunchCrossings[sample]);
00294                 rctEmCands->push_back(isolatedEGObjects.at(i));
00295                 rctEmCands->push_back(nonisolatedEGObjects.at(i));
00296               }
00297           }
00298       
00299       
00300         for (int i = 0; i < 18; i++)
00301           {
00302             std::vector<L1CaloRegion> regions = rct->getRegions(i);
00303             for (int j = 0; j < 22; j++)
00304               {
00305                 regions.at(j).setBx(bunchCrossings[sample]);
00306                 rctRegions->push_back(regions.at(j));
00307               }
00308           }
00309 
00310       }
00311 
00312   
00313   //putting stuff back into event
00314   event.put(rctEmCands);
00315   event.put(rctRegions);
00316   
00317 }