CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
RawDataConverter.h
Go to the documentation of this file.
7 
8 // Forward declarations
9 class TFile;
10 class TTree;
11 
13 public:
14  explicit RawDataConverter(const edm::ParameterSet&);
15  ~RawDataConverter() override;
16 
17 private:
19  void beginJob() override;
20  void beginRun(edm::Run const&, edm::EventSetup const&) override;
21  void analyze(const edm::Event&, const edm::EventSetup&) override;
22  void endJob() override;
23 
24  void fillDetectorId(void);
25  void ClearData(void);
27  const edm::Event&
28  iEvent); // Check what kind of file is being processed and get valid module and instance labels returns the type of Digis that was found
29 
30  template <class T>
31  void GetDigis(const edm::Event&); // Copy the Digis into the local container (theData)
32 
33  std::vector<std::string> theDigiModuleLabels;
34  std::vector<std::string> theProductInstanceLabels;
35 
38 
39  TFile* theOutputFile;
40  TTree* theOutputTree;
42 
43  int latency;
45  int runnumber;
46  int lumiBlock;
48 };
49 
50 // Copy the Digis into the local container (theData)
51 // Currently this has only been implemented and tested for SiStripDigis
52 // SiStripRawDigis and SiStripProcessedRawDigis will need some changes to work (this is the final goal)
53 template <class Digitype>
55  LogDebug("RawDataConverter") << "Fill ZeroSuppressed Digis into the Tree";
56 
57  // Get the DetSetVector for the SiStripDigis
58  // This is a vector with all the modules, each module containing zero or more strips with signal (Digis)
59  edm::Handle<edm::DetSetVector<Digitype> > detSetVector; // Handle for holding the DetSetVector
60  iEvent.getByLabel(CurrentModuleLabel, CurrentInstanceLabel, detSetVector);
61  if (!detSetVector.isValid())
62  throw std::runtime_error("Could not find the Digis");
63 
64  // set everything in the local container to zero
65  ClearData();
66 
67  // Fill the Digis into the Raw Data Container
68 
69  LASGlobalLoop loop; // loop helper
70  int det, ring, beam, disk, pos; // and its variables
71 
72  // loop over TEC+- (internal) modules
73  det = 0;
74  ring = 0;
75  beam = 0;
76  disk = 0;
77  do {
78  // Find the module in the DetSetVector and get a pointer (iterator) to it
80  detSetVector->find(detectorId.GetTECEntry(det, ring, beam, disk));
81 
82  if (theModule != detSetVector->end()) {
83  // loop over all the Digis in this Module
85  for (theDigi = theModule->data.begin(); theDigi != theModule->data.end(); ++theDigi) {
86  // fill the number of adc counts into the local container
87  if (theDigi->channel() < 512)
88  theData.GetTECEntry(det, ring, beam, disk).at(theDigi->channel()) = theDigi->adc();
89  }
90  }
91  } while (loop.TECLoop(det, ring, beam, disk));
92 
93  // loop TIB/TOB
94  det = 2;
95  beam = 0;
96  pos = 0; // <- set det = 2 (TIB)
97  do {
98  // Find the module in the DetSetVector and get a pointer (iterator) to it
100  detSetVector->find(detectorId.GetTIBTOBEntry(det, beam, pos));
101 
102  if (theModule != detSetVector->end()) {
103  // loop over all the Digis in this Module
104  typename edm::DetSet<Digitype>::const_iterator theDigi;
105  for (theDigi = theModule->data.begin(); theDigi != theModule->data.end(); ++theDigi) {
106  // fill the number of adc counts into the local container
107  if (theDigi->channel() < 512)
108  theData.GetTIBTOBEntry(det, beam, pos).at(theDigi->channel()) = theDigi->adc();
109  }
110  }
111  } while (loop.TIBTOBLoop(det, beam, pos));
112 
113  // loop TEC (AT)
114  det = 0;
115  beam = 0;
116  disk = 0;
117  do {
118  // Find the module in the DetSetVector and get a pointer (iterator) to it
120  detSetVector->find(detectorId.GetTEC2TECEntry(det, beam, disk));
121 
122  if (theModule != detSetVector->end()) {
123  // loop over all the Digis in this Module
124  typename edm::DetSet<Digitype>::const_iterator theDigi;
125  for (theDigi = theModule->data.begin(); theDigi != theModule->data.end(); ++theDigi) {
126  // fill the number of adc counts into the local container
127  if (theDigi->channel() < 512)
128  theData.GetTEC2TECEntry(det, beam, disk).at(theDigi->channel()) = theDigi->adc();
129  }
130  }
131  } while (loop.TEC2TECLoop(det, beam, disk));
132 }
void fillDetectorId(void)
DigiType GetValidLabels(const edm::Event &iEvent)
LASGlobalData< int > detectorId
void beginJob() override
void analyze(const edm::Event &, const edm::EventSetup &) override
int iEvent
Definition: GenABIO.cc:224
std::vector< std::string > theDigiModuleLabels
void beginRun(edm::Run const &, edm::EventSetup const &) override
bool TEC2TECLoop(int &, int &, int &) const
void GetDigis(const edm::Event &)
T & GetTIBTOBEntry(int subdetector, int beam, int tibTobPosition)
bool isValid() const
Definition: HandleBase.h:70
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:500
iterator end()
Return the off-the-end iterator.
Definition: DetSetVector.h:325
RawDataConverter(const edm::ParameterSet &)
~RawDataConverter() override
T & GetTEC2TECEntry(int subdetector, int beam, int tecDisk)
void endJob() override
T & GetTECEntry(int subdetector, int tecRing, int beam, int tecDisk)
Definition: LASGlobalData.h:84
std::string CurrentModuleLabel
LASGlobalData< std::vector< float > > theData
bool TECLoop(int &, int &, int &, int &) const
bool TIBTOBLoop(int &, int &, int &) const
iterator begin()
Return an iterator to the first DetSet.
Definition: DetSetVector.h:314
std::string CurrentInstanceLabel
collection_type::const_iterator const_iterator
Definition: DetSet.h:31
collection_type::const_iterator const_iterator
Definition: DetSetVector.h:102
std::vector< std::string > theProductInstanceLabels
Definition: Run.h:45
#define LogDebug(id)