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