Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-11-14 04:15:52

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   desc.add<uint32_t>("DontCountDetsAboveNClusters", 0);
0041 }
0042 
0043 ClusterChecker::~ClusterChecker() {}
0044 
0045 size_t ClusterChecker::tooManyClusters(const edm::Event& e) const {
0046   if (!doACheck_)
0047     return 0;
0048 
0049   // get special input for cosmic cluster multiplicity filter
0050   edm::Handle<edmNew::DetSetVector<SiStripCluster> > clusterDSV;
0051   e.getByToken(token_sc, clusterDSV);
0052   reco::utils::ClusterTotals totals;
0053   if (!clusterDSV.failedToGet()) {
0054     const edmNew::DetSetVector<SiStripCluster>& input = *clusterDSV;
0055 
0056     if (ignoreDetsAboveNClusters_ == 0) {
0057       totals.strip = input.dataSize();
0058       totals.stripdets = input.size();
0059     } else {
0060       //loop over detectors
0061       totals.strip = 0;
0062       totals.stripdets = 0;
0063       edmNew::DetSetVector<SiStripCluster>::const_iterator DSViter = input.begin(), DSViter_end = input.end();
0064       for (; DSViter != DSViter_end; DSViter++) {
0065         size_t siz = DSViter->size();
0066         if (siz > ignoreDetsAboveNClusters_)
0067           continue;
0068         totals.strip += siz;
0069         totals.stripdets++;
0070       }
0071     }
0072   }
0073   if (totals.strip > int(maxNrOfStripClusters_))
0074     return totals.strip;
0075 
0076   // get special input for pixel cluster multiplicity filter
0077   edm::Handle<edmNew::DetSetVector<SiPixelCluster> > pixelClusterDSV;
0078   e.getByToken(token_pc, pixelClusterDSV);
0079   if (!pixelClusterDSV.failedToGet()) {
0080     const edmNew::DetSetVector<SiPixelCluster>& input = *pixelClusterDSV;
0081 
0082     if (ignoreDetsAboveNClusters_ == 0) {
0083       totals.pixel = input.dataSize();
0084       totals.pixeldets = input.size();
0085     } else {
0086       //loop over detectors
0087       totals.pixel = 0;
0088       totals.pixeldets = 0;
0089       edmNew::DetSetVector<SiPixelCluster>::const_iterator DSViter = input.begin(), DSViter_end = input.end();
0090       for (; DSViter != DSViter_end; DSViter++) {
0091         size_t siz = DSViter->size();
0092         if (siz > ignoreDetsAboveNClusters_)
0093           continue;
0094         totals.pixel += siz;
0095         totals.pixeldets++;
0096       }
0097     }
0098   } else {
0099     //say something's wrong.
0100     edm::LogError("ClusterChecker")
0101         << "could not get any SiPixel cluster collections of type edm::DetSetVector<SiPixelCluster>  with label: "
0102         << pixelClusterCollectionInputTag_;
0103     totals.pixel = 999999;
0104   }
0105   if (totals.pixel > int(maxNrOfPixelClusters_))
0106     return totals.pixel;
0107 
0108   if (!selector_(totals))
0109     return totals.strip;
0110   return 0;
0111 }