Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:20:44

0001 #include "FWCore/Framework/interface/Frameworkfwd.h"
0002 #include "FWCore/Framework/interface/stream/EDProducer.h"
0003 #include "FWCore/Framework/interface/Event.h"
0004 #include "FWCore/Framework/interface/MakerMacros.h"
0005 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0006 
0007 #include "DataFormats/L1THGCal/interface/HGCalMulticluster.h"
0008 #include "CommonTools/Utils/interface/StringCutObjectSelector.h"
0009 #include "CommonTools/Utils/interface/StringObjectFunction.h"
0010 #include "CommonTools/MVAUtils/interface/TMVAZipReader.h"
0011 
0012 #include "TMVA/Factory.h"
0013 #include "TMVA/Reader.h"
0014 
0015 namespace l1t {
0016   class HGC3DClusterTMVASelector : public edm::stream::EDProducer<> {
0017   public:
0018     explicit HGC3DClusterTMVASelector(const edm::ParameterSet &);
0019     ~HGC3DClusterTMVASelector() override {}
0020 
0021   private:
0022     class Var {
0023     public:
0024       Var(const std::string &name, const std::string &expr) : name_(name), expr_(expr) {}
0025       void declare(TMVA::Reader &r) { r.AddVariable(name_, &val_); }
0026       void fill(const l1t::HGCalMulticluster &c) { val_ = expr_(c); }
0027 
0028     private:
0029       std::string name_;
0030       StringObjectFunction<l1t::HGCalMulticluster> expr_;
0031       float val_;
0032     };
0033 
0034     edm::EDGetTokenT<l1t::HGCalMulticlusterBxCollection> src_;
0035     StringCutObjectSelector<l1t::HGCalMulticluster> preselection_;
0036     std::vector<Var> variables_;
0037     std::string method_, weightsFile_;
0038     std::unique_ptr<TMVA::Reader> reader_;
0039     StringObjectFunction<l1t::HGCalMulticluster> wp_;
0040 
0041     void produce(edm::Event &, const edm::EventSetup &) override;
0042 
0043   };  // class
0044 }  // namespace l1t
0045 
0046 l1t::HGC3DClusterTMVASelector::HGC3DClusterTMVASelector(const edm::ParameterSet &iConfig)
0047     : src_(consumes<l1t::HGCalMulticlusterBxCollection>(iConfig.getParameter<edm::InputTag>("src"))),
0048       preselection_(iConfig.getParameter<std::string>("preselection")),
0049       method_(iConfig.getParameter<std::string>("method")),
0050       weightsFile_(iConfig.getParameter<std::string>("weightsFile")),
0051       reader_(new TMVA::Reader()),
0052       wp_(iConfig.getParameter<std::string>("wp")) {
0053   // first create all the variables
0054   for (const auto &psvar : iConfig.getParameter<std::vector<edm::ParameterSet>>("variables")) {
0055     variables_.emplace_back(psvar.getParameter<std::string>("name"), psvar.getParameter<std::string>("value"));
0056   }
0057   // then declare them
0058   for (auto &var : variables_)
0059     var.declare(*reader_);
0060   // then read the weights
0061   if (weightsFile_[0] != '/' && weightsFile_[0] != '.') {
0062     weightsFile_ = edm::FileInPath(weightsFile_).fullPath();
0063   }
0064   reco::details::loadTMVAWeights(&*reader_, method_, weightsFile_);
0065   // finally, declare outputs
0066   produces<l1t::HGCalMulticlusterBxCollection>();
0067   produces<l1t::HGCalMulticlusterBxCollection>("fail");
0068 }
0069 
0070 void l1t::HGC3DClusterTMVASelector::produce(edm::Event &iEvent, const edm::EventSetup &) {
0071   std::unique_ptr<l1t::HGCalMulticlusterBxCollection> out = std::make_unique<l1t::HGCalMulticlusterBxCollection>();
0072   std::unique_ptr<l1t::HGCalMulticlusterBxCollection> fail = std::make_unique<l1t::HGCalMulticlusterBxCollection>();
0073 
0074   edm::Handle<l1t::HGCalMulticlusterBxCollection> multiclusters;
0075   iEvent.getByToken(src_, multiclusters);
0076 
0077   for (int bx = multiclusters->getFirstBX(); bx <= multiclusters->getLastBX(); ++bx) {
0078     for (auto it = multiclusters->begin(bx), ed = multiclusters->end(bx); it != ed; ++it) {
0079       const auto &c = *it;
0080       if (preselection_(c)) {
0081         for (auto &var : variables_)
0082           var.fill(c);
0083         float mvaOut = reader_->EvaluateMVA(method_);
0084         if (mvaOut > wp_(c)) {
0085           out->push_back(bx, c);
0086         } else {
0087           fail->push_back(bx, c);
0088         }
0089       }
0090     }
0091   }
0092 
0093   iEvent.put(std::move(out));
0094   iEvent.put(std::move(fail), "fail");
0095 }
0096 using l1t::HGC3DClusterTMVASelector;
0097 DEFINE_FWK_MODULE(HGC3DClusterTMVASelector);