CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Member Functions | Private Attributes
HLTPixelAsymmetryFilter Class Reference

#include <HLTPixelAsymmetryFilter.h>

Inheritance diagram for HLTPixelAsymmetryFilter:
HLTFilter edm::EDFilter edm::ProducerBase edm::ProductRegistryHelper

Public Member Functions

 HLTPixelAsymmetryFilter (const edm::ParameterSet &)
 
 ~HLTPixelAsymmetryFilter ()
 
- Public Member Functions inherited from HLTFilter
 HLTFilter (const edm::ParameterSet &config)
 
int module () const
 
const std::string * moduleLabel () const
 
int path () const
 
const std::string * pathName () const
 
std::pair< int, int > pmid () const
 
bool saveTags () const
 
virtual ~HLTFilter ()
 
- Public Member Functions inherited from edm::EDFilter
 EDFilter ()
 
virtual ~EDFilter ()
 
- Public Member Functions inherited from edm::ProducerBase
 ProducerBase ()
 
void registerProducts (ProducerBase *, ProductRegistry *, ModuleDescription const &)
 
boost::function< void(const
BranchDescription &)> 
registrationCallback () const
 used by the fwk to register list of products More...
 
virtual ~ProducerBase ()
 

Private Member Functions

virtual bool hltFilter (edm::Event &, const edm::EventSetup &, trigger::TriggerFilterObjectWithRefs &filterproduct)
 

Private Attributes

double bmincharge_
 
double clus_thresh_
 
edm::InputTag inputTag_
 
double max_asym_
 
double min_asym_
 

Additional Inherited Members

- Public Types inherited from edm::EDFilter
typedef EDFilter ModuleType
 
typedef WorkerT< EDFilterWorkerType
 
- Public Types inherited from edm::ProducerBase
typedef
ProductRegistryHelper::TypeLabelList 
TypeLabelList
 
- Static Public Member Functions inherited from HLTFilter
static void makeHLTFilterDescription (edm::ParameterSetDescription &desc)
 
- Static Public Member Functions inherited from edm::EDFilter
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &)
 
- Protected Member Functions inherited from edm::EDFilter
CurrentProcessingContext const * currentContext () const
 
- Protected Member Functions inherited from edm::ProducerBase
template<class TProducer , class TMethod >
void callWhenNewProductsRegistered (TProducer *iProd, TMethod iMethod)
 

Detailed Description

Definition at line 50 of file HLTPixelAsymmetryFilter.h.

Constructor & Destructor Documentation

HLTPixelAsymmetryFilter::HLTPixelAsymmetryFilter ( const edm::ParameterSet config)
explicit

Definition at line 17 of file HLTPixelAsymmetryFilter.cc.

References bmincharge_, clus_thresh_, inputTag_, LogDebug, max_asym_, and min_asym_.

17  : HLTFilter(config),
18  inputTag_ (config.getParameter<edm::InputTag>("inputTag")),
19  min_asym_ (config.getParameter<double>("MinAsym")),
20  max_asym_ (config.getParameter<double>("MaxAsym")),
21  clus_thresh_ (config.getParameter<double>("MinCharge")),
22  bmincharge_ (config.getParameter<double>("MinBarrel"))
23 {
24  LogDebug("") << "Using the " << inputTag_ << " input collection";
25  LogDebug("") << "Requesting events with a charge repartition asymmetry between " << min_asym_ << " and " << max_asym_;
26  LogDebug("") << "Mean cluster charge in the barrel should be higher than" << bmincharge_ << " electrons ";
27  LogDebug("") << "Only clusters with a charge larger than " << clus_thresh_ << " electrons will be used ";
28 }
#define LogDebug(id)
T getParameter(std::string const &) const
HLTFilter(const edm::ParameterSet &config)
Definition: HLTFilter.cc:18
HLTPixelAsymmetryFilter::~HLTPixelAsymmetryFilter ( )

Definition at line 30 of file HLTPixelAsymmetryFilter.cc.

31 {
32 }

Member Function Documentation

bool HLTPixelAsymmetryFilter::hltFilter ( edm::Event event,
const edm::EventSetup iSetup,
trigger::TriggerFilterObjectWithRefs filterproduct 
)
privatevirtual

Implements HLTFilter.

Definition at line 39 of file HLTPixelAsymmetryFilter.cc.

References accept(), trigger::TriggerFilterObjectWithRefs::addCollectionTag(), Reference_intrackfit_cff::barrel, begin, bmincharge_, clus_thresh_, cond::rpcobgas::detid, end, Reference_intrackfit_cff::endcap, PixelEndcapName::halfCylinder(), i, inputTag_, LogDebug, PixelEndcapName::mI, min_asym_, PixelEndcapName::mO, PixelEndcapName::pI, PixelSubdetector::PixelBarrel, PixelSubdetector::PixelEndcap, PixelEndcapName::pO, position, HLTFilter::saveTags(), and DetId::subdetId().

40 {
41  // All HLT filters must create and fill an HLT filter object,
42  // recording any reconstructed physics objects satisfying (or not)
43  // this HLT filter, and place it in the Event.
44 
45  // The filter object
46  if (saveTags()) filterproduct.addCollectionTag(inputTag_);
47 
48  // get hold of products from Event
50  event.getByLabel(inputTag_, clusterColl);
51 
52  unsigned int clusterSize = clusterColl->dataSize();
53  LogDebug("") << "Number of clusters accepted: " << clusterSize;
54 
55  bool accept = (clusterSize >= 2); // Not necessary to go further in this case
56 
57  if (!accept) return accept;
58 
59  double asym_pix_1 = -1;
60  double asym_pix_2 = -1;
61 
62  int n_clus[3] = {0,0,0};
63  double e_clus[3] = {0.,0.,0.};
64 
65  for (edmNew::DetSetVector<SiPixelCluster>::const_iterator DSViter=clusterColl->begin(); DSViter!=clusterColl->end();DSViter++ )
66  {
69  uint32_t detid = DSViter->id();
70 
71  bool barrel = DetId(detid).subdetId() == static_cast<int>(PixelSubdetector::PixelBarrel);
72  bool endcap = DetId(detid).subdetId() == static_cast<int>(PixelSubdetector::PixelEndcap);
73 
74  int detpos = -1;
75 
76  // First we look if we are in the endcap or in the barrel PIXELS
77  // Asymmetry is computed with endcap pixels only
78 
79  if (endcap)
80  {
82 
83  if (position == PixelEndcapName::mI || position == PixelEndcapName::mO) // FPIX-
84  detpos = 0;
85 
86  if (position == PixelEndcapName::pI || position == PixelEndcapName::pO) // FPIX+
87  detpos = 2;
88  }
89 
90  if (barrel)
91  detpos = 1;
92 
93  if (detpos<0) continue;
94 
95  for(edmNew::DetSet<SiPixelCluster>::const_iterator iter=begin;iter!=end;++iter)
96  {
97  if (iter->charge()>clus_thresh_ )
98  {
99  ++n_clus[detpos];
100  e_clus[detpos]+=iter->charge();
101  }
102  }
103  } // End of cluster loop
104 
105  for (int i=0;i<3;++i)
106  {
107  if (n_clus[i])
108  e_clus[i] /= n_clus[i];
109  }
110 
111  if (e_clus[1] < bmincharge_) return false; // Reject event if the Barrel mean charge is too low
112 
113 
114  if (e_clus[0]+e_clus[2] != 0) // Compute asyms, if applicable
115  {
116  asym_pix_1 = e_clus[0]/(e_clus[0]+e_clus[2]);
117  asym_pix_2 = e_clus[2]/(e_clus[0]+e_clus[2]);
118  }
119  else // Otherwise reject event
120  {
121  return false;
122  }
123 
124  bool pkam_1 = (asym_pix_1 <= max_asym_ && asym_pix_1 >= min_asym_);
125  bool pkam_2 = (asym_pix_2 <= max_asym_ && asym_pix_2 >= min_asym_);
126 
127  if (pkam_1 || pkam_2) return accept; // Final selection
128 
129  return false;
130 }
#define LogDebug(id)
int i
Definition: DBlmapReader.cc:9
boost::transform_iterator< IterHelp, const_IdIter > const_iterator
data_type const * const_iterator
Definition: DetSetNew.h:25
bool accept(const edm::Event &event, const edm::TriggerResults &triggerTable, const std::string &triggerPath)
Definition: TopDQMHelpers.h:22
#define end
Definition: vmac.h:38
int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:39
Definition: DetId.h:20
void addCollectionTag(const edm::InputTag &collectionTag)
collectionTags
bool saveTags() const
Definition: HLTFilter.h:45
#define begin
Definition: vmac.h:31
static int position[264][3]
Definition: ReadPGInfo.cc:509
HalfCylinder halfCylinder() const

Member Data Documentation

double HLTPixelAsymmetryFilter::bmincharge_
private

Definition at line 62 of file HLTPixelAsymmetryFilter.h.

Referenced by hltFilter(), and HLTPixelAsymmetryFilter().

double HLTPixelAsymmetryFilter::clus_thresh_
private

Definition at line 61 of file HLTPixelAsymmetryFilter.h.

Referenced by hltFilter(), and HLTPixelAsymmetryFilter().

edm::InputTag HLTPixelAsymmetryFilter::inputTag_
private

Definition at line 58 of file HLTPixelAsymmetryFilter.h.

Referenced by hltFilter(), and HLTPixelAsymmetryFilter().

double HLTPixelAsymmetryFilter::max_asym_
private

Definition at line 60 of file HLTPixelAsymmetryFilter.h.

Referenced by HLTPixelAsymmetryFilter().

double HLTPixelAsymmetryFilter::min_asym_
private

Definition at line 59 of file HLTPixelAsymmetryFilter.h.

Referenced by hltFilter(), and HLTPixelAsymmetryFilter().