Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:26:27

0001 
0002 
0003 #include "FWCore/Framework/interface/MakerMacros.h"
0004 #include "FWCore/Framework/interface/Frameworkfwd.h"
0005 #include "FWCore/Framework/interface/stream/EDProducer.h"
0006 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0007 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0008 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0009 #include "FWCore/Utilities/interface/InputTag.h"
0010 #include "DataFormats/SiStripCluster/interface/SiStripApproximateCluster.h"
0011 #include "DataFormats/SiStripCluster/interface/SiStripCluster.h"
0012 #include "DataFormats/Common/interface/DetSetVectorNew.h"
0013 #include "DataFormats/Common/interface/DetSetVector.h"
0014 
0015 #include <vector>
0016 #include <memory>
0017 
0018 class SiStripApprox2ApproxClusters : public edm::stream::EDProducer<> {
0019 public:
0020   explicit SiStripApprox2ApproxClusters(const edm::ParameterSet& conf);
0021   void produce(edm::Event&, const edm::EventSetup&) override;
0022 
0023   static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0024 
0025 private:
0026   edm::InputTag inputApproxClusters;
0027   uint8_t approxVersion;
0028   std::string approxVersionS;
0029   edm::EDGetTokenT<edmNew::DetSetVector<SiStripApproximateCluster>> clusterToken;
0030 };
0031 
0032 SiStripApprox2ApproxClusters::SiStripApprox2ApproxClusters(const edm::ParameterSet& conf) {
0033   inputApproxClusters = conf.getParameter<edm::InputTag>("inputApproxClusters");
0034   approxVersionS = conf.getParameter<std::string>("approxVersion");
0035 
0036   approxVersion = -1;
0037 
0038   if (approxVersionS == "ORIGINAL")
0039     approxVersion = 0;
0040   else if (approxVersionS == "FULL_WIDTH")
0041     approxVersion = 1;
0042   else if (approxVersionS == "BARY_RES_0.1")
0043     approxVersion = 2;
0044   else if (approxVersionS == "BARY_CHARGE_RES_0.1")
0045     approxVersion = 3;
0046 
0047   clusterToken = consumes<edmNew::DetSetVector<SiStripApproximateCluster>>(inputApproxClusters);
0048   produces<edmNew::DetSetVector<SiStripApproximateCluster>>();
0049 }
0050 
0051 void SiStripApprox2ApproxClusters::produce(edm::Event& event, edm::EventSetup const&) {
0052   auto result = std::make_unique<edmNew::DetSetVector<SiStripApproximateCluster>>();
0053   const auto& clusterCollection = event.get(clusterToken);
0054 
0055   for (const auto& detClusters : clusterCollection) {
0056     edmNew::DetSetVector<SiStripApproximateCluster>::FastFiller ff{*result, detClusters.id()};
0057 
0058     for (const auto& cluster : detClusters) {
0059       float barycenter = cluster.barycenter();
0060       uint8_t width = cluster.width();
0061       float avgCharge = cluster.avgCharge();
0062       bool filter = cluster.filter();
0063       bool isSaturated = cluster.isSaturated();
0064 
0065       switch (approxVersion) {
0066         case 0:  //ORIGINAL
0067           barycenter = std::round(barycenter);
0068           if (width > 0x3F)
0069             width = 0x3F;
0070           avgCharge = std::round(avgCharge);
0071           break;
0072         case 1:  //FULL_WIDTH
0073           barycenter = std::round(barycenter);
0074           avgCharge = std::round(avgCharge);
0075           break;
0076         case 2:  //BARY_RES_0.1
0077           barycenter = std::round(barycenter * 10) / 10;
0078           if (width > 0x3F)
0079             width = 0x3F;
0080           avgCharge = std::round(avgCharge);
0081           break;
0082         case 3:  //BARY_CHARGE_RES_0.1
0083           barycenter = std::round(barycenter * 10) / 10;
0084           if (width > 0x3F)
0085             width = 0x3F;
0086           avgCharge = std::round(avgCharge * 10) / 10;
0087           break;
0088       }
0089 
0090       ff.push_back(SiStripApproximateCluster(barycenter, width, avgCharge, filter, isSaturated));
0091     }
0092   }
0093 
0094   event.put(std::move(result));
0095 }
0096 
0097 void SiStripApprox2ApproxClusters::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0098   edm::ParameterSetDescription desc;
0099   desc.add<edm::InputTag>("inputApproxClusters", edm::InputTag("siStripClusters"));
0100   desc.add<std::string>("approxVersion", std::string("ORIGINAL"));
0101 
0102   descriptions.add("SiStripApprox2ApproxClusters", desc);
0103 }
0104 
0105 DEFINE_FWK_MODULE(SiStripApprox2ApproxClusters);