Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2023-10-25 09:53:17

0001 // -*- C++ -*-
0002 //
0003 // Package:    HLTHcalCalibTypeFilter
0004 // Class:      HLTHcalCalibTypeFilter
0005 //
0006 /**\class HLTHcalCalibTypeFilter HLTHcalCalibTypeFilter.cc filter/HLTHcalCalibTypeFilter/src/HLTHcalCalibTypeFilter.cc
0007 
0008 Description: Filter to select HCAL abort gap events
0009 
0010 Implementation:
0011 <Notes on implementation>
0012 */
0013 //
0014 // Original Author:  Bryan DAHMES
0015 //         Created:  Tue Jan 22 13:55:00 CET 2008
0016 //
0017 //
0018 
0019 // system include files
0020 #include <string>
0021 #include <iostream>
0022 #include <memory>
0023 
0024 // CMSSW include files
0025 #include "DataFormats/FEDRawData/interface/FEDNumbering.h"
0026 #include "DataFormats/FEDRawData/interface/FEDRawData.h"
0027 #include "DataFormats/FEDRawData/interface/FEDRawDataCollection.h"
0028 #include "DataFormats/HcalDigi/interface/HcalCalibrationEventTypes.h"
0029 #include "EventFilter/HcalRawToDigi/interface/AMC13Header.h"
0030 #include "EventFilter/HcalRawToDigi/interface/HcalDCCHeader.h"
0031 #include "EventFilter/HcalRawToDigi/interface/HcalUHTRData.h"
0032 #include "FWCore/Framework/interface/Event.h"
0033 #include "FWCore/Framework/interface/global/EDFilter.h"
0034 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0035 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0036 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0037 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0038 #include "HLTrigger/HLTcore/interface/HLTFilter.h"
0039 
0040 //
0041 // class declaration
0042 //
0043 
0044 class HLTHcalCalibTypeFilter : public edm::global::EDFilter<> {
0045 public:
0046   explicit HLTHcalCalibTypeFilter(const edm::ParameterSet&);
0047   ~HLTHcalCalibTypeFilter() override = default;
0048   static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0049 
0050 private:
0051   bool filter(edm::StreamID, edm::Event&, const edm::EventSetup&) const override;
0052 
0053   // ----------member data ---------------------------
0054   const edm::EDGetTokenT<FEDRawDataCollection> inputToken_;
0055   const std::vector<int> calibTypes_;
0056 };
0057 
0058 //
0059 // constructors and destructor
0060 //
0061 HLTHcalCalibTypeFilter::HLTHcalCalibTypeFilter(const edm::ParameterSet& config)
0062     : inputToken_(consumes<FEDRawDataCollection>(config.getParameter<edm::InputTag>("InputTag"))),
0063       calibTypes_(config.getParameter<std::vector<int> >("CalibTypes")) {}
0064 
0065 void HLTHcalCalibTypeFilter::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0066   edm::ParameterSetDescription desc;
0067   desc.add<edm::InputTag>("InputTag", edm::InputTag("source"));
0068   desc.add<std::vector<int> >("CalibTypes", {1, 2, 3, 4, 5});
0069   descriptions.add("hltHcalCalibTypeFilter", desc);
0070 }
0071 
0072 //
0073 // member functions
0074 //
0075 
0076 // ------------ method called on each new Event  ------------
0077 bool HLTHcalCalibTypeFilter::filter(edm::StreamID, edm::Event& iEvent, const edm::EventSetup& iSetup) const {
0078   auto const& rawdata = iEvent.get(inputToken_);
0079 
0080   // some inits
0081   int numZeroes = 0, numPositives = 0;
0082 
0083   // loop over all HCAL FEDs
0084   for (int fed = FEDNumbering::MINHCALFEDID; fed <= FEDNumbering::MAXHCALuTCAFEDID; fed++) {
0085     // skip FEDs in between VME and uTCA
0086     if (fed > FEDNumbering::MAXHCALFEDID && fed < FEDNumbering::MINHCALuTCAFEDID)
0087       continue;
0088 
0089     // get raw data and check if there are empty feds
0090     const FEDRawData& fedData = rawdata.FEDData(fed);
0091     if (fedData.size() < 24)
0092       continue;
0093 
0094     if (fed <= FEDNumbering::MAXHCALFEDID) {
0095       // VME get event type
0096       int eventtype = ((const HcalDCCHeader*)(fedData.data()))->getCalibType();
0097       if (eventtype == 0)
0098         numZeroes++;
0099       else
0100         numPositives++;
0101     } else {
0102       // UTCA
0103       hcal::AMC13Header const* hamc13 = (hcal::AMC13Header const*)fedData.data();
0104       for (int iamc = 0; iamc < hamc13->NAMC(); iamc++) {
0105         HcalUHTRData uhtr(hamc13->AMCPayload(iamc), hamc13->AMCSize(iamc));
0106         int eventtype = uhtr.getEventType();
0107         if (eventtype == 0)
0108           numZeroes++;
0109         else
0110           numPositives++;
0111       }
0112     }
0113   }
0114 
0115   // if there are FEDs with Non-Collision event type, check what the majority is
0116   // if calibs - true
0117   // if 0s - false
0118   if (numPositives > 0) {
0119     if (numPositives > numZeroes)
0120       return true;
0121     else
0122       edm::LogWarning("HLTHcalCalibTypeFilter") << "Conflicting Calibration Types found";
0123   }
0124 
0125   // return false if there are no positives
0126   // and if the majority has 0 calib type
0127   return false;
0128 }
0129 
0130 // declare this class as a framework plugin
0131 #include "FWCore/Framework/interface/MakerMacros.h"
0132 DEFINE_FWK_MODULE(HLTHcalCalibTypeFilter);