Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-07-03 04:18:13

0001 // Authors: Olivie Franklova - olivie.abigail.franklova@cern.ch
0002 // Date: 03/2023
0003 // @file merge layer clusters which were produce by HGCalLayerClusterProducer
0004 
0005 #include "FWCore/Framework/interface/Frameworkfwd.h"
0006 #include "FWCore/Framework/interface/stream/EDProducer.h"
0007 
0008 #include "FWCore/Framework/interface/Event.h"
0009 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0010 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0011 
0012 #include "DataFormats/EgammaReco/interface/BasicCluster.h"
0013 #include "DataFormats/Common/interface/ValueMap.h"
0014 
0015 class MergeClusterProducer : public edm::stream::EDProducer<> {
0016 public:
0017   /**
0018    * @brief Constructor with parameter settings - which can be changed in  ...todo.
0019    * Constructor will set all variables by input param ps.
0020    *
0021    * @param[in] ps parametr set to set variables
0022   */
0023   MergeClusterProducer(const edm::ParameterSet &);
0024   ~MergeClusterProducer() override {}
0025   /**
0026    * @brief Method fill description which will be used in pyhton file.
0027    *
0028    * @param[out] description to be fill
0029   */
0030   static void fillDescriptions(edm::ConfigurationDescriptions &descriptions);
0031 
0032   /**
0033    * @brief Method will merge the producers and put them back to event
0034    *
0035    * @param[in, out] evt from get info and put result
0036    * @param[in] es to get event setup info
0037   */
0038   void produce(edm::Event &, const edm::EventSetup &) override;
0039 
0040 private:
0041   edm::EDGetTokenT<std::vector<reco::CaloCluster>> EEclusters_token_;
0042   edm::EDGetTokenT<std::vector<reco::CaloCluster>> HSiclusters_token_;
0043   edm::EDGetTokenT<std::vector<reco::CaloCluster>> HSciclusters_token_;
0044 
0045   std::string timeClname_;
0046   const edm::EDGetTokenT<edm::ValueMap<std::pair<float, float>>> clustersTimeEE_token_;
0047   const edm::EDGetTokenT<edm::ValueMap<std::pair<float, float>>> clustersTimeHSi_token_;
0048   const edm::EDGetTokenT<edm::ValueMap<std::pair<float, float>>> clustersTimeHSci_token_;
0049 
0050   /**
0051    * @brief method merge three vectors of reco::CaloCluster to one
0052    *
0053    * @param[out] merge the vector into which others vectors will be merge
0054    * @param[in] EE vector for Electromagnetic silicon
0055    * @param[in] HSi vector for Hardon silicon
0056    * @param[in] ESci vector for hadron scintillator
0057   */
0058   void mergeTogether(std::vector<reco::CaloCluster> &merge,
0059                      const std::vector<reco::CaloCluster> &EE,
0060                      const std::vector<reco::CaloCluster> &HSi,
0061                      const std::vector<reco::CaloCluster> &HSci);
0062 
0063   /**
0064    * @brief copy all values from vm to to
0065    *
0066    * @param[in] vm Value map with values
0067    * @param[out] to vector to will be copy value map
0068   */
0069   void addTo(std::vector<std::pair<float, float>> &to, const edm::ValueMap<std::pair<float, float>> &vm) {
0070     size_t size = vm.size();
0071     for (size_t i = 0; i < size; ++i) {
0072       to.push_back(vm.get(i));
0073     }
0074   }
0075   /**
0076    * @brief Merge value map of time for all parts of detector together  to vector times
0077    *
0078    * @param[in] evt Event to get time value maps
0079    * @param[in] size of all 3 value maps
0080    * @param[out] times vector of merged time vectors
0081   */
0082   void mergeTime(edm::Event &evt, size_t size, std::vector<std::pair<float, float>> &times) {
0083     edm::Handle<edm::ValueMap<std::pair<float, float>>> EE, HSi, HSci;
0084     // get values from all three part of detectors
0085     evt.getByToken(clustersTimeEE_token_, EE);
0086     evt.getByToken(clustersTimeHSi_token_, HSi);
0087     evt.getByToken(clustersTimeHSci_token_, HSci);
0088 
0089     times.reserve(size);
0090     addTo(times, *EE);
0091     addTo(times, *HSi);
0092     addTo(times, *HSci);
0093   }
0094   /**
0095    * @brief get info form event and then call merge
0096    *
0097    * it is used for merge and clusters and time
0098    *
0099    * @param[in] evt Event
0100    * @param[in] EE_token token for Electromagnetic silicon
0101    * @param[in] HSi_token token for Hardon silicon
0102    * @param[in] ESci_token token for hadron scintillator
0103    * @return merged result
0104   */
0105   template <typename T>
0106   void createMerge(edm::Event &evt,
0107                    const edm::EDGetTokenT<T> &EE_token,
0108                    const edm::EDGetTokenT<T> &HSi_token,
0109                    const edm::EDGetTokenT<T> &HSci_token,
0110                    T &merge) {
0111     edm::Handle<T> EE, HSi, HSci;
0112     // get values from all three part of detectors
0113     evt.getByToken(EE_token, EE);
0114     evt.getByToken(HSi_token, HSi);
0115     evt.getByToken(HSci_token, HSci);
0116     mergeTogether(merge, *EE, *HSi, *HSci);
0117   }
0118 };
0119 
0120 MergeClusterProducer::MergeClusterProducer(const edm::ParameterSet &ps)
0121     : timeClname_(ps.getParameter<std::string>("timeClname")),
0122       clustersTimeEE_token_(
0123           consumes<edm::ValueMap<std::pair<float, float>>>(ps.getParameter<edm::InputTag>("time_layerclustersEE"))),
0124       clustersTimeHSi_token_(
0125           consumes<edm::ValueMap<std::pair<float, float>>>(ps.getParameter<edm::InputTag>("time_layerclustersHSi"))),
0126       clustersTimeHSci_token_(
0127           consumes<edm::ValueMap<std::pair<float, float>>>(ps.getParameter<edm::InputTag>("time_layerclustersHSci"))) {
0128   EEclusters_token_ = consumes<std::vector<reco::CaloCluster>>(ps.getParameter<edm::InputTag>("layerClustersEE"));
0129   HSiclusters_token_ = consumes<std::vector<reco::CaloCluster>>(ps.getParameter<edm::InputTag>("layerClustersHSi"));
0130   HSciclusters_token_ = consumes<std::vector<reco::CaloCluster>>(ps.getParameter<edm::InputTag>("layerClustersHSci"));
0131 
0132   produces<std::vector<float>>("InitialLayerClustersMask");
0133   produces<std::vector<reco::BasicCluster>>();
0134   produces<std::vector<reco::BasicCluster>>("sharing");
0135   //time for layer clusters
0136   produces<edm::ValueMap<std::pair<float, float>>>(timeClname_);
0137 }
0138 
0139 void MergeClusterProducer::fillDescriptions(edm::ConfigurationDescriptions &descriptions) {
0140   // hgcalMergeLayerClusters
0141   edm::ParameterSetDescription desc;
0142   //layer clusters
0143   desc.add<edm::InputTag>("layerClustersEE", edm::InputTag("hgcalLayerClustersEE"));
0144   desc.add<edm::InputTag>("layerClustersHSi", edm::InputTag("hgcalLayerClustersHSi"));
0145   desc.add<edm::InputTag>("layerClustersHSci", edm::InputTag("hgcalLayerClustersHSci"));
0146 
0147   //time
0148   desc.add<edm::InputTag>("time_layerclustersEE", edm::InputTag("hgcalLayerClustersEE", "timeLayerCluster"));
0149   desc.add<edm::InputTag>("time_layerclustersHSi", edm::InputTag("hgcalLayerClustersHSi", "timeLayerCluster"));
0150   desc.add<edm::InputTag>("time_layerclustersHSci", edm::InputTag("hgcalLayerClustersHSci", "timeLayerCluster"));
0151 
0152   desc.add<std::string>("timeClname", "timeLayerCluster");
0153   descriptions.add("hgcalMergeLayerClusters", desc);
0154 }
0155 
0156 void MergeClusterProducer::produce(edm::Event &evt, const edm::EventSetup &es) {
0157   //merge clusters
0158   std::unique_ptr<std::vector<reco::BasicCluster>> clusters(new std::vector<reco::BasicCluster>);
0159   createMerge(evt, EEclusters_token_, HSiclusters_token_, HSciclusters_token_, *clusters);
0160   //put new clusters to event
0161   auto clusterHandle = evt.put(std::move(clusters));
0162 
0163   //create layer cluster mask
0164   std::unique_ptr<std::vector<float>> layerClustersMask(new std::vector<float>);
0165   layerClustersMask->resize(clusterHandle->size(), 1.0);
0166   //put it into event
0167   evt.put(std::move(layerClustersMask), "InitialLayerClustersMask");
0168 
0169   //time
0170   std::vector<std::pair<float, float>> times;
0171   mergeTime(evt, clusterHandle->size(), times);
0172 
0173   auto timeCl = std::make_unique<edm::ValueMap<std::pair<float, float>>>();
0174   edm::ValueMap<std::pair<float, float>>::Filler filler(*timeCl);
0175   filler.insert(clusterHandle, times.begin(), times.end());
0176   filler.fill();
0177   evt.put(std::move(timeCl), timeClname_);
0178 }
0179 
0180 void MergeClusterProducer::mergeTogether(std::vector<reco::CaloCluster> &merge,
0181                                          const std::vector<reco::CaloCluster> &EE,
0182                                          const std::vector<reco::CaloCluster> &HSi,
0183                                          const std::vector<reco::CaloCluster> &HSci) {
0184   auto clusterSize = EE.size() + HSi.size() + HSci.size();
0185   merge.reserve(clusterSize);
0186 
0187   merge.insert(merge.end(), EE.begin(), EE.end());
0188   merge.insert(merge.end(), HSi.begin(), HSi.end());
0189   merge.insert(merge.end(), HSci.begin(), HSci.end());
0190 }
0191 
0192 #include "FWCore/Framework/interface/MakerMacros.h"
0193 DEFINE_FWK_MODULE(MergeClusterProducer);