Validation
SiPixelPhase1DigisV
src
SiPixelPhase1DigisV.cc
Go to the documentation of this file.
1
// -*- C++ -*-
2
//
3
// Package: SiPixelPhase1DigisV
4
// Class: SiPixelPhase1DigisV
5
//
6
7
// Original Author: Marcel Schneider
8
9
#include "
Validation/SiPixelPhase1DigisV/interface/SiPixelPhase1DigisV.h
"
10
// Additional Authors: Alexander Morton - modifying code for validation use
11
12
// C++ stuff
13
#include <iostream>
14
15
// CMSSW stuff
16
#include "
FWCore/Framework/interface/Event.h
"
17
#include "
FWCore/Framework/interface/MakerMacros.h
"
18
#include "
FWCore/ParameterSet/interface/ParameterSet.h
"
19
20
// DQM Stuff
21
#include "
DQMServices/Core/interface/DQMStore.h
"
22
23
SiPixelPhase1DigisV::SiPixelPhase1DigisV
(
const
edm::ParameterSet
&iConfig) :
SiPixelPhase1Base
(iConfig) {
24
srcToken_
= consumes<edm::DetSetVector<PixelDigi>>(iConfig.
getParameter
<
edm::InputTag
>(
"src"
));
25
}
26
27
void
SiPixelPhase1DigisV::analyze
(
const
edm::Event
&
iEvent
,
const
edm::EventSetup
&iSetup) {
28
edm::Handle<edm::DetSetVector<PixelDigi>
>
input
;
29
iEvent
.getByToken(
srcToken_
,
input
);
30
if
(!
input
.isValid())
31
return
;
32
33
edm::DetSetVector<PixelDigi>::const_iterator
it;
34
for
(it =
input
->begin(); it !=
input
->end(); ++it) {
35
for
(
PixelDigi
const
&digi : *it) {
36
histo
[
ADC
].fill((
double
)digi.adc(),
DetId
(it->detId()), &
iEvent
);
37
histo
[
NDIGIS
].fill(
DetId
(it->detId()), &
iEvent
);
// count
38
histo
[
ROW
].fill((
double
)digi.row(),
DetId
(it->detId()), &
iEvent
);
39
histo
[
COLUMN
].fill((
double
)digi.column(),
DetId
(it->detId()), &
iEvent
);
40
}
41
}
42
histo
[
NDIGIS
].executePerEventHarvesting(&
iEvent
);
43
}
44
45
DEFINE_FWK_MODULE
(
SiPixelPhase1DigisV
);
input
static const std::string input
Definition:
EdmProvDump.cc:48
SiPixelPhase1DigisV::SiPixelPhase1DigisV
SiPixelPhase1DigisV(const edm::ParameterSet &conf)
Definition:
SiPixelPhase1DigisV.cc:23
PixelDigi
Definition:
PixelDigi.h:14
DQMStore.h
HistogramManagerHolder::histo
std::vector< HistogramManager > histo
Definition:
SiPixelPhase1Base.h:37
edm::Handle
Definition:
AssociativeIterator.h:50
SiPixelPhase1Base
Definition:
SiPixelPhase1Base.h:46
DetId
Definition:
DetId.h:17
MakerMacros.h
SiPixelPhase1DigisV::ROW
Definition:
SiPixelPhase1DigisV.h:24
DEFINE_FWK_MODULE
#define DEFINE_FWK_MODULE(type)
Definition:
MakerMacros.h:16
SiPixelPhase1DigisV::NDIGIS
Definition:
SiPixelPhase1DigisV.h:23
edm::ParameterSet
Definition:
ParameterSet.h:47
Event.h
edm::DetSetVector::const_iterator
collection_type::const_iterator const_iterator
Definition:
DetSetVector.h:102
iEvent
int iEvent
Definition:
GenABIO.cc:224
SiPixelPhase1DigisV
Definition:
SiPixelPhase1DigisV.h:19
edm::EventSetup
Definition:
EventSetup.h:57
SiPixelPhase1DigisV::ADC
Definition:
SiPixelPhase1DigisV.h:22
SiPixelPhase1DigisV::srcToken_
edm::EDGetTokenT< edm::DetSetVector< PixelDigi > > srcToken_
Definition:
SiPixelPhase1DigisV.h:36
SiPixelPhase1DigisV.h
SiPixelPhase1DigisV::COLUMN
Definition:
SiPixelPhase1DigisV.h:25
edm::ParameterSet::getParameter
T getParameter(std::string const &) const
Definition:
ParameterSet.h:303
ParameterSet.h
edm::Event
Definition:
Event.h:73
edm::InputTag
Definition:
InputTag.h:15
SiPixelPhase1DigisV::analyze
void analyze(const edm::Event &, const edm::EventSetup &) override
Definition:
SiPixelPhase1DigisV.cc:27
Generated for CMSSW Reference Manual by
1.8.16