CMS 3D CMS Logo

ClusterChecker.cc
Go to the documentation of this file.
2 
7 
9 
12  ClusterChecker(conf, iC)
13 {}
14 
17  doACheck_(conf.getParameter<bool>("doClusterCheck")),
18  selector_(conf.getParameter<bool>("doClusterCheck") && conf.existsAs<std::string>("cut") ?
19  conf.getParameter<std::string>("cut") :
20  "")
21 {
22  if (doACheck_){
23  clusterCollectionInputTag_ = conf.getParameter<edm::InputTag>("ClusterCollectionLabel");
24  pixelClusterCollectionInputTag_ = conf.getParameter<edm::InputTag>("PixelClusterCollectionLabel");
27  maxNrOfCosmicClusters_ = conf.getParameter<unsigned int>("MaxNumberOfCosmicClusters");
28  maxNrOfPixelClusters_ = conf.getParameter<unsigned int>("MaxNumberOfPixelClusters");
29  if (conf.existsAs<uint32_t>("DontCountDetsAboveNClusters")) {
30  ignoreDetsAboveNClusters_ = conf.getParameter<uint32_t>("DontCountDetsAboveNClusters");
31  } else {
33  }
34  }
35 }
36 
38  desc.add<bool>("doClusterCheck", true);
39  desc.add<unsigned>("MaxNumberOfCosmicClusters", 400000);
40  desc.add<edm::InputTag>("ClusterCollectionLabel", edm::InputTag("siStripClusters"));
41  desc.add<unsigned>("MaxNumberOfPixelClusters", 40000);
42  desc.add<edm::InputTag>("PixelClusterCollectionLabel", edm::InputTag("siPixelClusters"));
43  desc.add<std::string>("cut", "strip < 400000 && pixel < 40000 && (strip < 50000 + 10*pixel) && (pixel < 5000 + 0.1*strip)");
44 }
45 
46 
48 {
49 }
50 
52 {
53  if (!doACheck_) return 0;
54 
55  // get special input for cosmic cluster multiplicity filter
57  e.getByToken(token_sc, clusterDSV);
59  if (!clusterDSV.failedToGet()) {
60  const edmNew::DetSetVector<SiStripCluster> & input = *clusterDSV;
61 
62  if (ignoreDetsAboveNClusters_ == 0) {
63  totals.strip = input.dataSize();
64  totals.stripdets = input.size();
65  } else {
66  //loop over detectors
67  totals.strip = 0;
68  totals.stripdets = 0;
69  edmNew::DetSetVector<SiStripCluster>::const_iterator DSViter=input.begin(), DSViter_end=input.end();
70  for (; DSViter!=DSViter_end; DSViter++ ) {
71  size_t siz = DSViter->size();
72  if (siz > ignoreDetsAboveNClusters_) continue;
73  totals.strip += siz;
74  totals.stripdets++;
75  }
76  }
77  }
78  if (totals.strip > int(maxNrOfCosmicClusters_)) return totals.strip;
79 
80  // get special input for pixel cluster multiplicity filter
82  e.getByToken(token_pc, pixelClusterDSV);
83  if (!pixelClusterDSV.failedToGet()) {
84  const edmNew::DetSetVector<SiPixelCluster> & input = *pixelClusterDSV;
85 
86  if (ignoreDetsAboveNClusters_ == 0) {
87  totals.pixel = input.dataSize();
88  totals.pixeldets = input.size();
89  } else {
90  //loop over detectors
91  totals.pixel = 0;
92  totals.pixeldets = 0;
93  edmNew::DetSetVector<SiPixelCluster>::const_iterator DSViter=input.begin(), DSViter_end=input.end();
94  for (; DSViter!=DSViter_end; DSViter++ ) {
95  size_t siz = DSViter->size();
96  if (siz > ignoreDetsAboveNClusters_) continue;
97  totals.pixel += siz;
98  totals.pixeldets++;
99  }
100  }
101  }
102  else{
103  //say something's wrong.
104  edm::LogError("ClusterChecker")<<"could not get any SiPixel cluster collections of type edm::DetSetVector<SiPixelCluster> with label: "<<pixelClusterCollectionInputTag_;
105  totals.pixel = 999999;
106  }
107  if (totals.pixel > int(maxNrOfPixelClusters_)) return totals.pixel;
108 
109  if (!selector_(totals)) return totals.strip;
110  return 0;
111 }
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
T getParameter(std::string const &) const
int stripdets
number of pixel clusters
int pixeldets
number of strip detectors with at least one cluster
boost::transform_iterator< IterHelp, const_IdIter > const_iterator
const_iterator end(bool update=false) const
size_type dataSize() const
ClusterChecker()=delete
bool existsAs(std::string const &parameterName, bool trackiness=true) const
checks if a parameter exists as a given type
Definition: ParameterSet.h:186
edm::InputTag clusterCollectionInputTag_
StringCutObjectSelector< reco::utils::ClusterTotals > selector_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:508
unsigned int maxNrOfPixelClusters_
static void fillDescriptions(edm::ParameterSetDescription &description)
static std::string const input
Definition: EdmProvDump.cc:44
int pixel
number of strip clusters
unsigned int maxNrOfCosmicClusters_
size_t tooManyClusters(const edm::Event &e) const
unsigned int ignoreDetsAboveNClusters_
ParameterDescriptionBase * add(U const &iLabel, T const &value)
bool failedToGet() const
Definition: HandleBase.h:78
edm::InputTag pixelClusterCollectionInputTag_
edm::EDGetTokenT< edmNew::DetSetVector< SiPixelCluster > > token_pc
size_type size() const
edm::EDGetTokenT< edmNew::DetSetVector< SiStripCluster > > token_sc
const_iterator begin(bool update=false) const