CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
RawDataConverter.h
Go to the documentation of this file.
5 
6 // Forward declarations
7 class TFile;
8 class TTree;
9 
11 
12  public:
13  explicit RawDataConverter(const edm::ParameterSet&);
15 
16 
17  private:
19  virtual void beginJob() ;
20  virtual void beginRun(edm::Run const &, edm::EventSetup const &) ;
21  virtual void analyze(const edm::Event&, const edm::EventSetup&);
22  virtual void endJob() ;
23 
24  void fillDetectorId( void );
25  void ClearData( void );
26  DigiType GetValidLabels( const edm::Event& iEvent ); // Check what kind of file is being processed and get valid module and instance labels returns the type of Digis that was found
27 
28  template <class T>
29  void GetDigis(const edm::Event&); // Copy the Digis into the local container (theData)
30 
31  std::vector<std::string> theDigiModuleLabels;
32  std::vector<std::string> theProductInstanceLabels;
33 
36 
37  TFile* theOutputFile;
38  TTree* theOutputTree;
40 
41  int latency;
43  int runnumber;
44  int lumiBlock;
46 
47 };
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 {
56  LogDebug("RawDataConverter") << "Fill ZeroSuppressed Digis into the Tree";
57 
58  // Get the DetSetVector for the SiStripDigis
59  // This is a vector with all the modules, each module containing zero or more strips with signal (Digis)
60  edm::Handle< edm::DetSetVector< Digitype > > detSetVector; // Handle for holding the DetSetVector
61  iEvent.getByLabel( CurrentModuleLabel , CurrentInstanceLabel , detSetVector );
62  if( ! detSetVector.isValid() ) 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; ring = 0; beam = 0; disk = 0;
74  do {
75  // Find the module in the DetSetVector and get a pointer (iterator) to it
76  typename edm::DetSetVector< Digitype >::const_iterator theModule = detSetVector->find( detectorId.GetTECEntry( det, ring, beam, disk ) );
77 
78  if ( theModule != detSetVector->end() ) {
79  // loop over all the Digis in this Module
81  for (theDigi = theModule->data.begin(); theDigi != theModule->data.end(); ++theDigi ) {
82  // fill the number of adc counts into the local container
83  if ( theDigi->channel() < 512 ) theData.GetTECEntry( det, ring, beam, disk ).at( theDigi->channel() ) = theDigi->adc();
84  }
85  }
86  } while( loop.TECLoop( det, ring, beam, disk ) );
87 
88  // loop TIB/TOB
89  det = 2; beam = 0; pos = 0; // <- set det = 2 (TIB)
90  do {
91  // Find the module in the DetSetVector and get a pointer (iterator) to it
92  typename edm::DetSetVector< Digitype >::const_iterator theModule = detSetVector->find( detectorId.GetTIBTOBEntry( det, beam, pos ) );
93 
94  if ( theModule != detSetVector->end() ) {
95  // loop over all the Digis in this Module
97  for (theDigi = theModule->data.begin(); theDigi != theModule->data.end(); ++theDigi ) {
98  // fill the number of adc counts into the local container
99  if ( theDigi->channel() < 512 ) theData.GetTIBTOBEntry( det, beam, pos ).at( theDigi->channel() ) = theDigi->adc();
100  }
101  }
102  } while( loop.TIBTOBLoop( det, beam, pos ) );
103 
104 
105  // loop TEC (AT)
106  det = 0; beam = 0; disk = 0;
107  do {
108  // Find the module in the DetSetVector and get a pointer (iterator) to it
109  typename edm::DetSetVector< Digitype >::const_iterator theModule = detSetVector->find( detectorId.GetTEC2TECEntry( det, beam, disk ) );
110 
111  if ( theModule != detSetVector->end() ) {
112  // loop over all the Digis in this Module
114  for (theDigi = theModule->data.begin(); theDigi != theModule->data.end(); ++theDigi ) {
115  // fill the number of adc counts into the local container
116  if ( theDigi->channel() < 512 ) theData.GetTEC2TECEntry( det, beam, disk ).at( theDigi->channel() ) = theDigi->adc();
117  }
118  }
119  } while( loop.TEC2TECLoop( det, beam, disk ) );
120 }
121 
#define LogDebug(id)
virtual void endJob()
void fillDetectorId(void)
DigiType GetValidLabels(const edm::Event &iEvent)
LASGlobalData< int > detectorId
int loop
CMSSW
int iEvent
Definition: GenABIO.cc:230
std::vector< std::string > theDigiModuleLabels
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:76
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:390
iterator end()
Return the off-the-end iterator.
Definition: DetSetVector.h:363
RawDataConverter(const edm::ParameterSet &)
virtual void beginJob()
T & GetTEC2TECEntry(int subdetector, int beam, int tecDisk)
T & GetTECEntry(int subdetector, int tecRing, int beam, int tecDisk)
Definition: LASGlobalData.h:91
std::string CurrentModuleLabel
virtual void analyze(const edm::Event &, const edm::EventSetup &)
virtual void beginRun(edm::Run const &, edm::EventSetup const &)
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:348
std::string CurrentInstanceLabel
collection_type::const_iterator const_iterator
Definition: DetSet.h:33
collection_type::const_iterator const_iterator
Definition: DetSetVector.h:106
std::vector< std::string > theProductInstanceLabels
Definition: Run.h:41