DPGAnalysis
Skims
src
SelectHFMinBias.cc
Go to the documentation of this file.
1
// -*- C++ -*-
2
//
3
// Package: BeamSplash
4
// Class: BeamSPlash
5
//
6
//
7
// Original Author: Luca Malgeri
8
9
#include <memory>
10
#include <vector>
11
#include <map>
12
#include <set>
13
14
// user include files
15
16
#include "
FWCore/Utilities/interface/InputTag.h
"
17
#include "
FWCore/Framework/interface/Frameworkfwd.h
"
18
#include "
FWCore/Framework/interface/EDFilter.h
"
19
#include "
FWCore/Framework/interface/Event.h
"
20
#include "
FWCore/Framework/interface/MakerMacros.h
"
21
#include "
FWCore/ParameterSet/interface/ParameterSet.h
"
22
#include "
FWCore/Framework/interface/ESHandle.h
"
23
#include "
FWCore/MessageLogger/interface/MessageLogger.h
"
24
25
#include "
DPGAnalysis/Skims/interface/SelectHFMinBias.h
"
26
#include "
DataFormats/HcalDetId/interface/HcalSubdetector.h
"
27
#include "
DataFormats/CaloTowers/interface/CaloTowerCollection.h
"
28
29
using namespace
edm
;
30
using namespace
std
;
31
32
SelectHFMinBias::SelectHFMinBias
(
const
edm::ParameterSet
& iConfig) {}
33
34
SelectHFMinBias::~SelectHFMinBias
() {}
35
36
bool
SelectHFMinBias::filter
(
edm::Event
&
iEvent
,
const
edm::EventSetup
& iSetup) {
37
edm::Handle<CaloTowerCollection>
towers
;
38
iEvent
.getByLabel(
"towerMaker"
,
towers
);
39
40
int
negTowers = 0;
41
int
posTowers = 0;
42
for
(
CaloTowerCollection::const_iterator
cal =
towers
->begin(); cal !=
towers
->end(); ++cal) {
43
for
(
unsigned
int
i
= 0;
i
< cal->constituentsSize();
i
++) {
44
const
DetId
id
= cal->constituent(
i
);
45
if
(
id
.det() ==
DetId::Hcal
) {
46
HcalSubdetector
subdet = (
HcalSubdetector
(
id
.subdetId()));
47
if
(subdet ==
HcalForward
) {
48
if
(cal->energy() > 3. && cal->eta() < -3.)
49
negTowers++;
50
if
(cal->energy() > 3. && cal->eta() > 3.)
51
posTowers++;
52
}
53
}
54
}
55
}
56
if
(negTowers > 0 && posTowers > 0)
57
return
true
;
58
59
return
false
;
60
}
61
62
//define this as a plug-in
63
DEFINE_FWK_MODULE
(
SelectHFMinBias
);
mps_fire.i
i
Definition:
mps_fire.py:355
edm::SortedCollection< CaloTower >::const_iterator
std::vector< CaloTower >::const_iterator const_iterator
Definition:
SortedCollection.h:80
MessageLogger.h
ESHandle.h
edm
HLT enums.
Definition:
AlignableModifier.h:19
DetId::Hcal
Definition:
DetId.h:28
EDFilter.h
edm::Handle
Definition:
AssociativeIterator.h:50
DetId
Definition:
DetId.h:17
MakerMacros.h
DEFINE_FWK_MODULE
#define DEFINE_FWK_MODULE(type)
Definition:
MakerMacros.h:16
SelectHFMinBias::~SelectHFMinBias
~SelectHFMinBias() override
Definition:
SelectHFMinBias.cc:34
edm::ParameterSet
Definition:
ParameterSet.h:36
Event.h
SelectHFMinBias::SelectHFMinBias
SelectHFMinBias(const edm::ParameterSet &)
Definition:
SelectHFMinBias.cc:32
iEvent
int iEvent
Definition:
GenABIO.cc:224
edm::EventSetup
Definition:
EventSetup.h:57
HcalSubdetector.h
SelectHFMinBias.h
InputTag.h
CaloTowerCollection.h
HcalSubdetector
HcalSubdetector
Definition:
HcalAssistant.h:31
HcalForward
Definition:
HcalAssistant.h:36
SelectHFMinBias::filter
bool filter(edm::Event &, const edm::EventSetup &) override
Definition:
SelectHFMinBias.cc:36
HLT_2018_cff.towers
towers
Definition:
HLT_2018_cff.py:35030
std
Definition:
JetResolutionObject.h:76
Frameworkfwd.h
SelectHFMinBias
Definition:
SelectHFMinBias.h:33
ParameterSet.h
edm::Event
Definition:
Event.h:73
Generated for CMSSW Reference Manual by
1.8.16