HGCal Test Beam  03a93d6239a951948e06fb3ef8dae4cbdebfad30
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
FNALTelescopeRawToTracks.cc
Go to the documentation of this file.
1 #include "FWCore/Framework/interface/EDProducer.h"
2 #include "FWCore/Framework/interface/Event.h"
3 #include "DataFormats/Common/interface/Handle.h"
4 #include "DataFormats/FEDRawData/interface/FEDRawDataCollection.h"
5 
6 #include "FWCore/Framework/interface/EventSetup.h"
7 #include "FWCore/ParameterSet/interface/ParameterSet.h"
8 
9 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
10 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
11 
12 #include "HGCal/DataFormats/interface/HGCalTBTrackCollection.h"
13 
14 /**
15  * \class HGCal/RawToDigi/plugins/FNALTelescopeRawToTracks.h FNALTelescopeRawToTracks.h FNALTelescopeRawToTracks
16  * \brief Produces a track collection starting from FEDRawData
17  *
18  * Unpacking done according to what implemented in
19  * HGCalTBTextSource::parseAddTelescopeWords
20  */
21 
22 
23 class FNALTelescopeRawToTracks : public edm::EDProducer
24 {
25 public:
26  explicit FNALTelescopeRawToTracks(const edm::ParameterSet& ps);
27  virtual void produce(edm::Event& e, const edm::EventSetup& c);
28  virtual void beginJob();
29  void fillDescription(edm::ParameterSetDescription &desc);
30 
31 private:
32  edm::InputTag dataTag_;
33  std::vector<int> fedIds_;
34 
35 };
36 
37 FNALTelescopeRawToTracks::FNALTelescopeRawToTracks(edm::ParameterSet const& conf):
38  dataTag_(conf.getParameter<edm::InputTag>("InputLabel")),
39  fedIds_(conf.getUntrackedParameter<std::vector<int> >("fedIds"))
40 {
41  produces<HGCalTBTrackCollection>();
42  consumes<FEDRawDataCollection>(dataTag_);
43 }
44 
45 void FNALTelescopeRawToTracks::beginJob()
46 {
47 }
48 
49 void FNALTelescopeRawToTracks::produce(edm::Event& e, const edm::EventSetup& c)
50 {
51  edm::Handle<FEDRawDataCollection> rawraw;
52  e.getByLabel(dataTag_, rawraw);
53 
54  std::auto_ptr<HGCalTBTrackCollection> tracks = std::auto_ptr<HGCalTBTrackCollection>(new HGCalTBTrackCollection);
55 
56  for(auto fedId_ : fedIds_) {
57 
58  const FEDRawData& fed = rawraw->FEDData(fedId_);
59  if(fed.size() == 0) continue; // empty FEDs are allowed: for example when one FED is not in the readout for a particular run
60 
61  const float* pdata = (const float*)(fed.data());
62 
63  unsigned int size = fed.size() / HGCalTBTrack::getSizeof(); // returns the number of HGCalTBTrack
64  ///\todo missing check if the size is correct
65  for(unsigned int i = 0; i < size; ++i, pdata += HGCalTBTrack::getSizeof()) {
66  HGCalTBTrack track(pdata);
67  tracks->push_back(track);
68  }
69  }
70 
71 
72 
73 
74  // put it into the event
75  e.put(tracks);
76 }
77 
78 void FNALTelescopeRawToTracks::fillDescription(edm::ParameterSetDescription &desc)
79 {
80  desc.setComment("TEST");
81  desc.add<edm::InputTag>("InputLabel");
82  desc.addUntracked<std::vector<int> >("fedIds", std::vector<int>(99));
83 }
84 
85 #include "FWCore/PluginManager/interface/ModuleDef.h"
86 #include "FWCore/Framework/interface/MakerMacros.h"
87 
88 DEFINE_FWK_MODULE(FNALTelescopeRawToTracks);
89 
DEFINE_FWK_MODULE(Pedestals)
static size_t getSizeof(void)
Definition: HGCalTBTrack.h:28
std::vector< HGCalTBTrack > HGCalTBTrackCollection