Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:06:43

0001 #include "FWCore/Framework/interface/Event.h"
0002 #include "FWCore/Framework/interface/EventSetup.h"
0003 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0004 //#include "DPGAnalysis/SiStripTools/interface/APVLatency.h"
0005 //#include "DPGAnalysis/SiStripTools/interface/APVLatencyRcd.h"
0006 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0007 //#include "FWCore/Utilities/interface/Exception.h"
0008 
0009 #include "DPGAnalysis/SiStripTools/interface/EventWithHistoryFilter.h"
0010 
0011 EventWithHistoryFilter::EventWithHistoryFilter()
0012     : m_historyToken(),
0013       m_partition(),
0014       m_APVPhaseToken(),
0015       m_apvmodes(),
0016       m_dbxrange(),
0017       m_dbxrangelat(),
0018       m_bxrange(),
0019       m_bxrangelat(),
0020       m_bxcyclerange(),
0021       m_bxcyclerangelat(),
0022       m_dbxcyclerange(),
0023       m_dbxcyclerangelat(),
0024       m_dbxtrpltrange(),
0025       m_dbxgenericrange(),
0026       m_dbxgenericfirst(0),
0027       m_dbxgenericlast(1),
0028       m_noAPVPhase(true) {
0029   printConfig(edm::InputTag(), edm::InputTag());
0030 }
0031 
0032 EventWithHistoryFilter::EventWithHistoryFilter(const edm::ParameterSet& iConfig, edm::ConsumesCollector&& iC)
0033     : m_historyToken(iC.consumes<EventWithHistory>(
0034           iConfig.getUntrackedParameter<edm::InputTag>("historyProduct", edm::InputTag("consecutiveHEs")))),
0035       m_partition(iConfig.getUntrackedParameter<std::string>("partitionName", "Any")),
0036       m_APVPhaseToken(iC.consumes<APVCyclePhaseCollection>(
0037           edm::InputTag(iConfig.getUntrackedParameter<std::string>("APVPhaseLabel", "APVPhases")))),
0038       m_apvLatencyToken(iC.esConsumes<SiStripLatency, SiStripLatencyRcd>()),
0039       m_apvmodes(iConfig.getUntrackedParameter<std::vector<int> >("apvModes", std::vector<int>())),
0040       m_dbxrange(iConfig.getUntrackedParameter<std::vector<int> >("dbxRange", std::vector<int>())),
0041       m_dbxrangelat(iConfig.getUntrackedParameter<std::vector<int> >("dbxRangeLtcyAware", std::vector<int>())),
0042       m_bxrange(iConfig.getUntrackedParameter<std::vector<int> >("absBXRange", std::vector<int>())),
0043       m_bxrangelat(iConfig.getUntrackedParameter<std::vector<int> >("absBXRangeLtcyAware", std::vector<int>())),
0044       m_bxcyclerange(iConfig.getUntrackedParameter<std::vector<int> >("absBXInCycleRange", std::vector<int>())),
0045       m_bxcyclerangelat(
0046           iConfig.getUntrackedParameter<std::vector<int> >("absBXInCycleRangeLtcyAware", std::vector<int>())),
0047       m_dbxcyclerange(iConfig.getUntrackedParameter<std::vector<int> >("dbxInCycleRange", std::vector<int>())),
0048       m_dbxcyclerangelat(
0049           iConfig.getUntrackedParameter<std::vector<int> >("dbxInCycleRangeLtcyAware", std::vector<int>())),
0050       m_dbxtrpltrange(iConfig.getUntrackedParameter<std::vector<int> >("dbxTripletRange", std::vector<int>())),
0051       m_dbxgenericrange(iConfig.getUntrackedParameter<std::vector<int> >("dbxGenericRange", std::vector<int>())),
0052       m_dbxgenericfirst(iConfig.getUntrackedParameter<unsigned int>("dbxGenericFirst", 0)),
0053       m_dbxgenericlast(iConfig.getUntrackedParameter<unsigned int>("dbxGenericLast", 1))
0054 
0055 {
0056   m_noAPVPhase = isAPVPhaseNotNeeded();
0057   printConfig(iConfig.getUntrackedParameter<edm::InputTag>("historyProduct", edm::InputTag("consecutiveHEs")),
0058               edm::InputTag(iConfig.getUntrackedParameter<std::string>("APVPhaseLabel", "APVPhases")));
0059 }
0060 
0061 void EventWithHistoryFilter::set(const edm::ParameterSet& iConfig, edm::ConsumesCollector&& iC) {
0062   m_historyToken = iC.consumes<EventWithHistory>(
0063       iConfig.getUntrackedParameter<edm::InputTag>("historyProduct", edm::InputTag("consecutiveHEs")));
0064   m_partition = iConfig.getUntrackedParameter<std::string>("partitionName", "Any");
0065   m_APVPhaseToken = iC.consumes<APVCyclePhaseCollection>(
0066       edm::InputTag(iConfig.getUntrackedParameter<std::string>("APVPhaseLabel", "APVPhases")));
0067   m_apvLatencyToken = iC.esConsumes<SiStripLatency, SiStripLatencyRcd>();
0068   m_dbxrange = iConfig.getUntrackedParameter<std::vector<int> >("dbxRange", std::vector<int>());
0069   m_dbxrangelat = iConfig.getUntrackedParameter<std::vector<int> >("dbxRangeLtcyAware", std::vector<int>());
0070   m_bxrange = iConfig.getUntrackedParameter<std::vector<int> >("absBXRange", std::vector<int>());
0071   m_bxrangelat = iConfig.getUntrackedParameter<std::vector<int> >("absBXRangeLtcyAware", std::vector<int>());
0072   m_bxcyclerange = iConfig.getUntrackedParameter<std::vector<int> >("absBXInCycleRange", std::vector<int>());
0073   m_bxcyclerangelat =
0074       iConfig.getUntrackedParameter<std::vector<int> >("absBXInCycleRangeLtcyAware", std::vector<int>());
0075   m_dbxcyclerange = iConfig.getUntrackedParameter<std::vector<int> >("dbxInCycleRange", std::vector<int>());
0076   m_dbxcyclerangelat = iConfig.getUntrackedParameter<std::vector<int> >("dbxInCycleRangeLtcyAware", std::vector<int>());
0077   m_dbxtrpltrange = iConfig.getUntrackedParameter<std::vector<int> >("dbxTripletRange", std::vector<int>());
0078   m_dbxgenericrange = iConfig.getUntrackedParameter<std::vector<int> >("dbxGenericRange", std::vector<int>());
0079   m_dbxgenericfirst = iConfig.getUntrackedParameter<int>("dbxGenericFirst", 0);
0080   m_dbxgenericlast = iConfig.getUntrackedParameter<int>("dbxGenericLast", 1);
0081 
0082   m_noAPVPhase = isAPVPhaseNotNeeded();
0083   printConfig(iConfig.getUntrackedParameter<edm::InputTag>("historyProduct", edm::InputTag("consecutiveHEs")),
0084               edm::InputTag(iConfig.getUntrackedParameter<std::string>("APVPhaseLabel", "APVPhases")));
0085 }
0086 
0087 const bool EventWithHistoryFilter::selected(const EventWithHistory& he, const edm::EventSetup& iSetup) const {
0088   const std::vector<int> dummy;
0089   return is_selected(he, iSetup, dummy);
0090 }
0091 
0092 const bool EventWithHistoryFilter::selected(const EventWithHistory& he,
0093                                             const edm::Event& iEvent,
0094                                             const edm::EventSetup& iSetup) const {
0095   const std::vector<int> apvphases = getAPVPhase(iEvent);
0096   return is_selected(he, iSetup, apvphases);
0097 }
0098 
0099 const bool EventWithHistoryFilter::selected(const edm::Event& event, const edm::EventSetup& iSetup) const {
0100   const std::vector<int> apvphases = getAPVPhase(event);
0101 
0102   edm::Handle<EventWithHistory> hEvent;
0103   event.getByToken(m_historyToken, hEvent);
0104 
0105   return is_selected(*hEvent, iSetup, apvphases);
0106 }
0107 
0108 const bool EventWithHistoryFilter::is_selected(const EventWithHistory& he,
0109                                                const edm::EventSetup& iSetup,
0110                                                const std::vector<int>& _apvphases) const {
0111   const std::vector<int>& apvphases = _apvphases;
0112   const int latency = getAPVLatency(iSetup);
0113 
0114   bool selected = true;
0115 
0116   if (!isAPVModeNotNeeded()) {
0117     const int apvmode = getAPVMode(iSetup);
0118     bool modeok = false;
0119     for (std::vector<int>::const_iterator wantedmode = m_apvmodes.begin(); wantedmode != m_apvmodes.end();
0120          ++wantedmode) {
0121       modeok = modeok || (apvmode == *wantedmode);
0122     }
0123     if (!modeok)
0124       return false;
0125   }
0126 
0127   selected = selected && (isCutInactive(m_dbxrange) || isInRange(he.deltaBX(), m_dbxrange, he.depth() != 0));
0128 
0129   selected = selected && (isCutInactive(m_dbxrangelat) ||
0130                           isInRange(he.deltaBX() - latency, m_dbxrangelat, he.depth() != 0 && latency >= 0));
0131 
0132   selected = selected && (isCutInactive(m_bxrange) || isInRange(he.absoluteBX() % 70, m_bxrange, true));
0133 
0134   selected = selected &&
0135              (isCutInactive(m_bxrangelat) || isInRange((he.absoluteBX() - latency) % 70, m_bxrangelat, latency >= 0));
0136 
0137   // loop on all the phases and require that the cut is fulfilled for at least one of them
0138 
0139   bool phaseselected;
0140 
0141   phaseselected = isCutInactive(m_bxcyclerange);
0142   for (std::vector<int>::const_iterator phase = apvphases.begin(); phase != apvphases.end(); ++phase) {
0143     phaseselected = phaseselected || isInRange(he.absoluteBXinCycle(*phase) % 70, m_bxcyclerange, *phase >= 0);
0144   }
0145   selected = selected && phaseselected;
0146 
0147   phaseselected = isCutInactive(m_bxcyclerangelat);
0148   for (std::vector<int>::const_iterator phase = apvphases.begin(); phase != apvphases.end(); ++phase) {
0149     phaseselected =
0150         phaseselected ||
0151         isInRange((he.absoluteBXinCycle(*phase) - latency) % 70, m_bxcyclerangelat, *phase >= 0 && latency >= 0);
0152   }
0153   selected = selected && phaseselected;
0154 
0155   phaseselected = isCutInactive(m_dbxcyclerange);
0156   for (std::vector<int>::const_iterator phase = apvphases.begin(); phase != apvphases.end(); ++phase) {
0157     phaseselected =
0158         phaseselected || isInRange(he.deltaBXinCycle(*phase), m_dbxcyclerange, he.depth() != 0 && *phase >= 0);
0159   }
0160   selected = selected && phaseselected;
0161 
0162   phaseselected = isCutInactive(m_dbxcyclerangelat);
0163   for (std::vector<int>::const_iterator phase = apvphases.begin(); phase != apvphases.end(); ++phase) {
0164     phaseselected = phaseselected || isInRange(he.deltaBXinCycle(*phase) - latency,
0165                                                m_dbxcyclerangelat,
0166                                                he.depth() != 0 && *phase >= 0 && latency >= 0);
0167   }
0168   selected = selected && phaseselected;
0169 
0170   // end of phase-dependent cuts
0171 
0172   selected =
0173       selected && (isCutInactive(m_dbxtrpltrange) || isInRange(he.deltaBX(1, 2), m_dbxtrpltrange, he.depth() > 1));
0174 
0175   selected =
0176       selected &&
0177       (isCutInactive(m_dbxgenericrange) ||
0178        isInRange(he.deltaBX(m_dbxgenericfirst, m_dbxgenericlast), m_dbxgenericrange, he.depth() >= m_dbxgenericlast));
0179 
0180   return selected;
0181 }
0182 
0183 const int EventWithHistoryFilter::getAPVLatency(const edm::EventSetup& iSetup) const {
0184   if (isAPVLatencyNotNeeded())
0185     return -1;
0186 
0187   const auto& apvlat = iSetup.getData(m_apvLatencyToken);
0188   const int latency = apvlat.singleLatency() != 255 ? apvlat.singleLatency() : -1;
0189 
0190   // thrown an exception if latency value is invalid
0191   /*
0192   if(latency < 0  && !isAPVLatencyNotNeeded())
0193     throw cms::Exception("InvalidAPVLatency") << " invalid APVLatency found ";
0194   */
0195 
0196   return latency;
0197 }
0198 
0199 const int EventWithHistoryFilter::getAPVMode(const edm::EventSetup& iSetup) const {
0200   if (isAPVModeNotNeeded())
0201     return -1;
0202 
0203   const auto& apvlat = iSetup.getData(m_apvLatencyToken);
0204   int mode = -1;
0205   if (apvlat.singleReadOutMode() == 1)
0206     mode = 47;
0207   if (apvlat.singleReadOutMode() == 0)
0208     mode = 37;
0209 
0210   // thrown an exception if mode value is invalid
0211   /*
0212   if(mode < 0 && !isAPVModeNotNeeded())
0213     throw cms::Exception("InvalidAPVMode") << " invalid APVMode found ";
0214   */
0215 
0216   return mode;
0217 }
0218 
0219 const std::vector<int> EventWithHistoryFilter::getAPVPhase(const edm::Event& iEvent) const {
0220   if (m_noAPVPhase) {
0221     const std::vector<int> dummy;
0222     return dummy;
0223   }
0224 
0225   edm::Handle<APVCyclePhaseCollection> apvPhases;
0226   iEvent.getByToken(m_APVPhaseToken, apvPhases);
0227 
0228   const std::vector<int> phases = apvPhases->getPhases(m_partition);
0229 
0230   /*
0231   if(!m_noAPVPhase) {
0232     if(phases.size()==0) throw cms::Exception("NoPartitionAPVPhase")
0233       << " No APV phase for partition " << m_partition.c_str() << " : check if a proper partition has been chosen ";
0234   }
0235   */
0236 
0237   return phases;
0238 }
0239 
0240 const bool EventWithHistoryFilter::isAPVLatencyNotNeeded() const {
0241   return isCutInactive(m_bxrangelat) && isCutInactive(m_dbxrangelat) && isCutInactive(m_bxcyclerangelat) &&
0242          isCutInactive(m_dbxcyclerangelat);
0243 }
0244 
0245 const bool EventWithHistoryFilter::isAPVPhaseNotNeeded() const {
0246   return isCutInactive(m_bxcyclerange) && isCutInactive(m_dbxcyclerange) && isCutInactive(m_bxcyclerangelat) &&
0247          isCutInactive(m_dbxcyclerangelat);
0248 }
0249 
0250 const bool EventWithHistoryFilter::isAPVModeNotNeeded() const { return (m_apvmodes.empty()); }
0251 
0252 const bool EventWithHistoryFilter::isCutInactive(const std::vector<int>& range) const {
0253   return (
0254       (range.empty() || (range.size() == 1 && range[0] < 0) || (range.size() == 2 && range[0] < 0 && range[1] < 0)));
0255 }
0256 
0257 const bool EventWithHistoryFilter::isInRange(const long long bx,
0258                                              const std::vector<int>& range,
0259                                              const bool extra) const {
0260   bool cut1 = range.empty() || range[0] < 0 || (extra && bx >= range[0]);
0261   bool cut2 = range.size() < 2 || range[1] < 0 || (extra && bx <= range[1]);
0262 
0263   if (range.size() >= 2 && range[0] >= 0 && range[1] >= 0 && (range[0] > range[1])) {
0264     return cut1 || cut2;
0265   } else {
0266     return cut1 && cut2;
0267   }
0268 }
0269 
0270 void EventWithHistoryFilter::printConfig(const edm::InputTag& historyTag, const edm::InputTag& apvphaseTag) const {
0271   std::string msgcategory = "EventWithHistoryFilterConfiguration";
0272 
0273   if (!(isCutInactive(m_bxrange) && isCutInactive(m_bxrangelat) && isCutInactive(m_bxcyclerange) &&
0274         isCutInactive(m_bxcyclerangelat) && isCutInactive(m_dbxrange) && isCutInactive(m_dbxrangelat) &&
0275         isCutInactive(m_dbxcyclerange) && isCutInactive(m_dbxcyclerangelat) && isCutInactive(m_dbxtrpltrange) &&
0276         isCutInactive(m_dbxgenericrange))) {
0277     edm::LogInfo(msgcategory) << "historyProduct: " << historyTag << " APVCyclePhase: " << apvphaseTag;
0278 
0279     edm::LogVerbatim(msgcategory) << "-----------------------";
0280     edm::LogVerbatim(msgcategory) << "List of active cuts:";
0281     if (!isCutInactive(m_bxrange)) {
0282       edm::LogVerbatim(msgcategory) << "......................";
0283       if (!m_bxrange.empty())
0284         edm::LogVerbatim(msgcategory) << "absoluteBX lower limit " << m_bxrange[0];
0285       if (m_bxrange.size() >= 2)
0286         edm::LogVerbatim(msgcategory) << "absoluteBX upper limit " << m_bxrange[1];
0287       edm::LogVerbatim(msgcategory) << "......................";
0288     }
0289     if (!isCutInactive(m_bxrangelat)) {
0290       edm::LogVerbatim(msgcategory) << "......................";
0291       if (!m_bxrangelat.empty())
0292         edm::LogVerbatim(msgcategory) << "absoluteBXLtcyAware lower limit " << m_bxrangelat[0];
0293       if (m_bxrangelat.size() >= 2)
0294         edm::LogVerbatim(msgcategory) << "absoluteBXLtcyAware upper limit " << m_bxrangelat[1];
0295       edm::LogVerbatim(msgcategory) << "......................";
0296     }
0297     if (!isCutInactive(m_bxcyclerange)) {
0298       edm::LogVerbatim(msgcategory) << "......................";
0299       edm::LogVerbatim(msgcategory) << "absoluteBXinCycle partition: " << m_partition;
0300       if (!m_bxcyclerange.empty())
0301         edm::LogVerbatim(msgcategory) << "absoluteBXinCycle lower limit " << m_bxcyclerange[0];
0302       if (m_bxcyclerange.size() >= 2)
0303         edm::LogVerbatim(msgcategory) << "absoluteBXinCycle upper limit " << m_bxcyclerange[1];
0304       edm::LogVerbatim(msgcategory) << "......................";
0305     }
0306     if (!isCutInactive(m_bxcyclerangelat)) {
0307       edm::LogVerbatim(msgcategory) << "......................";
0308       edm::LogVerbatim(msgcategory) << "absoluteBXinCycleLtcyAware partition: " << m_partition;
0309       if (!m_bxcyclerangelat.empty())
0310         edm::LogVerbatim(msgcategory) << "absoluteBXinCycleLtcyAware lower limit " << m_bxcyclerangelat[0];
0311       if (m_bxcyclerangelat.size() >= 2)
0312         edm::LogVerbatim(msgcategory) << "absoluteBXinCycleLtcyAware upper limit " << m_bxcyclerangelat[1];
0313       edm::LogVerbatim(msgcategory) << "......................";
0314     }
0315     if (!isCutInactive(m_dbxrange)) {
0316       edm::LogVerbatim(msgcategory) << "......................";
0317       if (!m_dbxrange.empty())
0318         edm::LogVerbatim(msgcategory) << "deltaBX lower limit " << m_dbxrange[0];
0319       if (m_dbxrange.size() >= 2)
0320         edm::LogVerbatim(msgcategory) << "deltaBX upper limit " << m_dbxrange[1];
0321       edm::LogVerbatim(msgcategory) << "......................";
0322     }
0323     if (!isCutInactive(m_dbxrangelat)) {
0324       edm::LogVerbatim(msgcategory) << "......................";
0325       if (!m_dbxrangelat.empty())
0326         edm::LogVerbatim(msgcategory) << "deltaBXLtcyAware lower limit " << m_dbxrangelat[0];
0327       if (m_dbxrangelat.size() >= 2)
0328         edm::LogVerbatim(msgcategory) << "deltaBXLtcyAware upper limit " << m_dbxrangelat[1];
0329       edm::LogVerbatim(msgcategory) << "......................";
0330     }
0331     if (!isCutInactive(m_dbxcyclerange)) {
0332       edm::LogVerbatim(msgcategory) << "......................";
0333       edm::LogVerbatim(msgcategory) << "deltaBXinCycle partition: " << m_partition;
0334       if (!m_dbxcyclerange.empty())
0335         edm::LogVerbatim(msgcategory) << "deltaBXinCycle lower limit " << m_dbxcyclerange[0];
0336       if (m_dbxcyclerange.size() >= 2)
0337         edm::LogVerbatim(msgcategory) << "deltaBXinCycle upper limit " << m_dbxcyclerange[1];
0338       edm::LogVerbatim(msgcategory) << "......................";
0339     }
0340     if (!isCutInactive(m_dbxcyclerangelat)) {
0341       edm::LogVerbatim(msgcategory) << "......................";
0342       edm::LogVerbatim(msgcategory) << "deltaBXinCycleLtcyAware partition: " << m_partition;
0343       if (!m_dbxcyclerangelat.empty())
0344         edm::LogVerbatim(msgcategory) << "deltaBXinCycleLtcyAware lower limit " << m_dbxcyclerangelat[0];
0345       if (m_dbxcyclerangelat.size() >= 2)
0346         edm::LogVerbatim(msgcategory) << "deltaBXinCycleLtcyAware upper limit " << m_dbxcyclerangelat[1];
0347       edm::LogVerbatim(msgcategory) << "......................";
0348     }
0349     if (!isCutInactive(m_dbxtrpltrange)) {
0350       edm::LogVerbatim(msgcategory) << "......................";
0351       if (!m_dbxtrpltrange.empty())
0352         edm::LogVerbatim(msgcategory) << "TripletIsolation lower limit " << m_dbxtrpltrange[0];
0353       if (m_dbxtrpltrange.size() >= 2)
0354         edm::LogVerbatim(msgcategory) << "TripletIsolation upper limit " << m_dbxtrpltrange[1];
0355       edm::LogVerbatim(msgcategory) << "......................";
0356     }
0357     if (!isCutInactive(m_dbxgenericrange)) {
0358       edm::LogVerbatim(msgcategory) << "......................";
0359       edm::LogVerbatim(msgcategory) << "Generic DBX computed between n-" << m_dbxgenericfirst << " and n-"
0360                                     << m_dbxgenericlast << " trigger";
0361       if (!m_dbxgenericrange.empty())
0362         edm::LogVerbatim(msgcategory) << "Generic DBX cut lower limit " << m_dbxgenericrange[0];
0363       if (m_dbxgenericrange.size() >= 2)
0364         edm::LogVerbatim(msgcategory) << "Generic DBX upper limit " << m_dbxgenericrange[1];
0365       edm::LogVerbatim(msgcategory) << "......................";
0366     }
0367     edm::LogVerbatim(msgcategory) << "-----------------------";
0368   }
0369 }