Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:28:54

0001 #include "RecoTracker/TkSeedGenerator/interface/ClusterChecker.h"
0002 
0003 #include "FWCore/Framework/interface/Event.h"
0004 #include "DataFormats/Common/interface/Handle.h"
0005 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0006 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0007 
0008 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0009 
0010 ClusterChecker::ClusterChecker(const edm::ParameterSet& conf, edm::ConsumesCollector&& iC) : ClusterChecker(conf, iC) {}
0011 
0012 ClusterChecker::ClusterChecker(const edm::ParameterSet& conf, edm::ConsumesCollector& iC)
0013     : doACheck_(conf.getParameter<bool>("doClusterCheck")),
0014       selector_(conf.getParameter<bool>("doClusterCheck") && conf.existsAs<std::string>("cut")
0015                     ? conf.getParameter<std::string>("cut")
0016                     : "") {
0017   if (doACheck_) {
0018     clusterCollectionInputTag_ = conf.getParameter<edm::InputTag>("ClusterCollectionLabel");
0019     pixelClusterCollectionInputTag_ = conf.getParameter<edm::InputTag>("PixelClusterCollectionLabel");
0020     token_sc = iC.consumes<edmNew::DetSetVector<SiStripCluster> >(clusterCollectionInputTag_);
0021     token_pc = iC.consumes<edmNew::DetSetVector<SiPixelCluster> >(pixelClusterCollectionInputTag_);
0022     maxNrOfStripClusters_ = conf.getParameter<unsigned int>("MaxNumberOfStripClusters");
0023     maxNrOfPixelClusters_ = conf.getParameter<unsigned int>("MaxNumberOfPixelClusters");
0024     if (conf.existsAs<uint32_t>("DontCountDetsAboveNClusters")) {
0025       ignoreDetsAboveNClusters_ = conf.getParameter<uint32_t>("DontCountDetsAboveNClusters");
0026     } else {
0027       ignoreDetsAboveNClusters_ = 0;
0028     }
0029   }
0030 }
0031 
0032 void ClusterChecker::fillDescriptions(edm::ParameterSetDescription& desc) {
0033   desc.add<bool>("doClusterCheck", true);
0034   desc.add<unsigned>("MaxNumberOfStripClusters", 400000);
0035   desc.add<edm::InputTag>("ClusterCollectionLabel", edm::InputTag("siStripClusters"));
0036   desc.add<unsigned>("MaxNumberOfPixelClusters", 40000);
0037   desc.add<edm::InputTag>("PixelClusterCollectionLabel", edm::InputTag("siPixelClusters"));
0038   desc.add<std::string>("cut",
0039                         "strip < 400000 && pixel < 40000 && (strip < 50000 + 10*pixel) && (pixel < 5000 + 0.1*strip)");
0040 }
0041 
0042 ClusterChecker::~ClusterChecker() {}
0043 
0044 size_t ClusterChecker::tooManyClusters(const edm::Event& e) const {
0045   if (!doACheck_)
0046     return 0;
0047 
0048   // get special input for cosmic cluster multiplicity filter
0049   edm::Handle<edmNew::DetSetVector<SiStripCluster> > clusterDSV;
0050   e.getByToken(token_sc, clusterDSV);
0051   reco::utils::ClusterTotals totals;
0052   if (!clusterDSV.failedToGet()) {
0053     const edmNew::DetSetVector<SiStripCluster>& input = *clusterDSV;
0054 
0055     if (ignoreDetsAboveNClusters_ == 0) {
0056       totals.strip = input.dataSize();
0057       totals.stripdets = input.size();
0058     } else {
0059       //loop over detectors
0060       totals.strip = 0;
0061       totals.stripdets = 0;
0062       edmNew::DetSetVector<SiStripCluster>::const_iterator DSViter = input.begin(), DSViter_end = input.end();
0063       for (; DSViter != DSViter_end; DSViter++) {
0064         size_t siz = DSViter->size();
0065         if (siz > ignoreDetsAboveNClusters_)
0066           continue;
0067         totals.strip += siz;
0068         totals.stripdets++;
0069       }
0070     }
0071   }
0072   if (totals.strip > int(maxNrOfStripClusters_))
0073     return totals.strip;
0074 
0075   // get special input for pixel cluster multiplicity filter
0076   edm::Handle<edmNew::DetSetVector<SiPixelCluster> > pixelClusterDSV;
0077   e.getByToken(token_pc, pixelClusterDSV);
0078   if (!pixelClusterDSV.failedToGet()) {
0079     const edmNew::DetSetVector<SiPixelCluster>& input = *pixelClusterDSV;
0080 
0081     if (ignoreDetsAboveNClusters_ == 0) {
0082       totals.pixel = input.dataSize();
0083       totals.pixeldets = input.size();
0084     } else {
0085       //loop over detectors
0086       totals.pixel = 0;
0087       totals.pixeldets = 0;
0088       edmNew::DetSetVector<SiPixelCluster>::const_iterator DSViter = input.begin(), DSViter_end = input.end();
0089       for (; DSViter != DSViter_end; DSViter++) {
0090         size_t siz = DSViter->size();
0091         if (siz > ignoreDetsAboveNClusters_)
0092           continue;
0093         totals.pixel += siz;
0094         totals.pixeldets++;
0095       }
0096     }
0097   } else {
0098     //say something's wrong.
0099     edm::LogError("ClusterChecker")
0100         << "could not get any SiPixel cluster collections of type edm::DetSetVector<SiPixelCluster>  with label: "
0101         << pixelClusterCollectionInputTag_;
0102     totals.pixel = 999999;
0103   }
0104   if (totals.pixel > int(maxNrOfPixelClusters_))
0105     return totals.pixel;
0106 
0107   if (!selector_(totals))
0108     return totals.strip;
0109   return 0;
0110 }