Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Add layer counting capability to PixelActivityFilter #22952

Merged
merged 1 commit into from Apr 20, 2018
Merged
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Jump to
Jump to file
Failed to load files.
Diff view
Diff view
56 changes: 53 additions & 3 deletions HLTrigger/special/plugins/HLTPixelActivityFilter.cc
Expand Up @@ -5,6 +5,8 @@
#include "FWCore/ParameterSet/interface/ParameterSetDescription.h"

#include "DataFormats/SiPixelCluster/interface/SiPixelCluster.h"
#include "DataFormats/TrackerCommon/interface/TrackerTopology.h"
#include "Geometry/Records/interface/TrackerTopologyRcd.h"

//
// class declaration
Expand All @@ -18,10 +20,13 @@ class HLTPixelActivityFilter : public HLTFilter {

private:
bool hltFilter(edm::Event&, const edm::EventSetup&, trigger::TriggerFilterObjectWithRefs & filterproduct) const override;
// int countLayersWithClusters(edm::Handle<edmNew::DetSetVector<SiPixelCluster> > & clusterCol,const TrackerTopology& tTopo);

edm::InputTag inputTag_; // input tag identifying product containing pixel clusters
unsigned int min_clusters_; // minimum number of clusters
unsigned int max_clusters_; // maximum number of clusters
unsigned int min_layers_; // minimum number of clusters
unsigned int max_layers_; // maximum number of clusters
edm::EDGetTokenT<edmNew::DetSetVector<SiPixelCluster> > inputToken_;

};
Expand All @@ -32,21 +37,27 @@ class HLTPixelActivityFilter : public HLTFilter {
#include "DataFormats/Common/interface/Handle.h"
#include "DataFormats/HLTReco/interface/TriggerFilterObjectWithRefs.h"
#include "DataFormats/HLTReco/interface/TriggerTypeDefs.h"

//
// constructors and destructor
//

HLTPixelActivityFilter::HLTPixelActivityFilter(const edm::ParameterSet& config) : HLTFilter(config),
inputTag_ (config.getParameter<edm::InputTag>("inputTag")),
min_clusters_ (config.getParameter<unsigned int>("minClusters")),
max_clusters_ (config.getParameter<unsigned int>("maxClusters"))
max_clusters_ (config.getParameter<unsigned int>("maxClusters")),
min_layers_ (config.getParameter<unsigned int>("minLayers")),
max_layers_ (config.getParameter<unsigned int>("maxLayers"))
{
inputToken_ = consumes<edmNew::DetSetVector<SiPixelCluster> >(inputTag_);
LogDebug("") << "Using the " << inputTag_ << " input collection";
LogDebug("") << "Requesting at least " << min_clusters_ << " clusters";
if(max_clusters_ > 0)
LogDebug("") << "...but no more than " << max_clusters_ << " clusters";
if(min_layers_ > 0)
LogDebug("") << "Also requesting at least " << min_layers_ << " layers with clusters";
if(max_layers_ > 0)
LogDebug("") << "...but no more than " << max_layers_ << " layers with clusters";

}

HLTPixelActivityFilter::~HLTPixelActivityFilter() = default;
Expand All @@ -58,13 +69,14 @@ HLTPixelActivityFilter::fillDescriptions(edm::ConfigurationDescriptions& descrip
desc.add<edm::InputTag>("inputTag",edm::InputTag("hltSiPixelClusters"));
desc.add<unsigned int>("minClusters",3);
desc.add<unsigned int>("maxClusters",0);
desc.add<unsigned int>("minLayers",0);
desc.add<unsigned int>("maxLayers",0);
descriptions.add("hltPixelActivityFilter",desc);
}

//
// member functions
//

// ------------ method called to produce the data ------------
bool HLTPixelActivityFilter::hltFilter(edm::Event& event, const edm::EventSetup& iSetup, trigger::TriggerFilterObjectWithRefs & filterproduct) const
{
Expand All @@ -85,6 +97,44 @@ bool HLTPixelActivityFilter::hltFilter(edm::Event& event, const edm::EventSetup&
if(max_clusters_ > 0)
accept &= (clusterSize <= max_clusters_);

if (min_layers_ > 0 || max_layers_ > 0){

edm::ESHandle<TrackerTopology> tTopoHandle;
iSetup.get<TrackerTopologyRcd>().get(tTopoHandle);
const TrackerTopology& tTopo = *tTopoHandle;
unsigned int layerCount = 0;
const edmNew::DetSetVector<SiPixelCluster>& clusters = *clusterColl;

edmNew::DetSetVector<SiPixelCluster>::const_iterator DSViter=clusters.begin();


std::vector<int> foundLayersB;
std::vector<int> foundLayersEp;
std::vector<int> foundLayersEm;
for ( ; DSViter != clusters.end() ; DSViter++) {
unsigned int detid = DSViter->detId();
DetId detIdObject( detid );
const auto nCluster = DSViter->size();
const auto subdet = detIdObject.subdetId();
if (subdet == PixelSubdetector::PixelBarrel){
if(!(std::find(foundLayersB.begin(), foundLayersB.end(), tTopo.layer(detIdObject)) != foundLayersB.end()) && nCluster > 0) foundLayersB.push_back(tTopo.layer(detIdObject));
}
else if (subdet ==PixelSubdetector::PixelEndcap){
if (tTopo.side(detIdObject) == 2){
if(!(std::find(foundLayersEp.begin(), foundLayersEp.end(), tTopo.layer(detIdObject)) != foundLayersEp.end()) && nCluster > 0) foundLayersEp.push_back(tTopo.layer(detIdObject));
}
else if (tTopo.side(detIdObject) == 1){
if(!(std::find(foundLayersEm.begin(), foundLayersEm.end(), tTopo.layer(detIdObject)) != foundLayersEm.end()) && nCluster > 0) foundLayersEm.push_back(tTopo.layer(detIdObject));
}

}
}
layerCount = foundLayersB.size()+foundLayersEp.size()+foundLayersEm.size();
if (max_layers_ > 0) accept &= (layerCount <= max_layers_);
if (min_layers_ > 0) accept &= (layerCount >= min_layers_);


}
// return with final filter decision
return accept;
}
Expand Down