CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_1_8_patch12/src/RecoTracker/SpecialSeedGenerators/interface/ClusterChecker.h

Go to the documentation of this file.
00001 #ifndef RecoTracker_SpecialSeedGenerators_ClusterChecker_H
00002 #define RecoTracker_SpecialSeedGenerators_ClusterChecker_H
00003 
00004 #include "FWCore/Framework/interface/Event.h"
00005 #include "DataFormats/Common/interface/Handle.h"
00006 #include "FWCore/Utilities/interface/InputTag.h"
00007 
00008 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00009 #include "DataFormats/Common/interface/LazyGetter.h"
00010 #include "DataFormats/SiStripCluster/interface/SiStripCluster.h"
00011 #include "DataFormats/SiPixelCluster/interface/SiPixelCluster.h"
00012 #include "DataFormats/Common/interface/DetSetVectorNew.h"
00013 
00014 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00015 
00016 class ClusterChecker {
00017  public: 
00018   ClusterChecker(edm::ParameterSet conf) :
00019     doACheck_(conf.getParameter<bool>("doClusterCheck"))
00020     {
00021       if (doACheck_){
00022         clusterCollectionInputTag_ = conf.getParameter<edm::InputTag>("ClusterCollectionLabel");
00023         maxNrOfCosmicClusters_ = conf.getParameter<unsigned int>("MaxNumberOfCosmicClusters");
00024         pixelClusterCollectionInputTag_ = conf.getParameter<edm::InputTag>("PixelClusterCollectionLabel");
00025         maxNrOfPixelClusters_ = conf.getParameter<unsigned int>("MaxNumberOfPixelClusters");
00026         if (conf.existsAs<uint32_t>("DontCountDetsAboveNClusters")) {
00027             ignoreDetsAboveNClusters_ = conf.getParameter<uint32_t>("DontCountDetsAboveNClusters");
00028         } else {
00029             ignoreDetsAboveNClusters_ = 0;
00030         }
00031       }
00032     }
00033     
00034   size_t tooManyClusters(const edm::Event & e){
00035     if (!doACheck_) return 0;
00036 
00037     // get special input for cosmic cluster multiplicity filter
00038     edm::Handle<edmNew::DetSetVector<SiStripCluster> > clusterDSV;
00039     e.getByLabel(clusterCollectionInputTag_, clusterDSV);
00040     unsigned int totalClusters = 0;
00041     if (!clusterDSV.failedToGet()) {
00042       const edmNew::DetSetVector<SiStripCluster> & input = *clusterDSV;
00043 
00044       if (ignoreDetsAboveNClusters_ == 0) {
00045         totalClusters = input.dataSize();
00046       } else {
00047           //loop over detectors
00048           edmNew::DetSetVector<SiStripCluster>::const_iterator DSViter=input.begin(), DSViter_end=input.end();
00049           for (; DSViter!=DSViter_end; DSViter++ ) {
00050             size_t siz = DSViter->size();
00051             if (siz > ignoreDetsAboveNClusters_) continue;
00052             totalClusters += siz; 
00053           }
00054       }
00055     }
00056     else{
00057       edm::Handle<edm::LazyGetter<SiStripCluster> > lazyGH;
00058       e.getByLabel(clusterCollectionInputTag_, lazyGH);
00059       if (!lazyGH.failedToGet()){
00060         totalClusters = lazyGH->size();
00061       }else{
00062         //say something's wrong.
00063         edm::LogError("ClusterChecker")<<"could not get any SiStrip cluster collections of type edm::DetSetVector<SiStripCluster> or edm::LazyGetter<SiStripCluster, with label: "<<clusterCollectionInputTag_;
00064         totalClusters = 999999;
00065       }
00066     }
00067     if (totalClusters > maxNrOfCosmicClusters_) return totalClusters;
00068 
00069     // get special input for pixel cluster multiplicity filter
00070     edm::Handle<edmNew::DetSetVector<SiPixelCluster> > pixelClusterDSV;
00071     e.getByLabel(pixelClusterCollectionInputTag_, pixelClusterDSV);
00072     unsigned int totalPixelClusters = 0;
00073     if (!pixelClusterDSV.failedToGet()) {
00074       const edmNew::DetSetVector<SiPixelCluster> & input = *pixelClusterDSV;
00075 
00076       if (ignoreDetsAboveNClusters_ == 0) {
00077         totalPixelClusters = input.dataSize();
00078       } else {
00079           //loop over detectors
00080           edmNew::DetSetVector<SiPixelCluster>::const_iterator DSViter=input.begin(), DSViter_end=input.end();
00081           for (; DSViter!=DSViter_end; DSViter++ ) {
00082             size_t siz = DSViter->size();
00083             if (siz > ignoreDetsAboveNClusters_) continue;
00084             totalPixelClusters += siz; 
00085           }
00086       }
00087     }
00088     else{
00089       //say something's wrong.
00090       edm::LogError("ClusterChecker")<<"could not get any SiPixel cluster collections of type edm::DetSetVector<SiPixelCluster>  with label: "<<pixelClusterCollectionInputTag_;
00091       totalPixelClusters = 999999;
00092     }
00093     if (totalPixelClusters > maxNrOfPixelClusters_) return totalPixelClusters;
00094 
00095     return 0;
00096 
00097   }
00098 
00099  private: 
00100   bool doACheck_;
00101   edm::InputTag clusterCollectionInputTag_;
00102   edm::InputTag pixelClusterCollectionInputTag_;
00103   uint32_t maxNrOfCosmicClusters_;
00104   uint32_t maxNrOfPixelClusters_;
00105   uint32_t ignoreDetsAboveNClusters_;
00106 };
00107 
00108 #endif