CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
SiPixelDigiToRaw.cc
Go to the documentation of this file.
1 #include "SiPixelDigiToRaw.h"
2 
5 
7 
12 
15 
18 
20 #include "TH1D.h"
21 #include "TFile.h"
22 
23 using namespace std;
24 
26  cablingTree_(0),
27  frameReverter_(0),
28  config_(pset),
29  hCPU(0), hDigi(0), theTimer(0)
30 {
31 
32  // Define EDProduct type
33  produces<FEDRawDataCollection>();
34 
35  // start the counters
36  eventCounter = 0;
37  allDigiCounter = 0;
38  allWordCounter = 0;
39 
40  // Timing
41  bool timing = config_.getUntrackedParameter<bool>("Timing",false);
42  if (timing) {
43  theTimer = new R2DTimerObserver("**** MY TIMING REPORT ***");
44  hCPU = new TH1D ("hCPU","hCPU",100,0.,0.050);
45  hDigi = new TH1D("hDigi","hDigi",50,0.,15000.);
46  }
47 }
48 
49 // -----------------------------------------------------------------------------
51  delete cablingTree_;
52  delete frameReverter_;
53 
54  if (theTimer) {
55  TFile rootFile("analysis.root", "RECREATE", "my histograms");
56  hCPU->Write();
57  hDigi->Write();
58  delete theTimer;
59  }
60 }
61 
62 // -----------------------------------------------------------------------------
63 
64 // -----------------------------------------------------------------------------
66  const edm::EventSetup& es)
67 {
68  using namespace sipixelobjects;
69  eventCounter++;
70  edm::LogInfo("SiPixelDigiToRaw") << "[SiPixelDigiToRaw::produce] "
71  << "event number: " << eventCounter;
72 
74  label = config_.getParameter<edm::InputTag>("InputLabel");
76 
79  typedef vector< edm::DetSet<PixelDigi> >::const_iterator DI;
80 
81  int digiCounter = 0;
82  for (DI di=digiCollection->begin(); di != digiCollection->end(); di++) {
83  digiCounter += (di->data).size();
84  digis[ di->id] = di->data;
85  }
86  allDigiCounter += digiCounter;
87 
88  if (recordWatcher.check( es )) {
90  es.get<SiPixelFedCablingMapRcd>().get( cablingMap );
91  fedIds = cablingMap->fedIds();
92  if (cablingTree_) delete cablingTree_; cablingTree_= cablingMap->cablingTree();
93  if (frameReverter_) delete frameReverter_; frameReverter_ = new SiPixelFrameReverter( es, cablingMap.product() );
94  }
95 
97  if (debug) LogDebug("SiPixelDigiToRaw") << cablingTree_->version();
98 
101  if (theTimer) theTimer->start();
102 
103  // create product (raw data)
104  std::auto_ptr<FEDRawDataCollection> buffers( new FEDRawDataCollection );
105 
106  const vector<const PixelFEDCabling *> fedList = cablingTree_->fedList();
107 
108  // convert data to raw
109  formatter.formatRawData( ev.id().event(), rawdata, digis );
110 
111  // pack raw data into collection
112  typedef vector<const PixelFEDCabling *>::const_iterator FI;
113  for (FI it = fedList.begin(); it != fedList.end(); it++) {
114  LogDebug("SiPixelDigiToRaw")<<" PRODUCE DATA FOR FED_id: " << (**it).id();
115  FEDRawData& fedRawData = buffers->FEDData( (**it).id() );
116  PixelDataFormatter::RawData::iterator fedbuffer = rawdata.find( (**it).id() );
117  if( fedbuffer != rawdata.end() ) fedRawData = fedbuffer->second;
118  LogDebug("SiPixelDigiToRaw")<<"size of data in fedRawData: "<<fedRawData.size();
119  }
120  allWordCounter += formatter.nWords();
121  if (debug) LogDebug("SiPixelDigiToRaw")
122 
123  << "Words/Digis this ev: "<<digiCounter<<"(fm:"<<formatter.nDigis()<<")/"
124  <<formatter.nWords()
125  <<" all: "<< allDigiCounter <<"/"<<allWordCounter;
126 
127  if (theTimer) {
128  theTimer->stop();
129  LogDebug("SiPixelDigiToRaw") << "TIMING IS: (real)" << theTimer->lastMeasurement().real() ;
130  LogDebug("SiPixelDigiToRaw") << " (Words/Digis) this ev: "
131  <<formatter.nWords()<<"/"<<formatter.nDigis() << "--- all :"<<allWordCounter<<"/"<<allDigiCounter;
132  hCPU->Fill( theTimer->lastMeasurement().real() );
133  hDigi->Fill(formatter.nDigis());
134  }
135 
136  ev.put( buffers );
137 
138 }
139 
140 // -----------------------------------------------------------------------------
141 
#define LogDebug(id)
T getParameter(std::string const &) const
EventNumber_t event() const
Definition: EventID.h:44
T getUntrackedParameter(std::string const &, T const &) const
std::map< cms_uint32_t, DetDigis > Digis
SiPixelDigiToRaw(const edm::ParameterSet &)
ctor
void passFrameReverter(const SiPixelFrameReverter *reverter)
static MessageDrop * instance()
Definition: MessageDrop.cc:65
unsigned long eventCounter
static bool debugEnabled
Definition: MessageDrop.h:105
virtual void produce(edm::Event &, const edm::EventSetup &)
get data, convert to raw event, attach again to Event
size_t size() const
Lenght of the data buffer in bytes.
Definition: FEDRawData.h:49
std::vector< unsigned int > fedIds
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:94
virtual ~SiPixelDigiToRaw()
dtor
void formatRawData(unsigned int lvl1_ID, RawData &fedRawData, const Digis &digis)
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:361
R2DTimerObserver * theTimer
virtual std::string version() const
map version
const LastMeasurement & lastMeasurement()
const T & get() const
Definition: EventSetup.h:55
T const * product() const
Definition: ESHandle.h:62
std::map< int, FEDRawData > RawData
edm::ESWatcher< SiPixelFedCablingMapRcd > recordWatcher
edm::ParameterSet config_
bool check(const edm::EventSetup &iSetup)
Definition: ESWatcher.h:59
edm::EventID id() const
Definition: EventBase.h:56
SiPixelFrameReverter * frameReverter_
dictionary rawdata
Definition: lumiPlot.py:393
SiPixelFedCablingTree * cablingTree_
std::vector< const PixelFEDCabling * > fedList() const
tuple size
Write out results.
edm::InputTag label