Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2023-05-20 02:14:12

0001 #include "DataFormats/L1THGCal/interface/HGCalMulticluster.h"
0002 #include "Geometry/Records/interface/CaloGeometryRecord.h"
0003 #include "L1Trigger/L1THGCal/interface/HGCalTriggerGeometryBase.h"
0004 #include "L1Trigger/L1THGCalUtilities/interface/HGCalTriggerNtupleBase.h"
0005 #include "L1Trigger/L1THGCal/interface/backend/HGCalTriggerClusterIdentificationBase.h"
0006 
0007 class HGCalTriggerNtupleHGCMulticlusters : public HGCalTriggerNtupleBase {
0008 public:
0009   HGCalTriggerNtupleHGCMulticlusters(const edm::ParameterSet& conf);
0010   ~HGCalTriggerNtupleHGCMulticlusters() override{};
0011   void initialize(TTree&, const edm::ParameterSet&, edm::ConsumesCollector&&) final;
0012   void fill(const edm::Event& e, const HGCalTriggerNtupleEventSetup& es) final;
0013 
0014 private:
0015   void clear() final;
0016 
0017   edm::EDGetToken multiclusters_token_;
0018 
0019   bool fill_layer_info_;
0020   bool fill_interpretation_info_;
0021 
0022   std::unique_ptr<HGCalTriggerClusterIdentificationBase> id_;
0023 
0024   int cl3d_n_;
0025   std::vector<uint32_t> cl3d_id_;
0026   std::vector<float> cl3d_pt_;
0027   std::vector<float> cl3d_energy_;
0028   std::vector<float> cl3d_eta_;
0029   std::vector<float> cl3d_phi_;
0030   std::vector<int> cl3d_clusters_n_;
0031   std::vector<std::vector<uint32_t>> cl3d_clusters_id_;
0032   std::vector<std::vector<float>> cl3d_layer_pt_;
0033   // cluster shower shapes
0034   std::vector<int> cl3d_showerlength_;
0035   std::vector<int> cl3d_coreshowerlength_;
0036   std::vector<int> cl3d_firstlayer_;
0037   std::vector<int> cl3d_maxlayer_;
0038   std::vector<float> cl3d_seetot_;
0039   std::vector<float> cl3d_seemax_;
0040   std::vector<float> cl3d_spptot_;
0041   std::vector<float> cl3d_sppmax_;
0042   std::vector<float> cl3d_szz_;
0043   std::vector<float> cl3d_srrtot_;
0044   std::vector<float> cl3d_srrmax_;
0045   std::vector<float> cl3d_srrmean_;
0046   std::vector<float> cl3d_emaxe_;
0047   std::vector<float> cl3d_hoe_;
0048   std::vector<float> cl3d_meanz_;
0049   std::vector<float> cl3d_layer10_;
0050   std::vector<float> cl3d_layer50_;
0051   std::vector<float> cl3d_layer90_;
0052   std::vector<float> cl3d_ntc67_;
0053   std::vector<float> cl3d_ntc90_;
0054   std::vector<float> cl3d_bdteg_;
0055   std::vector<int> cl3d_quality_;
0056   std::vector<std::vector<float>> cl3d_ipt_;
0057   std::vector<std::vector<float>> cl3d_ienergy_;
0058 };
0059 
0060 DEFINE_EDM_PLUGIN(HGCalTriggerNtupleFactory, HGCalTriggerNtupleHGCMulticlusters, "HGCalTriggerNtupleHGCMulticlusters");
0061 
0062 HGCalTriggerNtupleHGCMulticlusters::HGCalTriggerNtupleHGCMulticlusters(const edm::ParameterSet& conf)
0063     : HGCalTriggerNtupleBase(conf),
0064       fill_layer_info_(conf.getParameter<bool>("FillLayerInfo")),
0065       fill_interpretation_info_(conf.getParameter<bool>("FillInterpretationInfo")) {
0066   accessEventSetup_ = false;
0067 }
0068 
0069 void HGCalTriggerNtupleHGCMulticlusters::initialize(TTree& tree,
0070                                                     const edm::ParameterSet& conf,
0071                                                     edm::ConsumesCollector&& collector) {
0072   multiclusters_token_ =
0073       collector.consumes<l1t::HGCalMulticlusterBxCollection>(conf.getParameter<edm::InputTag>("Multiclusters"));
0074   id_ = std::unique_ptr<HGCalTriggerClusterIdentificationBase>{
0075       HGCalTriggerClusterIdentificationFactory::get()->create("HGCalTriggerClusterIdentificationBDT")};
0076   id_->initialize(conf.getParameter<edm::ParameterSet>("EGIdentification"));
0077 
0078   std::string prefix(conf.getUntrackedParameter<std::string>("Prefix", "cl3d"));
0079 
0080   std::string bname;
0081   auto withPrefix([&prefix, &bname](char const* vname) -> char const* {
0082     bname = prefix + "_" + vname;
0083     return bname.c_str();
0084   });
0085 
0086   tree.Branch(withPrefix("n"), &cl3d_n_, (prefix + "_n/I").c_str());
0087   tree.Branch(withPrefix("id"), &cl3d_id_);
0088   tree.Branch(withPrefix("pt"), &cl3d_pt_);
0089   tree.Branch(withPrefix("energy"), &cl3d_energy_);
0090   tree.Branch(withPrefix("eta"), &cl3d_eta_);
0091   tree.Branch(withPrefix("phi"), &cl3d_phi_);
0092   tree.Branch(withPrefix("clusters_n"), &cl3d_clusters_n_);
0093   tree.Branch(withPrefix("clusters_id"), &cl3d_clusters_id_);
0094   if (fill_layer_info_)
0095     tree.Branch(withPrefix("layer_pt"), &cl3d_layer_pt_);
0096   tree.Branch(withPrefix("showerlength"), &cl3d_showerlength_);
0097   tree.Branch(withPrefix("coreshowerlength"), &cl3d_coreshowerlength_);
0098   tree.Branch(withPrefix("firstlayer"), &cl3d_firstlayer_);
0099   tree.Branch(withPrefix("maxlayer"), &cl3d_maxlayer_);
0100   tree.Branch(withPrefix("seetot"), &cl3d_seetot_);
0101   tree.Branch(withPrefix("seemax"), &cl3d_seemax_);
0102   tree.Branch(withPrefix("spptot"), &cl3d_spptot_);
0103   tree.Branch(withPrefix("sppmax"), &cl3d_sppmax_);
0104   tree.Branch(withPrefix("szz"), &cl3d_szz_);
0105   tree.Branch(withPrefix("srrtot"), &cl3d_srrtot_);
0106   tree.Branch(withPrefix("srrmax"), &cl3d_srrmax_);
0107   tree.Branch(withPrefix("srrmean"), &cl3d_srrmean_);
0108   tree.Branch(withPrefix("emaxe"), &cl3d_emaxe_);
0109   tree.Branch(withPrefix("hoe"), &cl3d_hoe_);
0110   tree.Branch(withPrefix("meanz"), &cl3d_meanz_);
0111   tree.Branch(withPrefix("layer10"), &cl3d_layer10_);
0112   tree.Branch(withPrefix("layer50"), &cl3d_layer50_);
0113   tree.Branch(withPrefix("layer90"), &cl3d_layer90_);
0114   tree.Branch(withPrefix("ntc67"), &cl3d_ntc67_);
0115   tree.Branch(withPrefix("ntc90"), &cl3d_ntc90_);
0116   tree.Branch(withPrefix("bdteg"), &cl3d_bdteg_);
0117   tree.Branch(withPrefix("quality"), &cl3d_quality_);
0118   if (fill_interpretation_info_) {
0119     tree.Branch(withPrefix("ipt"), &cl3d_ipt_);
0120     tree.Branch(withPrefix("ienergy"), &cl3d_ienergy_);
0121   }
0122 }
0123 
0124 void HGCalTriggerNtupleHGCMulticlusters::fill(const edm::Event& e, const HGCalTriggerNtupleEventSetup& es) {
0125   // retrieve clusters 3D
0126   edm::Handle<l1t::HGCalMulticlusterBxCollection> multiclusters_h;
0127   e.getByToken(multiclusters_token_, multiclusters_h);
0128   const l1t::HGCalMulticlusterBxCollection& multiclusters = *multiclusters_h;
0129 
0130   clear();
0131   for (auto cl3d_itr = multiclusters.begin(0); cl3d_itr != multiclusters.end(0); cl3d_itr++) {
0132     cl3d_n_++;
0133     cl3d_id_.emplace_back(cl3d_itr->detId());
0134     // physical values
0135     cl3d_pt_.emplace_back(cl3d_itr->pt());
0136     cl3d_energy_.emplace_back(cl3d_itr->energy());
0137     cl3d_eta_.emplace_back(cl3d_itr->eta());
0138     cl3d_phi_.emplace_back(cl3d_itr->phi());
0139     cl3d_clusters_n_.emplace_back(cl3d_itr->constituents().size());
0140     cl3d_showerlength_.emplace_back(cl3d_itr->showerLength());
0141     cl3d_coreshowerlength_.emplace_back(cl3d_itr->coreShowerLength());
0142     cl3d_firstlayer_.emplace_back(cl3d_itr->firstLayer());
0143     cl3d_maxlayer_.emplace_back(cl3d_itr->maxLayer());
0144     cl3d_seetot_.emplace_back(cl3d_itr->sigmaEtaEtaTot());
0145     cl3d_seemax_.emplace_back(cl3d_itr->sigmaEtaEtaMax());
0146     cl3d_spptot_.emplace_back(cl3d_itr->sigmaPhiPhiTot());
0147     cl3d_sppmax_.emplace_back(cl3d_itr->sigmaPhiPhiMax());
0148     cl3d_szz_.emplace_back(cl3d_itr->sigmaZZ());
0149     cl3d_srrtot_.emplace_back(cl3d_itr->sigmaRRTot());
0150     cl3d_srrmax_.emplace_back(cl3d_itr->sigmaRRMax());
0151     cl3d_srrmean_.emplace_back(cl3d_itr->sigmaRRMean());
0152     cl3d_emaxe_.emplace_back(cl3d_itr->eMax() / cl3d_itr->energy());
0153     cl3d_hoe_.emplace_back(cl3d_itr->hOverE());
0154     cl3d_meanz_.emplace_back(std::abs(cl3d_itr->zBarycenter()));
0155     cl3d_layer10_.emplace_back(cl3d_itr->layer10percent());
0156     cl3d_layer50_.emplace_back(cl3d_itr->layer50percent());
0157     cl3d_layer90_.emplace_back(cl3d_itr->layer90percent());
0158     cl3d_ntc67_.emplace_back(cl3d_itr->triggerCells67percent());
0159     cl3d_ntc90_.emplace_back(cl3d_itr->triggerCells90percent());
0160     cl3d_bdteg_.emplace_back(id_->value(*cl3d_itr));
0161     cl3d_quality_.emplace_back(cl3d_itr->hwQual());
0162     if (fill_interpretation_info_) {
0163       std::vector<float> iPts(cl3d_itr->interpretations_size());
0164       std::vector<float> iEnergies(cl3d_itr->interpretations_size());
0165       for (auto interp = cl3d_itr->interpretations_begin(); interp != cl3d_itr->interpretations_end(); ++interp) {
0166         iPts.emplace_back(cl3d_itr->iPt(*interp));
0167         iEnergies.emplace_back(cl3d_itr->iEnergy(*interp));
0168       }
0169       cl3d_ipt_.push_back(iPts);
0170       cl3d_ienergy_.push_back(iEnergies);
0171     }
0172 
0173     //Per layer cluster information
0174     if (fill_layer_info_) {
0175       const unsigned nlayers = es.geometry->lastTriggerLayer();
0176       std::vector<float> layer_pt(nlayers, 0.0);
0177       for (const auto& cl_ptr : cl3d_itr->constituents()) {
0178         unsigned layer = es.geometry->triggerLayer(cl_ptr.second->detId());
0179         layer_pt[layer] += cl_ptr.second->pt();
0180       }
0181       cl3d_layer_pt_.emplace_back(layer_pt);
0182     }
0183 
0184     // Retrieve indices of trigger cells inside cluster
0185     cl3d_clusters_id_.emplace_back(cl3d_itr->constituents().size());
0186     std::transform(cl3d_itr->constituents_begin(),
0187                    cl3d_itr->constituents_end(),
0188                    cl3d_clusters_id_.back().begin(),
0189                    [](const std::pair<uint32_t, edm::Ptr<l1t::HGCalCluster>>& id_cl) { return id_cl.second->detId(); });
0190   }
0191 }
0192 
0193 void HGCalTriggerNtupleHGCMulticlusters::clear() {
0194   cl3d_n_ = 0;
0195   cl3d_id_.clear();
0196   cl3d_pt_.clear();
0197   cl3d_energy_.clear();
0198   cl3d_eta_.clear();
0199   cl3d_phi_.clear();
0200   cl3d_clusters_n_.clear();
0201   cl3d_clusters_id_.clear();
0202   cl3d_layer_pt_.clear();
0203   cl3d_showerlength_.clear();
0204   cl3d_coreshowerlength_.clear();
0205   cl3d_firstlayer_.clear();
0206   cl3d_maxlayer_.clear();
0207   cl3d_seetot_.clear();
0208   cl3d_seemax_.clear();
0209   cl3d_spptot_.clear();
0210   cl3d_sppmax_.clear();
0211   cl3d_szz_.clear();
0212   cl3d_srrtot_.clear();
0213   cl3d_srrmax_.clear();
0214   cl3d_srrmean_.clear();
0215   cl3d_emaxe_.clear();
0216   cl3d_hoe_.clear();
0217   cl3d_meanz_.clear();
0218   cl3d_layer10_.clear();
0219   cl3d_layer50_.clear();
0220   cl3d_layer90_.clear();
0221   cl3d_ntc67_.clear();
0222   cl3d_ntc90_.clear();
0223   cl3d_bdteg_.clear();
0224   cl3d_quality_.clear();
0225   cl3d_ipt_.clear();
0226   cl3d_ienergy_.clear();
0227 }