CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_4_5_patch3/src/L1Trigger/RegionalCaloTrigger/plugins/L1RCTTPGProvider.cc

Go to the documentation of this file.
00001 #include "L1Trigger/RegionalCaloTrigger/interface/L1RCTTPGProvider.h"
00002 #include "DataFormats/EcalDigi/interface/EcalDigiCollections.h"
00003 #include "DataFormats/HcalDigi/interface/HcalDigiCollections.h"
00004 
00005 L1RCTTPGProvider::L1RCTTPGProvider(const edm::ParameterSet& iConfig):
00006   ecalTPG_(iConfig.getParameter<edm::InputTag>("ecalTPGs")),
00007   hcalTPG_(iConfig.getParameter<edm::InputTag>("hcalTPGs")),
00008   useHcalCosmicTiming(iConfig.getParameter<bool>("useHCALCosmicTiming")),
00009   useEcalCosmicTiming(iConfig.getParameter<bool>("useECALCosmicTiming")),
00010   preSamples(iConfig.getParameter<int>("preSamples")),
00011   postSamples(iConfig.getParameter<int>("postSamples")),
00012   hfShift(iConfig.getParameter<int>("HFShift")),
00013   hbShift(iConfig.getParameter<int>("HBShift"))
00014 {
00015 
00016   //Output :The new manipulated TPGs
00017   //make it smart - to name the collections
00018   //correctly
00019     char ecal_label[200];
00020     char hcal_label[200];
00021 
00022 
00023 
00024     for(int i=preSamples;i>0;--i)
00025       {
00026         sprintf(ecal_label,"ECALBxminus%d",i);
00027         sprintf(hcal_label,"HCALBxminus%d",i);
00028         produces<EcalTrigPrimDigiCollection>(ecal_label);
00029         produces<HcalTrigPrimDigiCollection>(hcal_label);
00030       }
00031     
00032     produces<EcalTrigPrimDigiCollection>("ECALBx0");
00033     produces<HcalTrigPrimDigiCollection>("HCALBx0");
00034 
00035     for(int i=0;i<postSamples;++i)
00036       {
00037        sprintf(ecal_label,"ECALBxplus%d",i+1);
00038        sprintf(hcal_label,"HCALBxplus%d",i+1);
00039        produces<EcalTrigPrimDigiCollection>(ecal_label);
00040        produces<HcalTrigPrimDigiCollection>(hcal_label);
00041       }
00042 }
00043 
00044 L1RCTTPGProvider::~L1RCTTPGProvider()
00045 {
00046    // do anything here that needs to be done at desctruction time
00047    // (e.g. close files, deallocate resources etc.)
00048 }
00049 
00050 
00051 //
00052 // member functions
00053 //
00054 
00055 // ------------ method called to produce the data  ------------
00056 void
00057 L1RCTTPGProvider::produce(edm::Event& iEvent, const edm::EventSetup& iSetup)
00058 {
00059    using namespace edm;
00060    //Declare handles
00061    Handle<EcalTrigPrimDigiCollection> ecal;
00062    Handle<HcalTrigPrimDigiCollection> hcal;
00063 
00064    //Declare vector of collection to send for output !
00065 
00066    std::vector<EcalTrigPrimDigiCollection> ecalColl(preSamples+1+postSamples);
00067    std::vector<HcalTrigPrimDigiCollection> hcalColl(preSamples+1+postSamples);
00068 
00069   unsigned nSamples = preSamples + postSamples + 1;
00070 
00071 
00072   if(iEvent.getByLabel(ecalTPG_,ecal))
00073     if (ecal.isValid()) 
00074       { 
00075       // loop through all ecal digis
00076       for (EcalTrigPrimDigiCollection::const_iterator ecal_it = ecal->begin(); ecal_it != ecal->end(); ecal_it++)
00077         {
00078           short zside             = ecal_it->id().zside();
00079           unsigned short ietaAbs  = ecal_it->id().ietaAbs();
00080           short iphi              = ecal_it->id().iphi();
00081           unsigned short digiSize = ecal_it->size();
00082           unsigned short nSOI     = (unsigned short) ( ecal_it->sampleOfInterest() );
00083           if (digiSize < nSamples || nSOI < preSamples || ((int)(digiSize - nSOI) < (int)(nSamples - preSamples)))
00084             {
00085               unsigned short preLoopsZero = (unsigned short) (preSamples)- nSOI;
00086 
00087               // fill extra bx's at beginning with zeros
00088               for (int sample = 0; sample < preLoopsZero; sample++)
00089                 {
00090                   // fill first few with zeros
00091                   EcalTriggerPrimitiveDigi
00092                     ecalDigi(EcalTrigTowerDetId((int) zside, EcalTriggerTower,
00093                                                 (int) ietaAbs, (int) iphi));
00094                   ecalDigi.setSize(1);
00095                   ecalDigi.setSample(0, EcalTriggerPrimitiveSample(0,false,0));
00096                   ecalColl[sample].push_back(ecalDigi);
00097                 }
00098 
00099               // loop through existing data
00100               for (int sample = preLoopsZero;sample < (preLoopsZero + digiSize); sample++)
00101                 {
00102                   // go through data
00103                   EcalTriggerPrimitiveDigi
00104                     ecalDigi(EcalTrigTowerDetId((int) zside, EcalTriggerTower,
00105                                                 (int) ietaAbs, (int) iphi));
00106                   ecalDigi.setSize(1);
00107                   if (useEcalCosmicTiming && iphi >= 1 && iphi <= 36)
00108                     {
00109                       if (nSOI < (preSamples + 1))
00110                         {
00111                           edm::LogWarning ("TooLittleData")
00112                             << "ECAL data needs at least one presample "
00113                             << "more than the number requested "
00114                             << "to use ecal cosmic timing mod!  "
00115                             << "reverting to useEcalCosmicTiming = false "
00116                             << "for rest of job.";
00117                           useEcalCosmicTiming = false;
00118                         }
00119                       else
00120                         {
00121                           ecalDigi.setSample(0, EcalTriggerPrimitiveSample
00122                                              (ecal_it->sample(nSOI + sample - 
00123                                                               preSamples - 
00124                                                               1).raw()));
00125                         }
00126                     }
00127                   //else
00128                   if ((!useEcalCosmicTiming) || (iphi >=37 && iphi <= 72))
00129                     {
00130                       ecalDigi.setSample(0, EcalTriggerPrimitiveSample
00131                                          (ecal_it->sample(nSOI + sample - 
00132                                                           preSamples).raw()));
00133                     }
00134                   ecalColl[sample].push_back(ecalDigi);
00135                 }
00136               
00137               // fill extra bx's at end with zeros
00138               for (unsigned int sample = (preLoopsZero + digiSize); 
00139                    sample < nSamples; sample++)
00140                 {
00141                   // fill zeros!
00142                   EcalTriggerPrimitiveDigi
00143                     ecalDigi(EcalTrigTowerDetId((int) zside, EcalTriggerTower,
00144                                                 (int) ietaAbs, (int) iphi));
00145                   ecalDigi.setSize(1);
00146                   ecalDigi.setSample(0, EcalTriggerPrimitiveSample(0,false,0));
00147                   ecalColl[sample].push_back(ecalDigi);
00148                 }
00149             }
00150           else
00151             {
00152               for (unsigned short sample = 0; sample < nSamples; sample++)
00153                 {
00154                   // put each time sample into its own digi
00155                   short zside = ecal_it->id().zside();
00156                   unsigned short ietaAbs = ecal_it->id().ietaAbs();
00157                   short iphi = ecal_it->id().iphi();
00158                   EcalTriggerPrimitiveDigi
00159                     ecalDigi(EcalTrigTowerDetId((int) zside, EcalTriggerTower,
00160                                                 (int) ietaAbs, (int) iphi));
00161                   ecalDigi.setSize(1);
00162 
00163                   if (useEcalCosmicTiming && iphi >= 1 && iphi <=36)
00164                     {
00165                       if (nSOI < (preSamples + 1))
00166                         {
00167                           edm::LogWarning ("TooLittleData")
00168                             << "ECAL data needs at least one presample "
00169                             << "more than the number requested "
00170                             << "to use ecal cosmic timing mod!  "
00171                             << "reverting to useEcalCosmicTiming = false "
00172                             << "for rest of job.";
00173                           useEcalCosmicTiming = false;
00174                         }
00175                       else
00176                         {
00177                           ecalDigi.setSample(0, EcalTriggerPrimitiveSample
00178                                              (ecal_it->sample
00179                                               (ecal_it->sampleOfInterest() + 
00180                                                sample - preSamples - 
00181                                                1).raw()));
00182                         }
00183                     }
00184                   //else
00185                   if ((!useEcalCosmicTiming) || (iphi >=37 && iphi <= 72))
00186                     {
00187                       ecalDigi.setSample(0, EcalTriggerPrimitiveSample
00188                                          (ecal_it->sample
00189                                           (ecal_it->sampleOfInterest() + 
00190                                            sample - preSamples).raw()));
00191                     }
00192                   // push back each digi into correct "time sample" of coll
00193                   ecalColl[sample].push_back(ecalDigi);
00194                 }
00195             }
00196         }
00197     }
00198 
00199   if(iEvent.getByLabel(hcalTPG_,hcal))
00200     if (hcal.isValid()) 
00201     { 
00202       // loop through all hcal digis
00203       for (HcalTrigPrimDigiCollection::const_iterator hcal_it = hcal->begin(); hcal_it != hcal->end(); hcal_it++)
00204         {
00205           short ieta = hcal_it->id().ieta();
00206           short iphi = hcal_it->id().iphi();
00207           // loop through time samples for each digi
00208           unsigned short digiSize = hcal_it->size();
00209           // (size of each digi must be no less than nSamples)
00210           unsigned short nSOI = (unsigned short) (hcal_it->presamples());
00211           if (digiSize < nSamples || nSOI < preSamples
00212               || ((int)(digiSize - nSOI) < (int)(nSamples - preSamples)))
00213             {
00214               unsigned short preLoopsZero = (unsigned short) (preSamples) 
00215                 - nSOI;
00216               // fill extra bx's at beginning with zeros
00217               for (int sample = 0; sample < preLoopsZero; sample++)
00218                 {
00219                   // fill first few with zeros
00220                   HcalTriggerPrimitiveDigi hcalDigi(HcalTrigTowerDetId((int) ieta, (int) iphi));
00221                   hcalDigi.setSize(1);
00222                   hcalDigi.setPresamples(0);
00223                   hcalDigi.setSample(0, HcalTriggerPrimitiveSample(0,false,0,0));
00224                   hcalColl[sample].push_back(hcalDigi);
00225                 }
00226 
00227               // loop through existing data
00228               for (int sample = preLoopsZero; 
00229                    sample < (preLoopsZero + digiSize); sample++)
00230                 {
00231                   // go through data
00232                   HcalTriggerPrimitiveDigi hcalDigi(HcalTrigTowerDetId((int) ieta, (int) iphi));
00233                   hcalDigi.setSize(1);
00234                   hcalDigi.setPresamples(0);
00235 
00236                   if (useHcalCosmicTiming && iphi >= 1 && iphi <= 36)
00237                     {
00238                       if (nSOI < (preSamples + 1))
00239                         {
00240                           edm::LogWarning ("TooLittleData")
00241                             << "HCAL data needs at least one presample "
00242                             << "more than the number requested "
00243                             << "to use hcal cosmic timing mod!  "
00244                             << "reverting to useHcalCosmicTiming = false "
00245                             << "for rest of job.";
00246                           useHcalCosmicTiming = false;
00247                         }
00248                       else
00249                         {
00250                           hcalDigi.setSample(0, HcalTriggerPrimitiveSample
00251                                              (hcal_it->sample(hcal_it->
00252                                                               presamples() + 
00253                                                               sample - 
00254                                                               preSamples - 
00255                                                               1).raw()));
00256                         }
00257                     }
00258                   //else
00259                   if ((!useHcalCosmicTiming) || (iphi >= 37 && iphi <= 72))
00260                     {
00261 
00262                           hcalDigi.setSample(0, HcalTriggerPrimitiveSample
00263                                              (hcal_it->sample(hcal_it->
00264                                                      presamples() + 
00265                                                       sample - 
00266                                                       preSamples).raw()));
00267 
00268                     }
00269                   hcalColl[sample].push_back(hcalDigi);
00270                 }
00271               
00272               // fill extra bx's at end with zeros
00273               for (unsigned int sample = (preLoopsZero + digiSize); 
00274                    sample < nSamples; sample++)
00275                 {
00276                   // fill zeros!
00277                   HcalTriggerPrimitiveDigi
00278                     hcalDigi(HcalTrigTowerDetId((int) ieta, (int) iphi));
00279                   hcalDigi.setSize(1);
00280                   hcalDigi.setPresamples(0);
00281                   hcalDigi.setSample(0, HcalTriggerPrimitiveSample(0,false,0,0));
00282                   hcalColl[sample].push_back(hcalDigi);
00283                 }
00284             }
00285           else
00286             {
00287               for (unsigned short sample = 0; sample < nSamples; sample++)
00288                 {
00289                   // put each (relevant) time sample into its own digi
00290                   HcalTriggerPrimitiveDigi hcalDigi(HcalTrigTowerDetId((int) ieta, (int) iphi));
00291                   hcalDigi.setSize(1);
00292                   hcalDigi.setPresamples(0);
00293 
00294                   if (useHcalCosmicTiming && iphi >= 1 && iphi <= 36)
00295                     {
00296                       if (nSOI < (preSamples + 1))
00297                         {
00298                           edm::LogWarning ("TooLittleData")
00299                             << "HCAL data needs at least one presample "
00300                             << "more than the number requested "
00301                             << "to use hcal cosmic timing mod!  "
00302                             << "reverting to useHcalCosmicTiming = false "
00303                             << "for rest of job.";
00304                           useHcalCosmicTiming = false;
00305                         }
00306                       else
00307                         {
00308                           hcalDigi.setSample(0, HcalTriggerPrimitiveSample
00309                                              (hcal_it->sample(hcal_it->
00310                                                               presamples() + 
00311                                                               sample - 
00312                                                               preSamples - 
00313                                                               1).raw()));
00314                         }
00315                     }
00316                   //else
00317                   if ((!useHcalCosmicTiming) || (iphi >= 37 && iphi <= 72))
00318                     {
00319                           if(ieta>-29 && ieta<29) 
00320                             hcalDigi.setSample(0, HcalTriggerPrimitiveSample
00321                                                (hcal_it->sample(hcal_it->
00322                                                               presamples() + 
00323                                                                 sample - 
00324                                                                 preSamples+hbShift).raw()));
00325                           if(ieta<=-29 || ieta>=29)
00326                             hcalDigi.setSample(0, HcalTriggerPrimitiveSample
00327                                                (hcal_it->sample(hcal_it->
00328                                                               presamples() + 
00329                                                                 sample - 
00330                                                                 preSamples+hfShift).raw()));
00331                     }
00332                   hcalColl[sample].push_back(hcalDigi);  
00333                 }
00334             }
00335         }
00336     }
00337 
00338 
00339    //Now put the events back to file
00340 
00341   for(int i=0;i<preSamples;++i)
00342     {
00343       char ecal_label[200];
00344       char hcal_label[200];
00345 
00346 
00347       sprintf(ecal_label,"ECALBxminus%d",preSamples-i);
00348       sprintf(hcal_label,"HCALBxminus%d",preSamples-i);
00349 
00350       std::auto_ptr<EcalTrigPrimDigiCollection> ecalIn(new EcalTrigPrimDigiCollection); 
00351       std::auto_ptr<HcalTrigPrimDigiCollection> hcalIn(new HcalTrigPrimDigiCollection);
00352           for(unsigned int j=0;j<ecalColl[i].size();++j)
00353             {
00354               ecalIn->push_back((ecalColl[i])[j]);
00355             }
00356           for(unsigned int j=0;j<hcalColl[i].size();++j)
00357             hcalIn->push_back((hcalColl[i])[j]);
00358 
00359          iEvent.put(ecalIn,ecal_label);
00360          iEvent.put(hcalIn,hcal_label);
00361 
00362     }
00363 
00364 
00365    std::auto_ptr<EcalTrigPrimDigiCollection> ecal0(new EcalTrigPrimDigiCollection); 
00366    std::auto_ptr<HcalTrigPrimDigiCollection> hcal0(new HcalTrigPrimDigiCollection);
00367    for(unsigned int j=0;j<ecalColl[preSamples].size();++j)
00368      ecal0->push_back((ecalColl[preSamples])[j]);
00369    for(unsigned int j=0;j<hcalColl[preSamples].size();++j)
00370      hcal0->push_back((hcalColl[preSamples])[j]);
00371 
00372    iEvent.put(ecal0,"ECALBx0");
00373    iEvent.put(hcal0,"HCALBx0");
00374 
00375 
00376   for(int i=preSamples+1;i<preSamples+postSamples+1;++i)
00377     {
00378       char ecal_label[200];
00379       char hcal_label[200];
00380 
00381       sprintf(ecal_label,"ECALBxplus%d",i-preSamples);
00382       sprintf(hcal_label,"HCALBxplus%d",i-preSamples);
00383 
00384       std::auto_ptr<EcalTrigPrimDigiCollection> ecalIn2(new EcalTrigPrimDigiCollection); 
00385       std::auto_ptr<HcalTrigPrimDigiCollection> hcalIn2(new HcalTrigPrimDigiCollection);
00386 
00387       for(unsigned int j=0;j<ecalColl[i].size();++j)
00388         ecalIn2->push_back((ecalColl[i])[j]);
00389 
00390       for(unsigned int j=0;j<hcalColl[i].size();++j)
00391         hcalIn2->push_back((hcalColl[i])[j]);
00392 
00393 
00394          iEvent.put(ecalIn2,ecal_label);
00395          iEvent.put(hcalIn2,hcal_label);
00396     }
00397 }
00398 
00399 // ------------ method called once each job just before starting event loop  ------------
00400 
00401 
00402 // ------------ method called once each job just after ending the event loop  ------------
00403 void 
00404 L1RCTTPGProvider::endJob() {
00405 }