Main Page
Namespaces
Classes
Package Documentation
HLTrigger
special
plugins
HLTPixelAsymmetryFilter.h
Go to the documentation of this file.
1
#ifndef HLTPixelAsymmetryFilter_h
2
#define HLTPixelAsymmetryFilter_h
3
4
6
//
7
// HLTPixelAsymmetryFilter
8
//
9
// Filter definition
10
//
11
// We perform a selection on PIXEL cluster repartition
12
//
13
// This filter is primarily used to select Beamgas (aka PKAM) events
14
//
15
// An asymmetry parameter, based on the pixel clusters, is computed as follows
16
//
17
// asym1 = fpix-/(fpix- + fpix+) for beam1
18
// asym2 = fpix+/(fpix- + fpix+) for beam2
19
//
20
// with:
21
//
22
// fpix- = mean cluster charge in FPIX-
23
// fpix+ = mean cluster charge in FPIX+
24
// bpix = mean cluster charge in BarrelPIX
25
//
26
// Usually for PKAM events, cluster repartition is quite uniform and asymmetry is around 0.5
27
//
28
//
29
// More details:
30
// http://sviret.web.cern.ch/sviret/Welcome.php?n=CMS.MIB
31
//
32
// S.Viret: 12/01/2011 (viret@in2p3.fr)
33
//
35
36
#include "
HLTrigger/HLTcore/interface/HLTFilter.h
"
37
38
#include "
FWCore/Framework/interface/Event.h
"
39
#include "
FWCore/Framework/interface/EventSetup.h
"
40
#include "
FWCore/ParameterSet/interface/ParameterSet.h
"
41
#include "
FWCore/ParameterSet/interface/ConfigurationDescriptions.h
"
42
#include "
FWCore/ParameterSet/interface/ParameterSetDescription.h
"
43
#include "
FWCore/MessageLogger/interface/MessageLogger.h
"
44
#include "
DataFormats/Common/interface/Handle.h
"
45
#include "
DataFormats/SiPixelCluster/interface/SiPixelCluster.h
"
46
#include "
DataFormats/SiPixelDetId/interface/PixelSubdetector.h
"
47
#include "
DataFormats/SiPixelDetId/interface/PixelBarrelName.h
"
48
#include "
DataFormats/SiPixelDetId/interface/PixelEndcapName.h
"
49
#include "
DataFormats/HLTReco/interface/TriggerFilterObjectWithRefs.h
"
50
#include "
DataFormats/HLTReco/interface/TriggerTypeDefs.h
"
51
52
class
HLTPixelAsymmetryFilter
:
public
HLTFilter
{
53
public
:
54
explicit
HLTPixelAsymmetryFilter
(
const
edm::ParameterSet
&);
55
~HLTPixelAsymmetryFilter
()
override
;
56
static
void
fillDescriptions
(
edm::ConfigurationDescriptions
& descriptions);
57
58
private
:
59
bool
hltFilter
(
edm::Event
&,
const
edm::EventSetup
&,
trigger::TriggerFilterObjectWithRefs
& filterproduct)
const override
;
60
61
edm::EDGetTokenT<edmNew::DetSetVector<SiPixelCluster>
>
inputToken_
;
62
edm::InputTag
inputTag_
;
// input tag identifying product containing pixel clusters
63
double
min_asym_
;
// minimum asymmetry
64
double
max_asym_
;
// maximum asymmetry
65
double
clus_thresh_
;
// minimum charge for a cluster to be selected (in e-)
66
double
bmincharge_
;
// minimum average charge in the barrel (bpix, in e-)
67
68
};
69
70
#endif
71
MessageLogger.h
HLTPixelAsymmetryFilter
Definition:
HLTPixelAsymmetryFilter.h:52
HLTPixelAsymmetryFilter::hltFilter
bool hltFilter(edm::Event &, const edm::EventSetup &, trigger::TriggerFilterObjectWithRefs &filterproduct) const override
Definition:
HLTPixelAsymmetryFilter.cc:50
HLTFilter.h
PixelEndcapName.h
HLTPixelAsymmetryFilter::max_asym_
double max_asym_
Definition:
HLTPixelAsymmetryFilter.h:64
HLTPixelAsymmetryFilter::bmincharge_
double bmincharge_
Definition:
HLTPixelAsymmetryFilter.h:66
HLTPixelAsymmetryFilter::inputToken_
edm::EDGetTokenT< edmNew::DetSetVector< SiPixelCluster > > inputToken_
Definition:
HLTPixelAsymmetryFilter.h:61
Event.h
HLTPixelAsymmetryFilter::HLTPixelAsymmetryFilter
HLTPixelAsymmetryFilter(const edm::ParameterSet &)
Definition:
HLTPixelAsymmetryFilter.cc:17
trigger::TriggerFilterObjectWithRefs
Definition:
TriggerFilterObjectWithRefs.h:36
EventSetup.h
HLTPixelAsymmetryFilter::min_asym_
double min_asym_
Definition:
HLTPixelAsymmetryFilter.h:63
TriggerFilterObjectWithRefs.h
edm::EDGetTokenT
Definition:
EDGetToken.h:33
ParameterSet.h
ParameterSetDescription.h
HLTFilter
Definition:
HLTFilter.h:28
HLTPixelAsymmetryFilter::inputTag_
edm::InputTag inputTag_
Definition:
HLTPixelAsymmetryFilter.h:62
edm::EventSetup
Definition:
EventSetup.h:57
PixelSubdetector.h
HLTPixelAsymmetryFilter::~HLTPixelAsymmetryFilter
~HLTPixelAsymmetryFilter() override
edm::InputTag
Definition:
InputTag.h:15
edm::ParameterSet
Definition:
ParameterSet.h:36
ConfigurationDescriptions.h
HLTPixelAsymmetryFilter::clus_thresh_
double clus_thresh_
Definition:
HLTPixelAsymmetryFilter.h:65
TriggerTypeDefs.h
edm::Event
Definition:
Event.h:71
SiPixelCluster.h
edm::ConfigurationDescriptions
Definition:
ConfigurationDescriptions.h:28
Handle.h
HLTPixelAsymmetryFilter::fillDescriptions
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
Definition:
HLTPixelAsymmetryFilter.cc:34
PixelBarrelName.h
Generated for CMSSW Reference Manual by
1.8.11