Skip to content

Commit

Permalink
split layer count into BPix and FPix
Browse files Browse the repository at this point in the history
  • Loading branch information
JanFSchulte committed May 25, 2018
1 parent e6ba533 commit 5ab05ab
Showing 1 changed file with 31 additions and 15 deletions.
46 changes: 31 additions & 15 deletions HLTrigger/special/plugins/HLTPixelActivityFilter.cc
Expand Up @@ -25,8 +25,10 @@ class HLTPixelActivityFilter : public HLTFilter {
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
unsigned int min_layersBPix_; // minimum number of clusters
unsigned int max_layersBPix_; // maximum number of clusters
unsigned int min_layersFPix_; // minimum number of clusters
unsigned int max_layersFPix_; // maximum number of clusters
edm::EDGetTokenT<edmNew::DetSetVector<SiPixelCluster> > inputToken_;

};
Expand All @@ -45,18 +47,25 @@ HLTPixelActivityFilter::HLTPixelActivityFilter(const edm::ParameterSet& config)
inputTag_ (config.getParameter<edm::InputTag>("inputTag")),
min_clusters_ (config.getParameter<unsigned int>("minClusters")),
max_clusters_ (config.getParameter<unsigned int>("maxClusters")),
min_layers_ (config.getParameter<unsigned int>("minLayers")),
max_layers_ (config.getParameter<unsigned int>("maxLayers"))
min_layersBPix_ (config.getParameter<unsigned int>("minLayersBPix")),
max_layersBPix_ (config.getParameter<unsigned int>("maxLayersBPix")),
min_layersFPix_ (config.getParameter<unsigned int>("minLayersFPix")),
max_layersFPix_ (config.getParameter<unsigned int>("maxLayersFPix"))
{
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";
if(min_layersBPix_ > 0)
LogDebug("") << "Also requesting at least " << min_layersBPix_ << " layers with clusters in BPix";
if(max_layersBPix_ > 0)
LogDebug("") << "...but no more than " << max_layersBPix_ << " layers with clusters in BPix";
if(min_layersFPix_ > 0)
LogDebug("") << "Also requesting at least " << min_layersFPix_ << " layers with clusters in FPix";
if(max_layersFPix_ > 0)
LogDebug("") << "...but no more than " << max_layersFPix_ << " layers with clusters in FPix";


}

Expand All @@ -69,9 +78,12 @@ 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);
desc.add<unsigned int>("minLayersBPix",0);
desc.add<unsigned int>("maxLayersBPix",0);
desc.add<unsigned int>("minLayersFPix",0);
desc.add<unsigned int>("maxLayersFPix",0);
descriptions.add("hltPixelActivityFilter",desc);

}

//
Expand All @@ -97,12 +109,13 @@ bool HLTPixelActivityFilter::hltFilter(edm::Event& event, const edm::EventSetup&
if(max_clusters_ > 0)
accept &= (clusterSize <= max_clusters_);

if (min_layers_ > 0 || max_layers_ > 0){
if (min_layersBPix_ > 0 || max_layersBPix_ > 0 || min_layersFPix_ > 0 || max_layersFPix_ > 0){

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

edmNew::DetSetVector<SiPixelCluster>::const_iterator DSViter=clusters.begin();
Expand All @@ -129,9 +142,12 @@ bool HLTPixelActivityFilter::hltFilter(edm::Event& event, const edm::EventSetup&

}
}
layerCount = foundLayersB.size()+foundLayersEp.size()+foundLayersEm.size();
if (max_layers_ > 0) accept &= (layerCount <= max_layers_);
if (min_layers_ > 0) accept &= (layerCount >= min_layers_);
layerCountBPix = foundLayersB.size();
layerCountFPix = foundLayersEp.size()+foundLayersEm.size();
if (max_layersBPix_ > 0) accept &= (layerCountBPix <= max_layersBPix_);
if (min_layersBPix_ > 0) accept &= (layerCountBPix >= min_layersBPix_);
if (max_layersFPix_ > 0) accept &= (layerCountFPix <= max_layersFPix_);
if (min_layersFPix_ > 0) accept &= (layerCountFPix >= min_layersFPix_);


}
Expand Down

0 comments on commit 5ab05ab

Please sign in to comment.