Back to home page

Project CMSSW displayed by LXR

 
 

    


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

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_varrr_;
0047   std::vector<float> cl3d_varzz_;
0048   std::vector<float> cl3d_varee_;
0049   std::vector<float> cl3d_varpp_;
0050   std::vector<float> cl3d_emaxe_;
0051   std::vector<float> cl3d_hoe_;
0052   std::vector<float> cl3d_meanz_;
0053   std::vector<float> cl3d_layer10_;
0054   std::vector<float> cl3d_layer50_;
0055   std::vector<float> cl3d_layer90_;
0056   std::vector<float> cl3d_first1layers_;
0057   std::vector<float> cl3d_first3layers_;
0058   std::vector<float> cl3d_first5layers_;
0059   std::vector<float> cl3d_firstHcal1layers_;
0060   std::vector<float> cl3d_firstHcal3layers_;
0061   std::vector<float> cl3d_firstHcal5layers_;
0062   std::vector<float> cl3d_last1layers_;
0063   std::vector<float> cl3d_last3layers_;
0064   std::vector<float> cl3d_last5layers_;
0065   std::vector<float> cl3d_emax1layers_;
0066   std::vector<float> cl3d_emax3layers_;
0067   std::vector<float> cl3d_emax5layers_;
0068   std::vector<float> cl3d_eot_;
0069   std::vector<int> cl3d_ebm0_;
0070   std::vector<int> cl3d_ebm1_;
0071   std::vector<int> cl3d_hbm_;
0072   std::vector<float> cl3d_ntc67_;
0073   std::vector<float> cl3d_ntc90_;
0074   std::vector<float> cl3d_bdteg_;
0075   std::vector<int> cl3d_quality_;
0076   std::vector<std::vector<float>> cl3d_ipt_;
0077   std::vector<std::vector<float>> cl3d_ienergy_;
0078 };
0079 
0080 DEFINE_EDM_PLUGIN(HGCalTriggerNtupleFactory, HGCalTriggerNtupleHGCMulticlusters, "HGCalTriggerNtupleHGCMulticlusters");
0081 
0082 HGCalTriggerNtupleHGCMulticlusters::HGCalTriggerNtupleHGCMulticlusters(const edm::ParameterSet& conf)
0083     : HGCalTriggerNtupleBase(conf),
0084       fill_layer_info_(conf.getParameter<bool>("FillLayerInfo")),
0085       fill_interpretation_info_(conf.getParameter<bool>("FillInterpretationInfo")) {
0086   accessEventSetup_ = false;
0087 }
0088 
0089 void HGCalTriggerNtupleHGCMulticlusters::initialize(TTree& tree,
0090                                                     const edm::ParameterSet& conf,
0091                                                     edm::ConsumesCollector&& collector) {
0092   multiclusters_token_ =
0093       collector.consumes<l1t::HGCalMulticlusterBxCollection>(conf.getParameter<edm::InputTag>("Multiclusters"));
0094   id_ = std::unique_ptr<HGCalTriggerClusterIdentificationBase>{
0095       HGCalTriggerClusterIdentificationFactory::get()->create("HGCalTriggerClusterIdentificationBDT")};
0096   id_->initialize(conf.getParameter<edm::ParameterSet>("EGIdentification"));
0097 
0098   std::string prefix(conf.getUntrackedParameter<std::string>("Prefix", "cl3d"));
0099 
0100   std::string bname;
0101   auto withPrefix([&prefix, &bname](char const* vname) -> char const* {
0102     bname = prefix + "_" + vname;
0103     return bname.c_str();
0104   });
0105 
0106   tree.Branch(withPrefix("n"), &cl3d_n_, (prefix + "_n/I").c_str());
0107   tree.Branch(withPrefix("id"), &cl3d_id_);
0108   tree.Branch(withPrefix("pt"), &cl3d_pt_);
0109   tree.Branch(withPrefix("energy"), &cl3d_energy_);
0110   tree.Branch(withPrefix("eta"), &cl3d_eta_);
0111   tree.Branch(withPrefix("phi"), &cl3d_phi_);
0112   tree.Branch(withPrefix("clusters_n"), &cl3d_clusters_n_);
0113   tree.Branch(withPrefix("clusters_id"), &cl3d_clusters_id_);
0114   if (fill_layer_info_)
0115     tree.Branch(withPrefix("layer_pt"), &cl3d_layer_pt_);
0116   tree.Branch(withPrefix("showerlength"), &cl3d_showerlength_);
0117   tree.Branch(withPrefix("coreshowerlength"), &cl3d_coreshowerlength_);
0118   tree.Branch(withPrefix("firstlayer"), &cl3d_firstlayer_);
0119   tree.Branch(withPrefix("maxlayer"), &cl3d_maxlayer_);
0120   tree.Branch(withPrefix("seetot"), &cl3d_seetot_);
0121   tree.Branch(withPrefix("seemax"), &cl3d_seemax_);
0122   tree.Branch(withPrefix("spptot"), &cl3d_spptot_);
0123   tree.Branch(withPrefix("sppmax"), &cl3d_sppmax_);
0124   tree.Branch(withPrefix("szz"), &cl3d_szz_);
0125   tree.Branch(withPrefix("srrtot"), &cl3d_srrtot_);
0126   tree.Branch(withPrefix("srrmax"), &cl3d_srrmax_);
0127   tree.Branch(withPrefix("srrmean"), &cl3d_srrmean_);
0128   tree.Branch(withPrefix("varrr"), &cl3d_varrr_);
0129   tree.Branch(withPrefix("varzz"), &cl3d_varzz_);
0130   tree.Branch(withPrefix("varee"), &cl3d_varee_);
0131   tree.Branch(withPrefix("varpp"), &cl3d_varpp_);
0132   tree.Branch(withPrefix("emaxe"), &cl3d_emaxe_);
0133   tree.Branch(withPrefix("hoe"), &cl3d_hoe_);
0134   tree.Branch(withPrefix("meanz"), &cl3d_meanz_);
0135   tree.Branch(withPrefix("layer10"), &cl3d_layer10_);
0136   tree.Branch(withPrefix("layer50"), &cl3d_layer50_);
0137   tree.Branch(withPrefix("layer90"), &cl3d_layer90_);
0138   tree.Branch(withPrefix("first1layers"), &cl3d_first1layers_);
0139   tree.Branch(withPrefix("first3layers"), &cl3d_first3layers_);
0140   tree.Branch(withPrefix("first5layers"), &cl3d_first5layers_);
0141   tree.Branch(withPrefix("firstHcal1layers"), &cl3d_firstHcal1layers_);
0142   tree.Branch(withPrefix("firstHcal3layers"), &cl3d_firstHcal3layers_);
0143   tree.Branch(withPrefix("firstHcal5layers"), &cl3d_firstHcal5layers_);
0144   tree.Branch(withPrefix("last1layers"), &cl3d_last1layers_);
0145   tree.Branch(withPrefix("last3layers"), &cl3d_last3layers_);
0146   tree.Branch(withPrefix("last5layers"), &cl3d_last5layers_);
0147   tree.Branch(withPrefix("emax1layers"), &cl3d_emax1layers_);
0148   tree.Branch(withPrefix("emax3layers"), &cl3d_emax3layers_);
0149   tree.Branch(withPrefix("emax5layers"), &cl3d_emax5layers_);
0150   tree.Branch(withPrefix("eot"), &cl3d_eot_);
0151   tree.Branch(withPrefix("ebm0"), &cl3d_ebm0_);
0152   tree.Branch(withPrefix("ebm1"), &cl3d_ebm1_);
0153   tree.Branch(withPrefix("hbm"), &cl3d_hbm_);
0154   tree.Branch(withPrefix("ntc67"), &cl3d_ntc67_);
0155   tree.Branch(withPrefix("ntc90"), &cl3d_ntc90_);
0156   tree.Branch(withPrefix("bdteg"), &cl3d_bdteg_);
0157   tree.Branch(withPrefix("quality"), &cl3d_quality_);
0158   if (fill_interpretation_info_) {
0159     tree.Branch(withPrefix("ipt"), &cl3d_ipt_);
0160     tree.Branch(withPrefix("ienergy"), &cl3d_ienergy_);
0161   }
0162 }
0163 
0164 void HGCalTriggerNtupleHGCMulticlusters::fill(const edm::Event& e, const HGCalTriggerNtupleEventSetup& es) {
0165   // retrieve clusters 3D
0166   edm::Handle<l1t::HGCalMulticlusterBxCollection> multiclusters_h;
0167   e.getByToken(multiclusters_token_, multiclusters_h);
0168   const l1t::HGCalMulticlusterBxCollection& multiclusters = *multiclusters_h;
0169 
0170   clear();
0171   for (auto cl3d_itr = multiclusters.begin(0); cl3d_itr != multiclusters.end(0); cl3d_itr++) {
0172     cl3d_n_++;
0173     cl3d_id_.emplace_back(cl3d_itr->detId());
0174     // physical values
0175     cl3d_pt_.emplace_back(cl3d_itr->pt());
0176     cl3d_energy_.emplace_back(cl3d_itr->energy());
0177     cl3d_eta_.emplace_back(cl3d_itr->eta());
0178     cl3d_phi_.emplace_back(cl3d_itr->phi());
0179     cl3d_clusters_n_.emplace_back(cl3d_itr->constituents().size());
0180     cl3d_showerlength_.emplace_back(cl3d_itr->showerLength());
0181     cl3d_coreshowerlength_.emplace_back(cl3d_itr->coreShowerLength());
0182     cl3d_firstlayer_.emplace_back(cl3d_itr->firstLayer());
0183     cl3d_maxlayer_.emplace_back(cl3d_itr->maxLayer());
0184     cl3d_seetot_.emplace_back(cl3d_itr->sigmaEtaEtaTot());
0185     cl3d_seemax_.emplace_back(cl3d_itr->sigmaEtaEtaMax());
0186     cl3d_spptot_.emplace_back(cl3d_itr->sigmaPhiPhiTot());
0187     cl3d_sppmax_.emplace_back(cl3d_itr->sigmaPhiPhiMax());
0188     cl3d_szz_.emplace_back(cl3d_itr->sigmaZZ());
0189     cl3d_srrtot_.emplace_back(cl3d_itr->sigmaRRTot());
0190     cl3d_srrmax_.emplace_back(cl3d_itr->sigmaRRMax());
0191     cl3d_srrmean_.emplace_back(cl3d_itr->sigmaRRMean());
0192     cl3d_varrr_.emplace_back(cl3d_itr->varRR());
0193     cl3d_varzz_.emplace_back(cl3d_itr->varZZ());
0194     cl3d_varee_.emplace_back(cl3d_itr->varEtaEta());
0195     cl3d_varpp_.emplace_back(cl3d_itr->varPhiPhi());
0196     cl3d_emaxe_.emplace_back(cl3d_itr->eMax() / cl3d_itr->energy());
0197     cl3d_hoe_.emplace_back(cl3d_itr->hOverE());
0198     cl3d_meanz_.emplace_back(std::abs(cl3d_itr->zBarycenter()));
0199     cl3d_layer10_.emplace_back(cl3d_itr->layer10percent());
0200     cl3d_layer50_.emplace_back(cl3d_itr->layer50percent());
0201     cl3d_layer90_.emplace_back(cl3d_itr->layer90percent());
0202     cl3d_first1layers_.emplace_back(cl3d_itr->first1layers());
0203     cl3d_first3layers_.emplace_back(cl3d_itr->first3layers());
0204     cl3d_first5layers_.emplace_back(cl3d_itr->first5layers());
0205     cl3d_firstHcal1layers_.emplace_back(cl3d_itr->firstHcal1layers());
0206     cl3d_firstHcal3layers_.emplace_back(cl3d_itr->firstHcal3layers());
0207     cl3d_firstHcal5layers_.emplace_back(cl3d_itr->firstHcal5layers());
0208     cl3d_last1layers_.emplace_back(cl3d_itr->last1layers());
0209     cl3d_last3layers_.emplace_back(cl3d_itr->last3layers());
0210     cl3d_last5layers_.emplace_back(cl3d_itr->last5layers());
0211     cl3d_emax1layers_.emplace_back(cl3d_itr->emax1layers());
0212     cl3d_emax3layers_.emplace_back(cl3d_itr->emax3layers());
0213     cl3d_emax5layers_.emplace_back(cl3d_itr->emax5layers());
0214     cl3d_eot_.emplace_back(cl3d_itr->eot());
0215     cl3d_ebm0_.emplace_back(cl3d_itr->ebm0());
0216     cl3d_ebm1_.emplace_back(cl3d_itr->ebm1());
0217     cl3d_hbm_.emplace_back(cl3d_itr->hbm());
0218     cl3d_ntc67_.emplace_back(cl3d_itr->triggerCells67percent());
0219     cl3d_ntc90_.emplace_back(cl3d_itr->triggerCells90percent());
0220     cl3d_bdteg_.emplace_back(id_->value(*cl3d_itr));
0221     cl3d_quality_.emplace_back(cl3d_itr->hwQual());
0222     if (fill_interpretation_info_) {
0223       std::vector<float> iPts(cl3d_itr->interpretations_size());
0224       std::vector<float> iEnergies(cl3d_itr->interpretations_size());
0225       for (auto interp = cl3d_itr->interpretations_begin(); interp != cl3d_itr->interpretations_end(); ++interp) {
0226         iPts.emplace_back(cl3d_itr->iPt(*interp));
0227         iEnergies.emplace_back(cl3d_itr->iEnergy(*interp));
0228       }
0229       cl3d_ipt_.push_back(iPts);
0230       cl3d_ienergy_.push_back(iEnergies);
0231     }
0232 
0233     //Per layer cluster information
0234     if (fill_layer_info_) {
0235       const unsigned nlayers = es.geometry->lastTriggerLayer();
0236       std::vector<float> layer_pt(nlayers, 0.0);
0237       for (const auto& cl_ptr : cl3d_itr->constituents()) {
0238         unsigned layer = es.geometry->triggerLayer(cl_ptr.second->detId());
0239         layer_pt[layer] += cl_ptr.second->pt();
0240       }
0241       cl3d_layer_pt_.emplace_back(layer_pt);
0242     }
0243 
0244     // Retrieve indices of trigger cells inside cluster
0245     cl3d_clusters_id_.emplace_back(cl3d_itr->constituents().size());
0246     std::transform(cl3d_itr->constituents_begin(),
0247                    cl3d_itr->constituents_end(),
0248                    cl3d_clusters_id_.back().begin(),
0249                    [](const std::pair<uint32_t, edm::Ptr<l1t::HGCalCluster>>& id_cl) { return id_cl.second->detId(); });
0250   }
0251 }
0252 
0253 void HGCalTriggerNtupleHGCMulticlusters::clear() {
0254   cl3d_n_ = 0;
0255   cl3d_id_.clear();
0256   cl3d_pt_.clear();
0257   cl3d_energy_.clear();
0258   cl3d_eta_.clear();
0259   cl3d_phi_.clear();
0260   cl3d_clusters_n_.clear();
0261   cl3d_clusters_id_.clear();
0262   cl3d_layer_pt_.clear();
0263   cl3d_showerlength_.clear();
0264   cl3d_coreshowerlength_.clear();
0265   cl3d_firstlayer_.clear();
0266   cl3d_maxlayer_.clear();
0267   cl3d_seetot_.clear();
0268   cl3d_seemax_.clear();
0269   cl3d_spptot_.clear();
0270   cl3d_sppmax_.clear();
0271   cl3d_szz_.clear();
0272   cl3d_srrtot_.clear();
0273   cl3d_srrmax_.clear();
0274   cl3d_srrmean_.clear();
0275   cl3d_varrr_.clear();
0276   cl3d_varzz_.clear();
0277   cl3d_varee_.clear();
0278   cl3d_varpp_.clear();
0279   cl3d_emaxe_.clear();
0280   cl3d_hoe_.clear();
0281   cl3d_meanz_.clear();
0282   cl3d_layer10_.clear();
0283   cl3d_layer50_.clear();
0284   cl3d_layer90_.clear();
0285   cl3d_first1layers_.clear();
0286   cl3d_first3layers_.clear();
0287   cl3d_first5layers_.clear();
0288   cl3d_firstHcal1layers_.clear();
0289   cl3d_firstHcal3layers_.clear();
0290   cl3d_firstHcal5layers_.clear();
0291   cl3d_last1layers_.clear();
0292   cl3d_last3layers_.clear();
0293   cl3d_last5layers_.clear();
0294   cl3d_emax1layers_.clear();
0295   cl3d_emax3layers_.clear();
0296   cl3d_emax5layers_.clear();
0297   cl3d_eot_.clear();
0298   cl3d_ebm0_.clear();
0299   cl3d_ebm1_.clear();
0300   cl3d_hbm_.clear();
0301   cl3d_ntc67_.clear();
0302   cl3d_ntc90_.clear();
0303   cl3d_bdteg_.clear();
0304   cl3d_quality_.clear();
0305   cl3d_ipt_.clear();
0306   cl3d_ienergy_.clear();
0307 }