Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-09-07 04:35:38

0001 #ifndef CondFormats_HcalObjects_HBHENegativeEFilter_h_
0002 #define CondFormats_HcalObjects_HBHENegativeEFilter_h_
0003 
0004 #include <vector>
0005 #include <utility>
0006 #include "FWCore/Utilities/interface/Exception.h"
0007 
0008 #include "boost/serialization/utility.hpp"
0009 #include "boost/serialization/access.hpp"
0010 #include "boost/serialization/split_member.hpp"
0011 
0012 #include "DataFormats/HcalDetId/interface/HcalDetId.h"
0013 #include "CondFormats/HcalObjects/interface/PiecewiseScalingPolynomial.h"
0014 #include <cstdint>
0015 
0016 class HBHENegativeEFilter {
0017 public:
0018   inline HBHENegativeEFilter() : minCharge_(0.), tFirst_(0), tLast_(0) {}
0019 
0020   // If the vector of cuts is empty, the filter will be disabled
0021   HBHENegativeEFilter(const std::vector<PiecewiseScalingPolynomial>& a1vec,
0022                       const std::vector<PiecewiseScalingPolynomial>& a2vec,
0023                       const std::vector<uint32_t>& iEtaLimits,
0024                       const std::vector<std::pair<double, double> >& cut,
0025                       double minCharge,
0026                       unsigned firstTimeSlice,
0027                       unsigned lastTimeSlice);
0028 
0029   // Does the sequence of time slices pass the filter?
0030   bool checkPassFilter(const HcalDetId& id, const double* ts, unsigned lenTS) const;
0031 
0032   // Examing various filter data elements
0033   inline const PiecewiseScalingPolynomial& getA1(const HcalDetId& id) const { return a1v_.at(getEtaIndex(id)); }
0034   inline const PiecewiseScalingPolynomial& getA2(const HcalDetId& id) const { return a2v_.at(getEtaIndex(id)); }
0035   inline const std::vector<uint32_t>& getEtaLimits() const { return iEtaLimits_; }
0036   inline const std::vector<std::pair<double, double> >& getCut() const { return cut_; }
0037   inline double getMinCharge() const { return minCharge_; }
0038   inline unsigned getFirstTimeSlice() const { return tFirst_; }
0039   inline unsigned getLastTimeSlice() const { return tLast_; }
0040   inline bool isEnabled() const { return !cut_.empty(); }
0041 
0042   // Comparison operators
0043   bool operator==(const HBHENegativeEFilter& r) const;
0044   inline bool operator!=(const HBHENegativeEFilter& r) const { return !(*this == r); }
0045 
0046 private:
0047   unsigned getEtaIndex(const HcalDetId& id) const;
0048   bool validate() const;
0049 
0050   std::vector<PiecewiseScalingPolynomial> a1v_;
0051   std::vector<PiecewiseScalingPolynomial> a2v_;
0052   std::vector<uint32_t> iEtaLimits_;
0053   std::vector<std::pair<double, double> > cut_;
0054   double minCharge_;
0055   uint32_t tFirst_;
0056   uint32_t tLast_;
0057 
0058   friend class boost::serialization::access;
0059 
0060   template <class Archive>
0061   inline void save(Archive& ar, const unsigned /* version */) const {
0062     if (!validate())
0063       throw cms::Exception("In HBHENegativeEFilter::save: invalid data");
0064     ar & a1v_ & a2v_ & iEtaLimits_ & cut_ & minCharge_ & tFirst_ & tLast_;
0065   }
0066 
0067   template <class Archive>
0068   inline void load(Archive& ar, const unsigned /* version */) {
0069     ar & a1v_ & a2v_ & iEtaLimits_ & cut_ & minCharge_ & tFirst_ & tLast_;
0070     if (!validate())
0071       throw cms::Exception("In HBHENegativeEFilter::load: invalid data");
0072   }
0073 
0074   BOOST_SERIALIZATION_SPLIT_MEMBER()
0075 };
0076 
0077 BOOST_CLASS_VERSION(HBHENegativeEFilter, 1)
0078 
0079 #endif  // CondFormats_HcalObjects_HBHENegativeEFilter_h_