Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 #include "DataFormats/EcalDigi/interface/EcalDigiCollections.h"
0002 #include "DataFormats/HcalDigi/interface/HcalDigiCollections.h"
0003 #include "L1Trigger/RegionalCaloTrigger/interface/L1RCTTPGProvider.h"
0004 
0005 L1RCTTPGProvider::L1RCTTPGProvider(const edm::ParameterSet &iConfig)
0006     : ecalTPG_(consumes(iConfig.getParameter<edm::InputTag>("ecalTPGs"))),
0007       hcalTPG_(consumes(iConfig.getParameter<edm::InputTag>("hcalTPGs"))),
0008       useHcalCosmicTiming(iConfig.getParameter<bool>("useHCALCosmicTiming")),
0009       useEcalCosmicTiming(iConfig.getParameter<bool>("useECALCosmicTiming")),
0010       preSamples(iConfig.getParameter<int>("preSamples")),
0011       postSamples(iConfig.getParameter<int>("postSamples")),
0012       hfShift(iConfig.getParameter<int>("HFShift")),
0013       hbShift(iConfig.getParameter<int>("HBShift")) {
0014   // Output :The new manipulated TPGs
0015   // make it smart - to name the collections
0016   // correctly
0017   char ecal_label[200];
0018   char hcal_label[200];
0019 
0020   for (int i = preSamples; i > 0; --i) {
0021     sprintf(ecal_label, "ECALBxminus%d", i);
0022     sprintf(hcal_label, "HCALBxminus%d", i);
0023     produces<EcalTrigPrimDigiCollection>(ecal_label);
0024     produces<HcalTrigPrimDigiCollection>(hcal_label);
0025   }
0026 
0027   produces<EcalTrigPrimDigiCollection>("ECALBx0");
0028   produces<HcalTrigPrimDigiCollection>("HCALBx0");
0029 
0030   for (int i = 0; i < postSamples; ++i) {
0031     sprintf(ecal_label, "ECALBxplus%d", i + 1);
0032     sprintf(hcal_label, "HCALBxplus%d", i + 1);
0033     produces<EcalTrigPrimDigiCollection>(ecal_label);
0034     produces<HcalTrigPrimDigiCollection>(hcal_label);
0035   }
0036 }
0037 
0038 L1RCTTPGProvider::~L1RCTTPGProvider() {
0039   // do anything here that needs to be done at desctruction time
0040   // (e.g. close files, deallocate resources etc.)
0041 }
0042 
0043 //
0044 // member functions
0045 //
0046 
0047 // ------------ method called to produce the data  ------------
0048 void L1RCTTPGProvider::produce(edm::StreamID, edm::Event &iEvent, const edm::EventSetup &iSetup) const {
0049   using namespace edm;
0050   // Declare handles
0051   Handle<EcalTrigPrimDigiCollection> ecal;
0052   Handle<HcalTrigPrimDigiCollection> hcal;
0053 
0054   // Declare vector of collection to send for output !
0055 
0056   std::vector<EcalTrigPrimDigiCollection> ecalColl(preSamples + 1 + postSamples);
0057   std::vector<HcalTrigPrimDigiCollection> hcalColl(preSamples + 1 + postSamples);
0058 
0059   unsigned nSamples = preSamples + postSamples + 1;
0060 
0061   if ((ecal = iEvent.getHandle(ecalTPG_)))
0062     if (ecal.isValid()) {
0063       // loop through all ecal digis
0064       for (EcalTrigPrimDigiCollection::const_iterator ecal_it = ecal->begin(); ecal_it != ecal->end(); ecal_it++) {
0065         short zside = ecal_it->id().zside();
0066         unsigned short ietaAbs = ecal_it->id().ietaAbs();
0067         short iphi = ecal_it->id().iphi();
0068         unsigned short digiSize = ecal_it->size();
0069         unsigned short nSOI = (unsigned short)(ecal_it->sampleOfInterest());
0070         if (digiSize < nSamples || nSOI < preSamples || ((int)(digiSize - nSOI) < (int)(nSamples - preSamples))) {
0071           unsigned short preLoopsZero = (unsigned short)(preSamples)-nSOI;
0072 
0073           // fill extra bx's at beginning with zeros
0074           for (int sample = 0; sample < preLoopsZero; sample++) {
0075             // fill first few with zeros
0076             EcalTriggerPrimitiveDigi ecalDigi(
0077                 EcalTrigTowerDetId((int)zside, EcalTriggerTower, (int)ietaAbs, (int)iphi));
0078             ecalDigi.setSize(1);
0079             ecalDigi.setSample(0, EcalTriggerPrimitiveSample(0, false, 0));
0080             ecalColl[sample].push_back(ecalDigi);
0081           }
0082 
0083           // loop through existing data
0084           for (int sample = preLoopsZero; sample < (preLoopsZero + digiSize); sample++) {
0085             // go through data
0086             EcalTriggerPrimitiveDigi ecalDigi(
0087                 EcalTrigTowerDetId((int)zside, EcalTriggerTower, (int)ietaAbs, (int)iphi));
0088             ecalDigi.setSize(1);
0089             if (useEcalCosmicTiming && iphi >= 1 && iphi <= 36) {
0090               if (nSOI < (preSamples + 1)) {
0091                 edm::LogWarning("TooLittleData") << "ECAL data needs at least one presample "
0092                                                  << "more than the number requested "
0093                                                  << "to use ecal cosmic timing mod!  "
0094                                                  << "reverting to useEcalCosmicTiming = false "
0095                                                  << "for rest of job.";
0096                 useEcalCosmicTiming = false;
0097               } else {
0098                 ecalDigi.setSample(0,
0099                                    EcalTriggerPrimitiveSample(ecal_it->sample(nSOI + sample - preSamples - 1).raw()));
0100               }
0101             }
0102             // else
0103             if ((!useEcalCosmicTiming) || (iphi >= 37 && iphi <= 72)) {
0104               ecalDigi.setSample(0, EcalTriggerPrimitiveSample(ecal_it->sample(nSOI + sample - preSamples).raw()));
0105             }
0106             ecalColl[sample].push_back(ecalDigi);
0107           }
0108 
0109           // fill extra bx's at end with zeros
0110           for (unsigned int sample = (preLoopsZero + digiSize); sample < nSamples; sample++) {
0111             // fill zeros!
0112             EcalTriggerPrimitiveDigi ecalDigi(
0113                 EcalTrigTowerDetId((int)zside, EcalTriggerTower, (int)ietaAbs, (int)iphi));
0114             ecalDigi.setSize(1);
0115             ecalDigi.setSample(0, EcalTriggerPrimitiveSample(0, false, 0));
0116             ecalColl[sample].push_back(ecalDigi);
0117           }
0118         } else {
0119           for (unsigned short sample = 0; sample < nSamples; sample++) {
0120             // put each time sample into its own digi
0121             short zside = ecal_it->id().zside();
0122             unsigned short ietaAbs = ecal_it->id().ietaAbs();
0123             short iphi = ecal_it->id().iphi();
0124             EcalTriggerPrimitiveDigi ecalDigi(
0125                 EcalTrigTowerDetId((int)zside, EcalTriggerTower, (int)ietaAbs, (int)iphi));
0126             ecalDigi.setSize(1);
0127 
0128             if (useEcalCosmicTiming && iphi >= 1 && iphi <= 36) {
0129               if (nSOI < (preSamples + 1)) {
0130                 edm::LogWarning("TooLittleData") << "ECAL data needs at least one presample "
0131                                                  << "more than the number requested "
0132                                                  << "to use ecal cosmic timing mod!  "
0133                                                  << "reverting to useEcalCosmicTiming = false "
0134                                                  << "for rest of job.";
0135                 useEcalCosmicTiming = false;
0136               } else {
0137                 ecalDigi.setSample(0,
0138                                    EcalTriggerPrimitiveSample(
0139                                        ecal_it->sample(ecal_it->sampleOfInterest() + sample - preSamples - 1).raw()));
0140               }
0141             }
0142             // else
0143             if ((!useEcalCosmicTiming) || (iphi >= 37 && iphi <= 72)) {
0144               ecalDigi.setSample(
0145                   0,
0146                   EcalTriggerPrimitiveSample(ecal_it->sample(ecal_it->sampleOfInterest() + sample - preSamples).raw()));
0147             }
0148             // push back each digi into correct "time sample" of coll
0149             ecalColl[sample].push_back(ecalDigi);
0150           }
0151         }
0152       }
0153     }
0154 
0155   if ((hcal = iEvent.getHandle(hcalTPG_)))
0156     if (hcal.isValid()) {
0157       // loop through all hcal digis
0158       for (HcalTrigPrimDigiCollection::const_iterator hcal_it = hcal->begin(); hcal_it != hcal->end(); hcal_it++) {
0159         short ieta = hcal_it->id().ieta();
0160         short iphi = hcal_it->id().iphi();
0161         // loop through time samples for each digi
0162         unsigned short digiSize = hcal_it->size();
0163         // (size of each digi must be no less than nSamples)
0164         unsigned short nSOI = (unsigned short)(hcal_it->presamples());
0165         if (digiSize < nSamples || nSOI < preSamples || ((int)(digiSize - nSOI) < (int)(nSamples - preSamples))) {
0166           unsigned short preLoopsZero = (unsigned short)(preSamples)-nSOI;
0167           // fill extra bx's at beginning with zeros
0168           for (int sample = 0; sample < preLoopsZero; sample++) {
0169             // fill first few with zeros
0170             HcalTriggerPrimitiveDigi hcalDigi(HcalTrigTowerDetId((int)ieta, (int)iphi));
0171             hcalDigi.setSize(1);
0172             hcalDigi.setPresamples(0);
0173             hcalDigi.setSample(0, HcalTriggerPrimitiveSample(0, false, 0, 0));
0174             hcalColl[sample].push_back(hcalDigi);
0175           }
0176 
0177           // loop through existing data
0178           for (int sample = preLoopsZero; sample < (preLoopsZero + digiSize); sample++) {
0179             // go through data
0180             HcalTriggerPrimitiveDigi hcalDigi(HcalTrigTowerDetId((int)ieta, (int)iphi));
0181             hcalDigi.setSize(1);
0182             hcalDigi.setPresamples(0);
0183 
0184             if (useHcalCosmicTiming && iphi >= 1 && iphi <= 36) {
0185               if (nSOI < (preSamples + 1)) {
0186                 edm::LogWarning("TooLittleData") << "HCAL data needs at least one presample "
0187                                                  << "more than the number requested "
0188                                                  << "to use hcal cosmic timing mod!  "
0189                                                  << "reverting to useHcalCosmicTiming = false "
0190                                                  << "for rest of job.";
0191                 useHcalCosmicTiming = false;
0192               } else {
0193                 hcalDigi.setSample(
0194                     0,
0195                     HcalTriggerPrimitiveSample(hcal_it->sample(hcal_it->presamples() + sample - preSamples - 1).raw()));
0196               }
0197             }
0198             // else
0199             if ((!useHcalCosmicTiming) || (iphi >= 37 && iphi <= 72)) {
0200               hcalDigi.setSample(
0201                   0, HcalTriggerPrimitiveSample(hcal_it->sample(hcal_it->presamples() + sample - preSamples).raw()));
0202             }
0203             hcalColl[sample].push_back(hcalDigi);
0204           }
0205 
0206           // fill extra bx's at end with zeros
0207           for (unsigned int sample = (preLoopsZero + digiSize); sample < nSamples; sample++) {
0208             // fill zeros!
0209             HcalTriggerPrimitiveDigi hcalDigi(HcalTrigTowerDetId((int)ieta, (int)iphi));
0210             hcalDigi.setSize(1);
0211             hcalDigi.setPresamples(0);
0212             hcalDigi.setSample(0, HcalTriggerPrimitiveSample(0, false, 0, 0));
0213             hcalColl[sample].push_back(hcalDigi);
0214           }
0215         } else {
0216           for (unsigned short sample = 0; sample < nSamples; sample++) {
0217             // put each (relevant) time sample into its own digi
0218             HcalTriggerPrimitiveDigi hcalDigi(HcalTrigTowerDetId((int)ieta, (int)iphi));
0219             hcalDigi.setSize(1);
0220             hcalDigi.setPresamples(0);
0221 
0222             if (useHcalCosmicTiming && iphi >= 1 && iphi <= 36) {
0223               if (nSOI < (preSamples + 1)) {
0224                 edm::LogWarning("TooLittleData") << "HCAL data needs at least one presample "
0225                                                  << "more than the number requested "
0226                                                  << "to use hcal cosmic timing mod!  "
0227                                                  << "reverting to useHcalCosmicTiming = false "
0228                                                  << "for rest of job.";
0229                 useHcalCosmicTiming = false;
0230               } else {
0231                 hcalDigi.setSample(
0232                     0,
0233                     HcalTriggerPrimitiveSample(hcal_it->sample(hcal_it->presamples() + sample - preSamples - 1).raw()));
0234               }
0235             }
0236             // else
0237             if ((!useHcalCosmicTiming) || (iphi >= 37 && iphi <= 72)) {
0238               if (ieta > -29 && ieta < 29)
0239                 hcalDigi.setSample(0,
0240                                    HcalTriggerPrimitiveSample(
0241                                        hcal_it->sample(hcal_it->presamples() + sample - preSamples + hbShift).raw()));
0242               if (ieta <= -29 || ieta >= 29)
0243                 hcalDigi.setSample(0,
0244                                    HcalTriggerPrimitiveSample(
0245                                        hcal_it->sample(hcal_it->presamples() + sample - preSamples + hfShift).raw()));
0246             }
0247             hcalColl[sample].push_back(hcalDigi);
0248           }
0249         }
0250       }
0251     }
0252 
0253   // Now put the events back to file
0254 
0255   for (int i = 0; i < preSamples; ++i) {
0256     char ecal_label[200];
0257     char hcal_label[200];
0258 
0259     sprintf(ecal_label, "ECALBxminus%d", preSamples - i);
0260     sprintf(hcal_label, "HCALBxminus%d", preSamples - i);
0261 
0262     std::unique_ptr<EcalTrigPrimDigiCollection> ecalIn(new EcalTrigPrimDigiCollection);
0263     std::unique_ptr<HcalTrigPrimDigiCollection> hcalIn(new HcalTrigPrimDigiCollection);
0264     for (unsigned int j = 0; j < ecalColl[i].size(); ++j) {
0265       ecalIn->push_back((ecalColl[i])[j]);
0266     }
0267     for (unsigned int j = 0; j < hcalColl[i].size(); ++j)
0268       hcalIn->push_back((hcalColl[i])[j]);
0269 
0270     iEvent.put(std::move(ecalIn), ecal_label);
0271     iEvent.put(std::move(hcalIn), hcal_label);
0272   }
0273 
0274   std::unique_ptr<EcalTrigPrimDigiCollection> ecal0(new EcalTrigPrimDigiCollection);
0275   std::unique_ptr<HcalTrigPrimDigiCollection> hcal0(new HcalTrigPrimDigiCollection);
0276   for (unsigned int j = 0; j < ecalColl[preSamples].size(); ++j)
0277     ecal0->push_back((ecalColl[preSamples])[j]);
0278   for (unsigned int j = 0; j < hcalColl[preSamples].size(); ++j)
0279     hcal0->push_back((hcalColl[preSamples])[j]);
0280 
0281   iEvent.put(std::move(ecal0), "ECALBx0");
0282   iEvent.put(std::move(hcal0), "HCALBx0");
0283 
0284   for (int i = preSamples + 1; i < preSamples + postSamples + 1; ++i) {
0285     char ecal_label[200];
0286     char hcal_label[200];
0287 
0288     sprintf(ecal_label, "ECALBxplus%d", i - preSamples);
0289     sprintf(hcal_label, "HCALBxplus%d", i - preSamples);
0290 
0291     std::unique_ptr<EcalTrigPrimDigiCollection> ecalIn2(new EcalTrigPrimDigiCollection);
0292     std::unique_ptr<HcalTrigPrimDigiCollection> hcalIn2(new HcalTrigPrimDigiCollection);
0293 
0294     for (unsigned int j = 0; j < ecalColl[i].size(); ++j)
0295       ecalIn2->push_back((ecalColl[i])[j]);
0296 
0297     for (unsigned int j = 0; j < hcalColl[i].size(); ++j)
0298       hcalIn2->push_back((hcalColl[i])[j]);
0299 
0300     iEvent.put(std::move(ecalIn2), ecal_label);
0301     iEvent.put(std::move(hcalIn2), hcal_label);
0302   }
0303 }