CMS 3D CMS Logo

myRawAna.cc
Go to the documentation of this file.
1 // myRawAna.cc
2 // Description: Access raw data, specifically event fragment size per FED
3 // Author: Jim Hirschauer
4 // Date: 28 - July - 2009
5 //
7 
11 
13 
15 #include <TROOT.h>
16 #include <TSystem.h>
17 #include <TFile.h>
18 #include <TCanvas.h>
19 #include <cmath>
20 
21 using namespace edm;
22 using namespace reco;
23 using namespace std;
24 
25 #define DEBUG 1
26 #define MAXJETS 100
27 
28 // Nothing passed from .cfg yet, but I leave the syntax here for
29 // future reference.
30 
32 
33 // ************************
34 // ************************
35 
38 
39  fedSize = fs->make<TH2F>("fedSize", "fedSize", 901, -0.5, 900.5, 20000, 0., 20000.);
40  totFedSize = fs->make<TH1F>("totFedSize", "totFedSize", 200, 0., 20000.);
41 }
42 
43 // ************************
44 // ************************
45 void myRawAna::analyze(const edm::Event& evt, const edm::EventSetup& es) {
46  // **************************************************************
47  // ** Access FED Information
48  // **************************************************************
49 
51  bool getFed = evt.getByLabel("source", theRaw);
52 
53  if (!getFed) {
54  std::cout << "fedRawData not available" << std::endl;
55  } else { // got the fed raw data
56  unsigned int totalFEDsize = 0;
57 
58  // HCAL FEDs are 700-730
59  unsigned int fedStart_ = 0;
60  unsigned int fedStop_ = 900;
61 
62  for (unsigned int i = fedStart_; i <= fedStop_; ++i) {
63  fedSize->Fill(i, theRaw->FEDData(i).size());
64  totalFEDsize += theRaw->FEDData(i).size();
65  }
66  totFedSize->Fill(totalFEDsize);
67  }
68 }
69 
70 // ***********************************
71 // ***********************************
static const std::string kSharedResource
Definition: TFileService.h:76
myRawAna(const edm::ParameterSet &)
Definition: myRawAna.cc:31
void beginJob() override
Definition: myRawAna.cc:36
size_t size() const
Lenght of the data buffer in bytes.
Definition: FEDRawData.h:45
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
const FEDRawData & FEDData(int fedid) const
retrieve data for fed
void endJob() override
Definition: myRawAna.cc:72
fixed size matrix
HLT enums.
void analyze(const edm::Event &, const edm::EventSetup &) override
Definition: myRawAna.cc:45
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:501