CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
SeedingLayerSetsBuilder.cc
Go to the documentation of this file.
2 
6 
9 
13 
16 
17 #include "HitExtractorPIX.h"
18 #include "HitExtractorSTRP.h"
19 
20 #include <iostream>
21 #include <sstream>
22 #include <ostream>
23 #include <fstream>
24 #include <map>
25 
26 
27 using namespace ctfseeding;
28 using namespace std;
29 
30 namespace {
33  int> nameToEnumId(const std::string& name) {
34  GeomDetEnumerators::SubDetector subdet = GeomDetEnumerators::invalidDet;
35  SeedingLayer::Side side = SeedingLayer::Barrel;
36  int idLayer = 0;
37 
38  size_t index;
39  //
40  // BPIX
41  //
42  if ((index = name.find("BPix")) != string::npos) {
44  side = SeedingLayer::Barrel;
45  idLayer = atoi(name.substr(index+4,1).c_str());
46  }
47  //
48  // FPIX
49  //
50  else if ((index = name.find("FPix")) != string::npos) {
52  idLayer = atoi(name.substr(index+4,1).c_str());
53  if ( name.find("pos") != string::npos ) {
55  } else {
57  }
58  }
59  //
60  // TIB
61  //
62  else if ((index = name.find("TIB")) != string::npos) {
63  subdet = GeomDetEnumerators::TIB;
64  side = SeedingLayer::Barrel;
65  idLayer = atoi(name.substr(index+3,1).c_str());
66  }
67  //
68  // TID
69  //
70  else if ((index = name.find("TID")) != string::npos) {
71  subdet = GeomDetEnumerators::TID;
72  idLayer = atoi(name.substr(index+3,1).c_str());
73  if ( name.find("pos") != string::npos ) {
75  } else {
77  }
78  }
79  //
80  // TOB
81  //
82  else if ((index = name.find("TOB")) != string::npos) {
83  subdet = GeomDetEnumerators::TOB;
84  side = SeedingLayer::Barrel;
85  idLayer = atoi(name.substr(index+3,1).c_str());
86  }
87  //
88  // TEC
89  //
90  else if ((index = name.find("TEC")) != string::npos) {
91  subdet = GeomDetEnumerators::TEC;
92  idLayer = atoi(name.substr(index+3,1).c_str());
93  if ( name.find("pos") != string::npos ) {
95  } else {
97  }
98  }
99  return std::make_tuple(subdet, side, idLayer);
100  }
101 }
102 
103 SeedingLayerSetsBuilder::LayerSpec::LayerSpec(unsigned short index, const std::string& layerName, const edm::ParameterSet& cfgLayer, edm::ConsumesCollector& iC):
104  nameIndex(index),
105  hitBuilder(cfgLayer.getParameter<string>("TTRHBuilder"))
106 {
107  usePixelHitProducer = false;
108  if (cfgLayer.exists("HitProducer")) {
109  pixelHitProducer = cfgLayer.getParameter<string>("HitProducer");
110  usePixelHitProducer = true;
111  }
112 
113  bool skipClusters = cfgLayer.exists("skipClusters");
114  if (skipClusters) {
115  LogDebug("SeedingLayerSetsBuilder")<<layerName<<" ready for skipping";
116  }
117  else{
118  LogDebug("SeedingLayerSetsBuilder")<<layerName<<" not skipping ";
119  }
120 
121  auto subdetData = nameToEnumId(layerName);
122  subdet = std::get<0>(subdetData);
123  side = std::get<1>(subdetData);
124  idLayer = std::get<2>(subdetData);
125  if(subdet == GeomDetEnumerators::PixelBarrel ||
127  extractor = std::make_shared<HitExtractorPIX>(side, idLayer, pixelHitProducer, iC);
128  }
129  else if(subdet != GeomDetEnumerators::invalidDet) { // strip
130  std::shared_ptr<HitExtractorSTRP> extr = std::make_shared<HitExtractorSTRP>(subdet, side, idLayer);
131  if (cfgLayer.exists("matchedRecHits")) {
132  extr->useMatchedHits(cfgLayer.getParameter<edm::InputTag>("matchedRecHits"), iC);
133  }
134  if (cfgLayer.exists("rphiRecHits")) {
135  extr->useRPhiHits(cfgLayer.getParameter<edm::InputTag>("rphiRecHits"), iC);
136  }
137  if (cfgLayer.exists("stereoRecHits")) {
138  extr->useStereoHits(cfgLayer.getParameter<edm::InputTag>("stereoRecHits"), iC);
139  }
140  if (cfgLayer.exists("useRingSlector") && cfgLayer.getParameter<bool>("useRingSlector")) {
141  extr->useRingSelector(cfgLayer.getParameter<int>("minRing"),
142  cfgLayer.getParameter<int>("maxRing"));
143  }
144  bool useSimpleRphiHitsCleaner = cfgLayer.exists("useSimpleRphiHitsCleaner") ? cfgLayer.getParameter<bool>("useSimpleRphiHitsCleaner") : true;
145  extr->useSimpleRphiHitsCleaner(useSimpleRphiHitsCleaner);
146 
147  double minAbsZ = cfgLayer.exists("MinAbsZ") ? cfgLayer.getParameter<double>("MinAbsZ") : 0.;
148  if(minAbsZ > 0.) {
149  extr->setMinAbsZ(minAbsZ);
150  }
151  if(skipClusters) {
152  bool useProjection = cfgLayer.exists("useProjection") ? cfgLayer.getParameter<bool>("useProjection") : false;
153  if(useProjection) {
154  LogDebug("SeedingLayerSetsBuilder")<<layerName<<" will project partially masked matched rechit";
155  }
156  else {
157  extr->setNoProjection();
158  }
159  }
160  extractor = std::move(extr);
161  }
162  if(extractor && skipClusters) {
163  extractor->useSkipClusters(cfgLayer.getParameter<edm::InputTag>("skipClusters"), iC);
164  }
165 }
167 
168 std::string SeedingLayerSetsBuilder::LayerSpec::print(const std::vector<std::string>& names) const
169 {
170  std::ostringstream str;
171  str << "Layer="<<names[nameIndex]<<", hitBldr: "<<hitBuilder;
172 
173  str << ", useRingSelector: ";
174  HitExtractorSTRP *ext = nullptr;
175  if((ext = dynamic_cast<HitExtractorSTRP *>(extractor.get())) &&
176  ext->useRingSelector()) {
177  auto minMaxRing = ext->getMinMaxRing();
178  str <<"true,"<<" Rings: ("<< std::get<0>(minMaxRing) <<","<< std::get<1>(minMaxRing) <<")";
179  } else str<<"false";
180 
181  return str.str();
182 }
183 
186  SeedingLayerSetsBuilder(cfg, iC)
187 {}
189 {
190  std::vector<std::string> namesPset = cfg.getParameter<std::vector<std::string> >("layerList");
191  std::vector<std::vector<std::string> > layerNamesInSets = this->layerNamesInSets(namesPset);
192 
193  // debug printout of layers
194  typedef std::vector<std::string>::const_iterator IS;
195  typedef std::vector<std::vector<std::string> >::const_iterator IT;
196  std::ostringstream str;
197  // The following should not be set to cout
198 // for (IT it = layerNamesInSets.begin(); it != layerNamesInSets.end(); it++) {
199 // str << "SET: ";
200 // for (IS is = it->begin(); is != it->end(); is++) str << *is <<" ";
201 // str << std::endl;
202 // }
203 // std::cout << str.str() << std::endl;
204  if(layerNamesInSets.size() == 0)
206  else
207  theNumberOfLayersInSet = layerNamesInSets[0].size();
208 
209 
210  for (IT it = layerNamesInSets.begin(); it != layerNamesInSets.end(); it++) {
211  if(it->size() != theNumberOfLayersInSet)
212  throw cms::Exception("Configuration") << "Assuming all SeedingLayerSets to have same number of layers. LayerSet " << (it-layerNamesInSets.begin()) << " has " << it->size() << " while 0th has " << theNumberOfLayersInSet;
213  for(const std::string& layerName: *it) {
214  auto found = std::find(theLayerNames.begin(), theLayerNames.end(), layerName);
215  unsigned short layerIndex = 0;
216  if(found != theLayerNames.end()) {
217  layerIndex = found-theLayerNames.begin();
218  }
219  else {
221  throw cms::Exception("Assert") << "Too many layers in " << __FILE__ << ":" << __LINE__ << ", we may have to enlarge the index type from unsigned short to unsigned int";
222  }
223 
224  layerIndex = theLayers.size();
225  theLayers.emplace_back(theLayerNames.size(), layerName, layerConfig(layerName, cfg), iC);
226  theLayerNames.push_back(layerName);
227  }
228  theLayerSetIndices.push_back(layerIndex);
229  }
230  }
231  theLayerDets.resize(theLayers.size());
232  theTTRHBuilders.resize(theLayers.size());
233 
234  // debug printout
235  // The following should not be set to cout
236  //for(const LayerSpec& layer: theLayers) {
237  // std::cout << layer.print(theLayerNames) << std::endl;
238  //}
239 }
240 
242 
244 {
246 
247  for (string::size_type iEnd=nameLayer.size(); iEnd > 0; --iEnd) {
248  string name = nameLayer.substr(0,iEnd);
249  if (cfg.exists(name)) return cfg.getParameter<edm::ParameterSet>(name);
250  }
251  edm::LogError("SeedingLayerSetsBuilder") <<"configuration for layer: "<<nameLayer<<" not found, job will probably crash!";
252  return result;
253 }
254 
255 vector<vector<string> > SeedingLayerSetsBuilder::layerNamesInSets( const vector<string> & namesPSet)
256 {
257  std::vector<std::vector<std::string> > result;
258  for (std::vector<std::string>::const_iterator is=namesPSet.begin(); is < namesPSet.end(); ++is) {
259  vector<std::string> layersInSet;
260  string line = *is;
261  string::size_type pos=0;
262  while (pos != string::npos ) {
263  pos=line.find("+");
264  string layer = line.substr(0,pos);
265  layersInSet.push_back(layer);
266  line=line.substr(pos+1,string::npos);
267  }
268  result.push_back(layersInSet);
269  }
270  return result;
271 }
272 
275  es.get<TrackerRecoGeometryRecord>().get( htracker );
276  const GeometricSearchTracker& tracker = *htracker;
277 
278  const std::vector<BarrelDetLayer const*>& bpx = tracker.barrelLayers();
279  const std::vector<BarrelDetLayer const*>& tib = tracker.tibLayers();
280  const std::vector<BarrelDetLayer const*>& tob = tracker.tobLayers();
281 
282  const std::vector<ForwardDetLayer const*>& fpx_pos = tracker.posForwardLayers();
283  const std::vector<ForwardDetLayer const*>& tid_pos = tracker.posTidLayers();
284  const std::vector<ForwardDetLayer const*>& tec_pos = tracker.posTecLayers();
285 
286  const std::vector<ForwardDetLayer const*>& fpx_neg = tracker.negForwardLayers();
287  const std::vector<ForwardDetLayer const*>& tid_neg = tracker.negTidLayers();
288  const std::vector<ForwardDetLayer const*>& tec_neg = tracker.negTecLayers();
289 
290  for(size_t i=0, n=theLayers.size(); i<n; ++i) {
291  const LayerSpec& layer = theLayers[i];
292  const DetLayer * detLayer = nullptr;
293  int index = layer.idLayer-1;
294 
296  detLayer = bpx[index];
297  }
298  else if (layer.subdet == GeomDetEnumerators::PixelEndcap) {
299  if (layer.side == SeedingLayer::PosEndcap) {
300  detLayer = fpx_pos[index];
301  } else {
302  detLayer = fpx_neg[index];
303  }
304  }
305  else if (layer.subdet == GeomDetEnumerators::TIB) {
306  detLayer = tib[index];
307  }
308  else if (layer.subdet == GeomDetEnumerators::TID) {
309  if (layer.side == SeedingLayer::PosEndcap) {
310  detLayer = tid_pos[index];
311  } else {
312  detLayer = tid_neg[index];
313  }
314  }
315  else if (layer.subdet == GeomDetEnumerators::TOB) {
316  detLayer = tob[index];
317  }
318  else if (layer.subdet == GeomDetEnumerators::TEC) {
319  if (layer.side == SeedingLayer::PosEndcap) {
320  detLayer = tec_pos[index];
321  } else {
322  detLayer = tec_neg[index];
323  }
324  }
325  else {
326  throw cms::Exception("Configuration") << "Did not find DetLayer for layer " << theLayerNames[layer.nameIndex];
327  }
328 
330  es.get<TransientRecHitRecord>().get(layer.hitBuilder, builder);
331 
332  theLayerDets[i] = detLayer;
333  theTTRHBuilders[i] = builder.product();
334  }
335 }
336 
338 {
339  updateEventSetup(es);
340 
341  typedef std::vector<SeedingLayer> Set;
343 
344  for(size_t i=0, n=theLayerSetIndices.size(); i<n; i += theNumberOfLayersInSet) {
345  Set set;
346  for(size_t j=0; j<theNumberOfLayersInSet; ++j) {
347  const unsigned short layerIndex = theLayerSetIndices[i+j];
348  const LayerSpec& layer = theLayers[layerIndex];
349  const DetLayer *detLayer = theLayerDets[layerIndex];
350 
351  set.push_back( SeedingLayer( theLayerNames[layerIndex], layerIndex, detLayer, theTTRHBuilders[layerIndex], layer.extractor.get()));
352  }
353  result.push_back(set);
354  }
355  return result;
356 }
357 
359  // We want to evaluate both in the first invocation (to properly
360  // initialize ESWatcher), and this way we avoid one branch compared
361  // to || (should be tiny effect)
362  return geometryWatcher_.check(es) | trhWatcher_.check(es);
363 }
364 
366 void
368  std::vector<unsigned int> & indices, ctfseeding::SeedingLayer::Hits & hits) const {
369  indices.reserve(theLayers.size());
370  for(unsigned int i=0; i<theLayers.size(); ++i) {
371  // The index of the first hit of this layer
372  indices.push_back(hits.size());
373 
374  // Obtain and copy the hits
375  ctfseeding::SeedingLayer::Hits && tmp = theLayers[i].extractor->hits((const TkTransientTrackingRecHitBuilder &)(*theTTRHBuilders[i]), ev, es);
376  std::move(tmp.begin(), tmp.end(), std::back_inserter(hits));
377  }
378 }
#define LogDebug(id)
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
ctfseeding::SeedingLayerSets layers(const edm::EventSetup &es)
static const HistoName names[]
std::vector< ForwardDetLayer const * > const & posForwardLayers() const
std::vector< LayerSpec > theLayers
edm::ESWatcher< TrackerRecoGeometryRecord > geometryWatcher_
LayerSpec(unsigned short index, const std::string &layerName, const edm::ParameterSet &cfgLayer, edm::ConsumesCollector &iC)
std::vector< std::string > theLayerNames
bool exists(std::string const &parameterName) const
checks if a parameter exists
std::vector< HitPointer > Hits
Definition: SeedingLayer.h:28
std::vector< BarrelDetLayer const * > const & tobLayers() const
bool check(const edm::EventSetup &es)
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:7
uint16_t size_type
std::vector< const TransientTrackingRecHitBuilder * > theTTRHBuilders
GeomDetEnumerators::SubDetector subdet
const T & max(const T &a, const T &b)
void updateEventSetup(const edm::EventSetup &es)
tuple result
Definition: query.py:137
std::vector< ForwardDetLayer const * > const & negForwardLayers() const
std::shared_ptr< ctfseeding::HitExtractor > extractor
int j
Definition: DBlmapReader.cc:9
std::vector< LinkConnSpec >::const_iterator IT
edm::ESWatcher< TransientRecHitRecord > trhWatcher_
std::vector< LayerSetIndex > theLayerSetIndices
void useRingSelector(int minRing, int maxRing)
std::vector< BarrelDetLayer const * > const & tibLayers() const
std::vector< ForwardDetLayer const * > const & posTecLayers() const
const T & get() const
Definition: EventSetup.h:55
T const * product() const
Definition: ESHandle.h:86
std::vector< ForwardDetLayer const * > const & negTidLayers() const
bool check(const edm::EventSetup &iSetup)
Definition: ESWatcher.h:57
std::vector< std::vector< double > > tmp
Definition: MVATrainer.cc:100
std::vector< ForwardDetLayer const * > const & posTidLayers() const
std::string print(const std::vector< std::string > &names) const
void hits(const edm::Event &ev, const edm::EventSetup &es, std::vector< unsigned int > &indices, ctfseeding::SeedingLayer::Hits &hits) const
std::vector< ForwardDetLayer const * > const & negTecLayers() const
std::tuple< int, int > getMinMaxRing() const
edm::ParameterSet layerConfig(const std::string &nameLayer, const edm::ParameterSet &cfg) const
std::vector< const DetLayer * > theLayerDets
std::vector< BarrelDetLayer const * > const & barrelLayers() const
ctfseeding::SeedingLayer::Side side
std::vector< std::vector< SeedingLayer > > SeedingLayerSets
std::vector< std::vector< std::string > > layerNamesInSets(const std::vector< std::string > &namesPSet)