Main Page
Namespaces
Classes
Package Documentation
All
Classes
Namespaces
Files
Functions
Variables
Typedefs
Enumerations
Enumerator
Properties
Friends
Macros
Pages
src
RecoJets
JetAnalyzers
src
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
//
6
#include "
RecoJets/JetAnalyzers/interface/myRawAna.h
"
7
8
#include "
DataFormats/Common/interface/Handle.h
"
9
#include "
DataFormats/Candidate/interface/Candidate.h
"
10
#include "
FWCore/Framework/interface/Event.h
"
11
12
#include "
DataFormats/FEDRawData/interface/FEDRawDataCollection.h
"
13
14
#include "
FWCore/ParameterSet/interface/ParameterSet.h
"
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
31
myRawAna::myRawAna
(
const
ParameterSet
& cfg ) {
32
}
33
34
35
// ************************
36
// ************************
37
38
void
myRawAna::beginJob
( ) {
39
40
edm::Service<TFileService>
fs;
41
42
fedSize = fs->
make
<TH2F>(
"fedSize"
,
"fedSize"
, 901, -0.5, 900.5, 20000, 0., 20000. );
43
totFedSize = fs->
make
<TH1F>(
"totFedSize"
,
"totFedSize"
, 200, 0., 20000. );
44
45
46
}
47
48
// ************************
49
// ************************
50
void
myRawAna::analyze
(
const
edm::Event
& evt,
const
edm::EventSetup
& es ) {
51
52
// **************************************************************
53
// ** Access FED Information
54
// **************************************************************
55
56
57
edm::Handle<FEDRawDataCollection>
theRaw;
58
bool
getFed = evt.
getByLabel
(
"source"
, theRaw);
59
60
if
( ! getFed ) {
61
std::cout
<<
"fedRawData not available"
<< std::endl;
62
}
else
{
// got the fed raw data
63
unsigned
int
totalFEDsize = 0 ;
64
65
// HCAL FEDs are 700-730
66
unsigned
int
fedStart_ = 0;
67
unsigned
int
fedStop_ = 900;
68
69
for
(
unsigned
int
i
=fedStart_;
i
<=fedStop_; ++
i
) {
70
fedSize->Fill(
i
,theRaw->FEDData(
i
).size());
71
totalFEDsize += theRaw->FEDData(
i
).size() ;
72
}
73
totFedSize->Fill(totalFEDsize);
74
}
75
76
77
}
78
79
// ***********************************
80
// ***********************************
81
void
myRawAna::endJob
() {
82
83
}
84
#include "
FWCore/Framework/interface/MakerMacros.h
"
85
DEFINE_FWK_MODULE
(
myRawAna
);
i
int i
Definition:
DBlmapReader.cc:9
myRawAna::myRawAna
myRawAna(const edm::ParameterSet &)
Definition:
myRawAna.cc:31
edm::Service< TFileService >
Event.h
MakerMacros.h
myRawAna.h
FEDRawDataCollection.h
TFileService::make
T * make(const Args &...args) const
make new ROOT object
Definition:
TFileService.h:64
edm::DEFINE_FWK_MODULE
DEFINE_FWK_MODULE(HiMixingModule)
Handle.h
edm::Handle
Definition:
AssociativeIterator.h:47
dt_dqm_sourceclient_common_cff.reco
tuple reco
Definition:
dt_dqm_sourceclient_common_cff.py:105
myRawAna
Definition:
myRawAna.h:21
ParameterSet.h
Candidate.h
myRawAna::beginJob
void beginJob()
Definition:
myRawAna.cc:38
edm::EventSetup
Definition:
EventSetup.h:44
myRawAna::analyze
void analyze(const edm::Event &, const edm::EventSetup &)
Definition:
myRawAna.cc:50
edm::Event::getByLabel
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition:
Event.h:390
myRawAna::endJob
void endJob()
Definition:
myRawAna.cc:81
edm::ParameterSet
Definition:
ParameterSet.h:35
gather_cfg.cout
tuple cout
Definition:
gather_cfg.py:121
edm::Event
Definition:
Event.h:62
Generated for CMSSW Reference Manual by
1.8.5