Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:03:19

0001 /****************************************************************************
0002  * Authors:
0003  *   Jan Kašpar
0004  ****************************************************************************/
0005 
0006 #include "FWCore/Framework/interface/MakerMacros.h"
0007 #include "FWCore/Utilities/interface/Exception.h"
0008 #include "FWCore/Framework/interface/stream/EDFilter.h"
0009 #include "FWCore/Framework/interface/Event.h"
0010 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0011 #include "FWCore/Framework/interface/EventSetup.h"
0012 #include "CondTools/RunInfo/interface/LHCInfoCombined.h"
0013 
0014 //----------------------------------------------------------------------------------------------------
0015 
0016 class XangleBetaStarFilter : public edm::stream::EDFilter<> {
0017 public:
0018   explicit XangleBetaStarFilter(const edm::ParameterSet &);
0019 
0020   static void fillDescriptions(edm::ConfigurationDescriptions &descriptions);
0021 
0022 private:
0023   edm::ESGetToken<LHCInfo, LHCInfoRcd> lhcInfoToken_;
0024   edm::ESGetToken<LHCInfoPerLS, LHCInfoPerLSRcd> lhcInfoPerLSToken_;
0025   edm::ESGetToken<LHCInfoPerFill, LHCInfoPerFillRcd> lhcInfoPerFillToken_;
0026 
0027   bool useNewLHCInfo_;
0028 
0029   double xangle_min_;
0030   double xangle_max_;
0031 
0032   double beta_star_min_;
0033   double beta_star_max_;
0034 
0035   bool filter(edm::Event &, const edm::EventSetup &) override;
0036 };
0037 
0038 //----------------------------------------------------------------------------------------------------
0039 
0040 XangleBetaStarFilter::XangleBetaStarFilter(const edm::ParameterSet &iConfig)
0041     : lhcInfoToken_(
0042           esConsumes<LHCInfo, LHCInfoRcd>(edm::ESInputTag{"", iConfig.getParameter<std::string>("lhcInfoLabel")})),
0043       lhcInfoPerLSToken_(esConsumes<LHCInfoPerLS, LHCInfoPerLSRcd>(
0044           edm::ESInputTag{"", iConfig.getParameter<std::string>("lhcInfoPerLSLabel")})),
0045       lhcInfoPerFillToken_(esConsumes<LHCInfoPerFill, LHCInfoPerFillRcd>(
0046           edm::ESInputTag{"", iConfig.getParameter<std::string>("lhcInfoPerFillLabel")})),
0047 
0048       useNewLHCInfo_(iConfig.getParameter<bool>("useNewLHCInfo")),
0049 
0050       xangle_min_(iConfig.getParameter<double>("xangle_min")),
0051       xangle_max_(iConfig.getParameter<double>("xangle_max")),
0052       beta_star_min_(iConfig.getParameter<double>("beta_star_min")),
0053       beta_star_max_(iConfig.getParameter<double>("beta_star_max")) {}
0054 
0055 //----------------------------------------------------------------------------------------------------
0056 
0057 void XangleBetaStarFilter::fillDescriptions(edm::ConfigurationDescriptions &descriptions) {
0058   edm::ParameterSetDescription desc;
0059 
0060   desc.add<std::string>("lhcInfoLabel", "")->setComment("label of the LHCInfo record");
0061   desc.add<std::string>("lhcInfoPerLSLabel", "")->setComment("label of the LHCInfoPerLS record");
0062   desc.add<std::string>("lhcInfoPerFillLabel", "")->setComment("label of the LHCInfoPerFill record");
0063 
0064   desc.add<bool>("useNewLHCInfo", false)->setComment("flag whether to use new LHCInfoPer* records or old LHCInfo");
0065 
0066   desc.add<double>("xangle_min", 0.);
0067   desc.add<double>("xangle_max", 1000.);
0068 
0069   desc.add<double>("beta_star_min", 0.);
0070   desc.add<double>("beta_star_max", 1000.);
0071 
0072   descriptions.add("xangleBetaStarFilter", desc);
0073 }
0074 
0075 //----------------------------------------------------------------------------------------------------
0076 
0077 bool XangleBetaStarFilter::filter(edm::Event & /*iEvent*/, const edm::EventSetup &iSetup) {
0078   LHCInfoCombined lhcInfoCombined(iSetup, lhcInfoPerLSToken_, lhcInfoPerFillToken_, lhcInfoToken_, useNewLHCInfo_);
0079 
0080   return (xangle_min_ <= lhcInfoCombined.crossingAngle() && lhcInfoCombined.crossingAngle() < xangle_max_) &&
0081          (beta_star_min_ <= lhcInfoCombined.betaStarX && lhcInfoCombined.betaStarX < beta_star_max_);
0082 }
0083 
0084 //----------------------------------------------------------------------------------------------------
0085 
0086 DEFINE_FWK_MODULE(XangleBetaStarFilter);