CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_1/src/RecoTracker/TkSeedGenerator/src/ClusterChecker.cc

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