CMS 3D CMS Logo

Public Member Functions | Private Attributes

ClusterChecker Class Reference

#include <ClusterChecker.h>

List of all members.

Public Member Functions

 ClusterChecker (edm::ParameterSet conf)
size_t tooManyClusters (const edm::Event &e)

Private Attributes

edm::InputTag clusterCollectionInputTag_
bool doACheck_
uint32_t ignoreDetsAboveNClusters_
uint32_t maxNrOfCosmicClusters_
uint32_t maxNrOfPixelClusters_
edm::InputTag pixelClusterCollectionInputTag_

Detailed Description

Definition at line 16 of file ClusterChecker.h.


Constructor & Destructor Documentation

ClusterChecker::ClusterChecker ( edm::ParameterSet  conf) [inline]

Definition at line 18 of file ClusterChecker.h.

References clusterCollectionInputTag_, doACheck_, edm::ParameterSet::existsAs(), edm::ParameterSet::getParameter(), ignoreDetsAboveNClusters_, maxNrOfCosmicClusters_, maxNrOfPixelClusters_, and pixelClusterCollectionInputTag_.

                                       :
    doACheck_(conf.getParameter<bool>("doClusterCheck"))
    {
      if (doACheck_){
        clusterCollectionInputTag_ = conf.getParameter<edm::InputTag>("ClusterCollectionLabel");
        maxNrOfCosmicClusters_ = conf.getParameter<unsigned int>("MaxNumberOfCosmicClusters");
        pixelClusterCollectionInputTag_ = conf.getParameter<edm::InputTag>("PixelClusterCollectionLabel");
        maxNrOfPixelClusters_ = conf.getParameter<unsigned int>("MaxNumberOfPixelClusters");
        if (conf.existsAs<uint32_t>("DontCountDetsAboveNClusters")) {
            ignoreDetsAboveNClusters_ = conf.getParameter<uint32_t>("DontCountDetsAboveNClusters");
        } else {
            ignoreDetsAboveNClusters_ = 0;
        }
      }
    }

Member Function Documentation

size_t ClusterChecker::tooManyClusters ( const edm::Event e) [inline]

Definition at line 34 of file ClusterChecker.h.

References edmNew::DetSetVector< T >::begin(), clusterCollectionInputTag_, edmNew::DetSetVector< T >::dataSize(), doACheck_, edmNew::DetSetVector< T >::end(), edm::HandleBase::failedToGet(), edm::Event::getByLabel(), ignoreDetsAboveNClusters_, collect_tpl::input, maxNrOfCosmicClusters_, maxNrOfPixelClusters_, pixelClusterCollectionInputTag_, and edmNew::DetSetVector< T >::size().

Referenced by CosmicSeedGenerator::produce(), SeedGeneratorFromRegionHitsEDProducer::produce(), CtfSpecialSeedGenerator::produce(), SimpleCosmicBONSeeder::produce(), and RoadSearchSeedFinder::produce().

                                            {
    if (!doACheck_) return 0;

    // get special input for cosmic cluster multiplicity filter
    edm::Handle<edmNew::DetSetVector<SiStripCluster> > clusterDSV;
    e.getByLabel(clusterCollectionInputTag_, clusterDSV);
    unsigned int totalClusters = 0;
    if (!clusterDSV.failedToGet()) {
      const edmNew::DetSetVector<SiStripCluster> & input = *clusterDSV;

      if (ignoreDetsAboveNClusters_ == 0) {
        totalClusters = input.dataSize();
      } else {
          //loop over detectors
          edmNew::DetSetVector<SiStripCluster>::const_iterator DSViter=input.begin(), DSViter_end=input.end();
          for (; DSViter!=DSViter_end; DSViter++ ) {
            size_t siz = DSViter->size();
            if (siz > ignoreDetsAboveNClusters_) continue;
            totalClusters += siz; 
          }
      }
    }
    else{
      edm::Handle<edm::LazyGetter<SiStripCluster> > lazyGH;
      e.getByLabel(clusterCollectionInputTag_, lazyGH);
      if (!lazyGH.failedToGet()){
        totalClusters = lazyGH->size();
      }else{
        //say something's wrong.
        edm::LogError("ClusterChecker")<<"could not get any SiStrip cluster collections of type edm::DetSetVector<SiStripCluster> or edm::LazyGetter<SiStripCluster, with label: "<<clusterCollectionInputTag_;
        totalClusters = 999999;
      }
    }
    if (totalClusters > maxNrOfCosmicClusters_) return totalClusters;

    // get special input for pixel cluster multiplicity filter
    edm::Handle<edmNew::DetSetVector<SiPixelCluster> > pixelClusterDSV;
    e.getByLabel(pixelClusterCollectionInputTag_, pixelClusterDSV);
    unsigned int totalPixelClusters = 0;
    if (!pixelClusterDSV.failedToGet()) {
      const edmNew::DetSetVector<SiPixelCluster> & input = *pixelClusterDSV;

      if (ignoreDetsAboveNClusters_ == 0) {
        totalPixelClusters = input.dataSize();
      } else {
          //loop over detectors
          edmNew::DetSetVector<SiPixelCluster>::const_iterator DSViter=input.begin(), DSViter_end=input.end();
          for (; DSViter!=DSViter_end; DSViter++ ) {
            size_t siz = DSViter->size();
            if (siz > ignoreDetsAboveNClusters_) continue;
            totalPixelClusters += siz; 
          }
      }
    }
    else{
      //say something's wrong.
      edm::LogError("ClusterChecker")<<"could not get any SiPixel cluster collections of type edm::DetSetVector<SiPixelCluster>  with label: "<<pixelClusterCollectionInputTag_;
      totalPixelClusters = 999999;
    }
    if (totalPixelClusters > maxNrOfPixelClusters_) return totalPixelClusters;

    return 0;

  }

Member Data Documentation

Definition at line 101 of file ClusterChecker.h.

Referenced by ClusterChecker(), and tooManyClusters().

bool ClusterChecker::doACheck_ [private]

Definition at line 100 of file ClusterChecker.h.

Referenced by ClusterChecker(), and tooManyClusters().

Definition at line 105 of file ClusterChecker.h.

Referenced by ClusterChecker(), and tooManyClusters().

Definition at line 103 of file ClusterChecker.h.

Referenced by ClusterChecker(), and tooManyClusters().

Definition at line 104 of file ClusterChecker.h.

Referenced by ClusterChecker(), and tooManyClusters().

Definition at line 102 of file ClusterChecker.h.

Referenced by ClusterChecker(), and tooManyClusters().