File indexing completed on 2024-04-06 12:26:28
0001
0002 #include <memory>
0003 #include <iostream>
0004
0005
0006 #include "CommonTools/UtilAlgos/interface/TFileService.h"
0007 #include "CommonTools/Utils/interface/TFileDirectory.h"
0008 #include "DataFormats/Common/interface/DetSet.h"
0009 #include "DataFormats/Common/interface/DetSetVector.h"
0010 #include "DataFormats/Common/interface/DetSetVectorNew.h"
0011 #include "DataFormats/SiStripCluster/interface/SiStripApproximateCluster.h"
0012 #include "FWCore/Framework/interface/Event.h"
0013 #include "FWCore/Framework/interface/Frameworkfwd.h"
0014 #include "FWCore/Framework/interface/MakerMacros.h"
0015 #include "FWCore/Framework/interface/one/EDAnalyzer.h"
0016 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0017 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0018 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0019 #include "FWCore/ServiceRegistry/interface/Service.h"
0020
0021
0022 #include "TROOT.h"
0023 #include "TFile.h"
0024 #include "TNtuple.h"
0025 #include "TTree.h"
0026 #include "TMath.h"
0027 #include "TList.h"
0028 #include "TString.h"
0029
0030
0031
0032
0033
0034 class SiStripApproximatedClustersDump : public edm::one::EDAnalyzer<edm::one::SharedResources> {
0035 public:
0036 explicit SiStripApproximatedClustersDump(const edm::ParameterSet&);
0037 ~SiStripApproximatedClustersDump() override;
0038
0039 static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0040
0041 private:
0042 void analyze(const edm::Event&, const edm::EventSetup&) override;
0043
0044 edm::InputTag inputTagClusters;
0045 edm::EDGetTokenT<edmNew::DetSetVector<SiStripApproximateCluster>> clusterToken;
0046
0047 TTree* outNtuple;
0048 edm::Service<TFileService> fs;
0049
0050 uint32_t detId;
0051 uint16_t barycenter;
0052 uint16_t width;
0053 uint8_t avCharge;
0054 edm::EventNumber_t eventN;
0055 };
0056
0057 SiStripApproximatedClustersDump::SiStripApproximatedClustersDump(const edm::ParameterSet& conf) {
0058 inputTagClusters = conf.getParameter<edm::InputTag>("approxSiStripClustersTag");
0059 clusterToken = consumes<edmNew::DetSetVector<SiStripApproximateCluster>>(inputTagClusters);
0060
0061 usesResource("TFileService");
0062
0063 outNtuple = fs->make<TTree>("ApproxClusters", "ApproxClusters");
0064 outNtuple->Branch("event", &eventN, "event/i");
0065 outNtuple->Branch("detId", &detId, "detId/i");
0066 outNtuple->Branch("barycenter", &barycenter, "barycenter/F");
0067 outNtuple->Branch("width", &width, "width/b");
0068 outNtuple->Branch("charge", &avCharge, "charge/b");
0069 }
0070
0071 SiStripApproximatedClustersDump::~SiStripApproximatedClustersDump() = default;
0072
0073 void SiStripApproximatedClustersDump::analyze(const edm::Event& event, const edm::EventSetup& es) {
0074 edm::Handle<edmNew::DetSetVector<SiStripApproximateCluster>> clusterCollection = event.getHandle(clusterToken);
0075
0076 for (const auto& detClusters : *clusterCollection) {
0077 detId = detClusters.detId();
0078 eventN = event.id().event();
0079
0080 for (const auto& cluster : detClusters) {
0081 barycenter = cluster.barycenter();
0082 width = cluster.width();
0083 avCharge = cluster.avgCharge();
0084 outNtuple->Fill();
0085 }
0086 }
0087 }
0088
0089 void SiStripApproximatedClustersDump::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0090 edm::ParameterSetDescription desc;
0091 desc.add<edm::InputTag>("approxSiStripClustersTag", edm::InputTag("SiStripClusters2ApproxClusters"));
0092 descriptions.add("SiStripApproximatedClustersDump", desc);
0093 }
0094
0095 #include "FWCore/Framework/interface/MakerMacros.h"
0096 DEFINE_FWK_MODULE(SiStripApproximatedClustersDump);