CMS 3D CMS Logo

StubKiller.cc
Go to the documentation of this file.
2 
3 using namespace std;
4 
6  : killScenario_(0),
7  trackerTopology_(nullptr),
8  trackerGeometry_(nullptr),
9  layersToKill_(vector<int>()),
10  minPhiToKill_(0),
11  maxPhiToKill_(0),
12  minZToKill_(0),
13  maxZToKill_(0),
14  minRToKill_(0),
15  maxRToKill_(0),
16  fractionOfStubsToKillInLayers_(0),
17  fractionOfStubsToKillEverywhere_(0),
18  fractionOfModulesToKillEverywhere_(0) {}
19 
20 void StubKiller::initialise(unsigned int killScenario,
21  const TrackerTopology* trackerTopology,
22  const TrackerGeometry* trackerGeometry) {
23  killScenario_ = killScenario;
24  trackerTopology_ = trackerTopology;
25  trackerGeometry_ = trackerGeometry;
26 
27  switch (killScenario_) {
28  // kill layer 5 in one quadrant + 5% random module loss to connect to what was done before
29  case 1:
30  layersToKill_ = {5};
31  minPhiToKill_ = 0.0;
32  maxPhiToKill_ = TMath::PiOver2();
33  minZToKill_ = -1000.0;
34  maxZToKill_ = 0.0;
35  minRToKill_ = 0.0;
36  maxRToKill_ = 1000.0;
40  break;
41 
42  // kill layer 1 in one quadrant + 5% random module loss
43  case 2:
44  layersToKill_ = {1};
45  minPhiToKill_ = 0.0;
46  maxPhiToKill_ = TMath::PiOver2();
47  minZToKill_ = -1000.0;
48  maxZToKill_ = 0.0;
49  minRToKill_ = 0.0;
50  maxRToKill_ = 1000.0;
54  break;
55 
56  // kill layer 1 + layer 2, both in same quadrant
57  case 3:
58  layersToKill_ = {1, 2};
59  minPhiToKill_ = 0.0;
60  maxPhiToKill_ = TMath::PiOver2();
61  minZToKill_ = -1000.0;
62  maxZToKill_ = 0.0;
63  minRToKill_ = 0.0;
64  maxRToKill_ = 1000.0;
68  break;
69 
70  // kill layer 1 and disk 1, both in same quadrant
71  case 4:
72  layersToKill_ = {1, 11};
73  minPhiToKill_ = 0.0;
74  maxPhiToKill_ = TMath::PiOver2();
75  minZToKill_ = -1000.0;
76  maxZToKill_ = 0.0;
77  minRToKill_ = 0.0;
78  maxRToKill_ = 66.5;
82  break;
83 
84  // 5% random module loss throughout tracker
85  case 5:
86  layersToKill_ = {};
90  break;
91 
92  // 1% random module loss throughout tracker
93  case 6:
94  layersToKill_ = {};
98  break;
99 
100  // kill layer 5 in one quadrant + 1% random module loss
101  case 7:
102  layersToKill_ = {5};
103  minPhiToKill_ = 0.0;
104  maxPhiToKill_ = TMath::PiOver2();
105  minZToKill_ = -1000.0;
106  maxZToKill_ = 0.0;
107  minRToKill_ = 0.0;
108  maxRToKill_ = 1000.0;
112  break;
113 
114  // kill layer 1 in one quadrant +1 % random module loss
115  case 8:
116  layersToKill_ = {1};
117  minPhiToKill_ = 0.0;
118  maxPhiToKill_ = TMath::PiOver2();
119  minZToKill_ = -1000.0;
120  maxZToKill_ = 0.0;
121  minRToKill_ = 0.0;
122  maxRToKill_ = 1000.0;
126  break;
127 
128  // 10% random module loss throughout tracker
129  case 9:
130  layersToKill_ = {};
134  break;
135  }
136 
137  deadModules_.clear();
139  this->chooseModulesToKill();
140  }
142 }
143 
145  TRandom3* randomGenerator = new TRandom3();
146  randomGenerator->SetSeed(0);
147 
148  for (const GeomDetUnit* gd : trackerGeometry_->detUnits()) {
149  if (!trackerTopology_->isLower(gd->geographicalId()))
150  continue;
151  if (randomGenerator->Uniform(0.0, 1.0) < fractionOfModulesToKillEverywhere_) {
152  deadModules_[gd->geographicalId()] = 1;
153  }
154  }
155 }
156 
158  for (const GeomDetUnit* gd : trackerGeometry_->detUnits()) {
159  float moduleR = gd->position().perp();
160  float moduleZ = gd->position().z();
161  float modulePhi = gd->position().phi();
162  DetId geoDetId = gd->geographicalId();
163  bool isInBarrel = geoDetId.subdetId() == StripSubdetector::TOB || geoDetId.subdetId() == StripSubdetector::TIB;
164 
165  int layerID = 0;
166  if (isInBarrel) {
167  layerID = trackerTopology_->layer(geoDetId);
168  } else {
169  layerID = 10 * trackerTopology_->side(geoDetId) + trackerTopology_->tidWheel(geoDetId);
170  }
171 
172  if (find(layersToKill_.begin(), layersToKill_.end(), layerID) != layersToKill_.end()) {
173  if (modulePhi < -1.0 * TMath::Pi())
174  modulePhi += 2.0 * TMath::Pi();
175  else if (modulePhi > TMath::Pi())
176  modulePhi -= 2.0 * TMath::Pi();
177 
178  if (modulePhi > minPhiToKill_ && modulePhi < maxPhiToKill_ && moduleZ > minZToKill_ && moduleZ < maxZToKill_ &&
179  moduleR > minRToKill_ && moduleR < maxRToKill_) {
180  if (deadModules_.find(gd->geographicalId()) == deadModules_.end()) {
181  deadModules_[gd->geographicalId()] = fractionOfStubsToKillInLayers_;
182  }
183  }
184  }
185  }
186 }
187 
189  if (killScenario_ == 0)
190  return false;
191  else {
192  bool killStubRandomly = killStub(stub,
196  minZToKill_,
197  maxZToKill_,
198  minRToKill_,
199  maxRToKill_,
202  bool killStubInDeadModules = killStubInDeadModule(stub);
203  return killStubRandomly || killStubInDeadModules;
204  }
205 }
206 
207 // layersToKill - a vector stating the layers we are killing stubs in. Can be an empty vector.
208 // Barrel layers are encoded as 1-6. The endcap layers are encoded as 11-15 (-z) and 21-25 (+z)
209 // min/max Phi/Z/R - stubs within the region specified by these boundaries and layersToKill are flagged for killing
210 // fractionOfStubsToKillInLayers - The fraction of stubs to kill in the specified layers/region.
211 // fractionOfStubsToKillEverywhere - The fraction of stubs to kill throughout the tracker
213  const vector<int> layersToKill,
214  const double minPhiToKill,
215  const double maxPhiToKill,
216  const double minZToKill,
217  const double maxZToKill,
218  const double minRToKill,
219  const double maxRToKill,
220  const double fractionOfStubsToKillInLayers,
221  const double fractionOfStubsToKillEverywhere) {
222  // Only kill stubs in specified layers
223  if (!layersToKill.empty()) {
224  // Get the layer the stub is in, and check if it's in the layer you want to kill
225  DetId stackDetid = stub->getDetId();
226  DetId geoDetId(stackDetid.rawId() + 1);
227 
228  bool isInBarrel = geoDetId.subdetId() == StripSubdetector::TOB || geoDetId.subdetId() == StripSubdetector::TIB;
229 
230  int layerID = 0;
231  if (isInBarrel) {
232  layerID = trackerTopology_->layer(geoDetId);
233  } else {
234  layerID = 10 * trackerTopology_->side(geoDetId) + trackerTopology_->tidWheel(geoDetId);
235  }
236 
237  if (find(layersToKill.begin(), layersToKill.end(), layerID) != layersToKill.end()) {
238  // Get the phi and z of stub, and check if it's in the region you want to kill
239  const GeomDetUnit* det0 = trackerGeometry_->idToDetUnit(geoDetId);
240  const PixelGeomDetUnit* theGeomDet = dynamic_cast<const PixelGeomDetUnit*>(det0);
241  const PixelTopology* topol = dynamic_cast<const PixelTopology*>(&(theGeomDet->specificTopology()));
242  MeasurementPoint measurementPoint = stub->clusterRef(0)->findAverageLocalCoordinatesCentered();
243  LocalPoint clustlp = topol->localPosition(measurementPoint);
244  GlobalPoint pos = theGeomDet->surface().toGlobal(clustlp);
245 
246  // Just in case phi is outside of -pi -> pi
247  double stubPhi = pos.phi();
248  if (stubPhi < -1.0 * TMath::Pi())
249  stubPhi += 2.0 * TMath::Pi();
250  else if (stubPhi > TMath::Pi())
251  stubPhi -= 2.0 * TMath::Pi();
252 
253  if (stubPhi > minPhiToKill && stubPhi < maxPhiToKill && pos.z() > minZToKill && pos.z() < maxZToKill &&
254  pos.perp() > minRToKill && pos.perp() < maxRToKill) {
255  // Kill fraction of stubs
256  if (fractionOfStubsToKillInLayers == 1) {
257  return true;
258  } else {
259  static TRandom randomGenerator;
260  if (randomGenerator.Rndm() < fractionOfStubsToKillInLayers) {
261  return true;
262  }
263  }
264  }
265  }
266  }
267 
268  // Kill fraction of stubs throughout tracker
269  if (fractionOfStubsToKillEverywhere > 0) {
270  static TRandom randomGenerator;
271  if (randomGenerator.Rndm() < fractionOfStubsToKillEverywhere) {
272  return true;
273  }
274  }
275 
276  return false;
277 }
278 
280  if (!deadModules_.empty()) {
281  DetId stackDetid = stub->getDetId();
282  DetId geoDetId(stackDetid.rawId() + 1);
283  if (deadModules_.find(geoDetId) != deadModules_.end())
284  return true;
285  }
286 
287  return false;
288 }
const edm::Ref< edmNew::DetSetVector< TTCluster< T > >, TTCluster< T > > & clusterRef(unsigned int hitStackMember) const
Clusters composing the Stub – see https://twiki.cern.ch/twiki/bin/viewauth/CMS/SLHCTrackerTriggerSWT...
Definition: TTStub.h:150
const double Pi
std::vector< int > layersToKill_
Definition: StubKiller.h:46
const TrackerGeomDet * idToDetUnit(DetId) const override
Return the pointer to the GeomDetUnit corresponding to a given DetId.
virtual LocalPoint localPosition(const MeasurementPoint &) const =0
bool killStubInDeadModule(const TTStub< Ref_Phase2TrackerDigi_ > *stub)
Definition: StubKiller.cc:279
std::map< DetId, float > deadModules_
Definition: StubKiller.h:57
unsigned int tidWheel(const DetId &id) const
double minPhiToKill_
Definition: StubKiller.h:47
const DetContainer & detUnits() const override
Returm a vector of all GeomDet.
unsigned int side(const DetId &id) const
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
DetId getDetId() const
Detector element.
Definition: TTStub.h:44
double fractionOfStubsToKillInLayers_
Definition: StubKiller.h:53
unsigned int layer(const DetId &id) const
double minRToKill_
Definition: StubKiller.h:51
void initialise(unsigned int killScenario, const TrackerTopology *trackerTopology, const TrackerGeometry *trackerGeometry)
Definition: StubKiller.cc:20
bool killStub(const TTStub< Ref_Phase2TrackerDigi_ > *stub, const std::vector< int > layersToKill, const double minPhiToKill, const double maxPhiToKill, const double minZToKill, const double maxZToKill, const double minRToKill, const double maxRToKill, const double fractionOfStubsToKillInLayers, const double fractionOfStubsToKillEverywhere)
Definition: StubKiller.cc:212
const TrackerGeometry * trackerGeometry_
Definition: StubKiller.h:44
constexpr int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:48
static constexpr auto TOB
Class to store the L1 Track Trigger stubs.
Definition: TTStub.h:22
double minZToKill_
Definition: StubKiller.h:49
GlobalPoint toGlobal(const Point2DBase< Scalar, LocalTag > lp) const
Definition: Surface.h:79
Definition: DetId.h:17
bool isLower(const DetId &id) const
static constexpr auto TIB
const Plane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:37
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
unsigned int killScenario_
Definition: StubKiller.h:42
double maxRToKill_
Definition: StubKiller.h:52
void chooseModulesToKill()
Definition: StubKiller.cc:144
virtual const PixelTopology & specificTopology() const
Returns a reference to the pixel proxy topology.
double fractionOfStubsToKillEverywhere_
Definition: StubKiller.h:54
void addDeadLayerModulesToDeadModuleList()
Definition: StubKiller.cc:157
double maxPhiToKill_
Definition: StubKiller.h:48
double maxZToKill_
Definition: StubKiller.h:50
double fractionOfModulesToKillEverywhere_
Definition: StubKiller.h:55
const TrackerTopology * trackerTopology_
Definition: StubKiller.h:43