CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
List of all members | Public Types | Public Member Functions | Private Member Functions | Private Attributes
tmtt::StubKiller Class Reference

#include <StubKiller.h>

Public Types

enum  KillOptions {
  KillOptions::none = 0, KillOptions::layer5 = 1, KillOptions::layer1 = 2, KillOptions::layer1layer2 = 3,
  KillOptions::layer1disk1 = 4, KillOptions::random = 5
}
 

Public Member Functions

bool killStub (const TTStub< Ref_Phase2TrackerDigi_ > *stub) const
 
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) const
 
bool killStubInDeadModule (const TTStub< Ref_Phase2TrackerDigi_ > *stub) const
 
const std::map< DetId, float > & listOfDeadModules () const
 
 StubKiller (KillOptions killScenario, const TrackerTopology *trackerTopology, const TrackerGeometry *trackerGeometry, const edm::Event &iEvent)
 

Private Member Functions

void addDeadLayerModulesToDeadModuleList ()
 
void chooseModulesToKill ()
 

Private Attributes

std::map< DetId, float > deadModules_
 
double fractionOfModulesToKillEverywhere_
 
double fractionOfStubsToKillEverywhere_
 
double fractionOfStubsToKillInLayers_
 
KillOptions killScenario_
 
std::vector< int > layersToKill_
 
double maxPhiToKill_
 
double maxRToKill_
 
double maxZToKill_
 
double minPhiToKill_
 
double minRToKill_
 
double minZToKill_
 
CLHEP::HepRandomEngine * rndmEngine_
 
edm::Service
< edm::RandomNumberGenerator
rndmService_
 
const TrackerGeometrytrackerGeometry_
 
const TrackerTopologytrackerTopology_
 

Detailed Description

Definition at line 20 of file StubKiller.h.

Member Enumeration Documentation

Enumerator
none 
layer5 
layer1 
layer1layer2 
layer1disk1 
random 

Definition at line 22 of file StubKiller.h.

22 { none = 0, layer5 = 1, layer1 = 2, layer1layer2 = 3, layer1disk1 = 4, random = 5 };

Constructor & Destructor Documentation

tmtt::StubKiller::StubKiller ( StubKiller::KillOptions  killScenario,
const TrackerTopology trackerTopology,
const TrackerGeometry trackerGeometry,
const edm::Event iEvent 
)

Definition at line 9 of file StubKiller.cc.

References addDeadLayerModulesToDeadModuleList(), chooseModulesToKill(), deadModules_, Exception, fractionOfModulesToKillEverywhere_, fractionOfStubsToKillEverywhere_, fractionOfStubsToKillInLayers_, edm::RandomNumberGenerator::getEngine(), edm::Service< T >::isAvailable(), killScenario_, layer1, layer1disk1, layer1layer2, layer5, layersToKill_, M_PI, maxPhiToKill_, maxRToKill_, maxZToKill_, minPhiToKill_, minRToKill_, minZToKill_, random, rndmEngine_, rndmService_, and edm::Event::streamID().

13  : killScenario_(killScenario),
14  trackerTopology_(trackerTopology),
15  trackerGeometry_(trackerGeometry),
16  minPhiToKill_(0),
17  maxPhiToKill_(0),
18  minZToKill_(0),
19  maxZToKill_(0),
20  minRToKill_(0),
21  maxRToKill_(0),
25  if (rndmService_.isAvailable()) {
26  rndmEngine_ = &(rndmService_->getEngine(event.streamID()));
27  } else {
28  throw cms::Exception("BadConfig")
29  << "StubKiller: requires RandomNumberGeneratorService, not present in cfg file, namely:" << endl
30  << "process.RandomNumberGeneratorService=cms.Service('RandomNumberGeneratorService',TMTrackProducer=cms.PSet("
31  "initialSeed=cms.untracked.uint32(12345)))";
32  }
33 
34  // These scenarios correspond to slide 12 of https://indico.cern.ch/event/719985/contributions/2970687/attachments/1634587/2607365/StressTestTF-Acosta-Apr18.pdf
35  // Scenario 1
36 
37  // kill layer 5 in one quadrant +5 % random module loss to connect to what was done before
39  layersToKill_ = {5};
40  minPhiToKill_ = 0;
41  maxPhiToKill_ = 0.5 * M_PI;
42  minZToKill_ = -1000;
43  maxZToKill_ = 0;
44  minRToKill_ = 0;
45  maxRToKill_ = 1000;
49  }
50  // Scenario 2
51  // kill layer 1 in one quadrant +5 % random module loss
52  else if (killScenario_ == KillOptions::layer1) {
53  layersToKill_ = {1};
54  minPhiToKill_ = 0;
55  maxPhiToKill_ = 0.5 * M_PI;
56  minZToKill_ = -1000;
57  maxZToKill_ = 0;
58  minRToKill_ = 0;
59  maxRToKill_ = 1000;
63  }
64  // Scenario 3
65  // kill layer 1 + layer 2, both in same quadrant
67  layersToKill_ = {1, 2};
68  minPhiToKill_ = 0;
69  maxPhiToKill_ = 0.5 * M_PI;
70  minZToKill_ = -1000;
71  maxZToKill_ = 0;
72  minRToKill_ = 0;
73  maxRToKill_ = 1000;
77  }
78  // Scenario 4
79  // kill layer 1 and disk 1, both in same quadrant
81  layersToKill_ = {1, 11};
82  minPhiToKill_ = 0;
83  maxPhiToKill_ = 0.5 * M_PI;
84  minZToKill_ = -1000;
85  maxZToKill_ = 0;
86  minRToKill_ = 0;
87  maxRToKill_ = 66.5;
91  }
92  // An extra scenario not listed in the slides
93  // 5% random module loss throughout tracker
94  else if (killScenario_ == KillOptions::random) {
95  layersToKill_ = {};
99  }
100 
101  deadModules_.clear();
103  this->chooseModulesToKill();
104  }
106  }
double minPhiToKill_
Definition: StubKiller.h:62
const TrackerGeometry * trackerGeometry_
Definition: StubKiller.h:59
double fractionOfModulesToKillEverywhere_
Definition: StubKiller.h:70
double maxPhiToKill_
Definition: StubKiller.h:63
KillOptions killScenario_
Definition: StubKiller.h:57
CLHEP::HepRandomEngine * rndmEngine_
Definition: StubKiller.h:75
double fractionOfStubsToKillInLayers_
Definition: StubKiller.h:68
virtual CLHEP::HepRandomEngine & getEngine(StreamID const &)=0
Use this engine in event methods.
std::map< DetId, float > deadModules_
Definition: StubKiller.h:72
void addDeadLayerModulesToDeadModuleList()
Definition: StubKiller.cc:239
std::vector< int > layersToKill_
Definition: StubKiller.h:61
bool isAvailable() const
Definition: Service.h:40
double minRToKill_
Definition: StubKiller.h:66
#define M_PI
const TrackerTopology * trackerTopology_
Definition: StubKiller.h:58
double fractionOfStubsToKillEverywhere_
Definition: StubKiller.h:69
void chooseModulesToKill()
Definition: StubKiller.cc:227
double minZToKill_
Definition: StubKiller.h:64
edm::Service< edm::RandomNumberGenerator > rndmService_
Definition: StubKiller.h:74
double maxZToKill_
Definition: StubKiller.h:65
double maxRToKill_
Definition: StubKiller.h:67

Member Function Documentation

void tmtt::StubKiller::addDeadLayerModulesToDeadModuleList ( )
private

Definition at line 239 of file StubKiller.cc.

References deadModules_, reco::deltaPhi(), TrackerGeometry::detUnits(), spr::find(), fractionOfStubsToKillInLayers_, TrackerTopology::layer(), layersToKill_, maxRToKill_, maxZToKill_, minPhiToKill_, minRToKill_, minZToKill_, TrackerTopology::side(), DetId::subdetId(), StripSubdetector::TIB, TrackerTopology::tidWheel(), StripSubdetector::TOB, trackerGeometry_, and trackerTopology_.

Referenced by StubKiller().

239  {
240  for (const GeomDetUnit* gd : trackerGeometry_->detUnits()) {
241  float moduleR = gd->position().perp();
242  float moduleZ = gd->position().z();
243  float modulePhi = reco::deltaPhi(gd->position().phi(), 0.);
244  DetId geoDetId = gd->geographicalId();
245  bool isInBarrel = geoDetId.subdetId() == StripSubdetector::TOB || geoDetId.subdetId() == StripSubdetector::TIB;
246 
247  int layerID = 0;
248  if (isInBarrel) {
249  layerID = trackerTopology_->layer(geoDetId);
250  } else {
251  layerID = 10 * trackerTopology_->side(geoDetId) + trackerTopology_->tidWheel(geoDetId);
252  }
253  if (find(layersToKill_.begin(), layersToKill_.end(), layerID) != layersToKill_.end()) {
254  if (modulePhi > minPhiToKill_ && modulePhi < maxPhiToKill_ && moduleZ > minZToKill_ && moduleZ < maxZToKill_ &&
255  moduleR > minRToKill_ && moduleR < maxRToKill_) {
256  if (deadModules_.find(gd->geographicalId()) == deadModules_.end()) {
257  deadModules_[gd->geographicalId()] = fractionOfStubsToKillInLayers_;
258  }
259  }
260  }
261  }
262  }
constexpr double deltaPhi(double phi1, double phi2)
Definition: deltaPhi.h:26
double minPhiToKill_
Definition: StubKiller.h:62
const TrackerGeometry * trackerGeometry_
Definition: StubKiller.h:59
unsigned int side(const DetId &id) const
unsigned int tidWheel(const DetId &id) const
const DetContainer & detUnits() const override
Returm a vector of all GeomDet.
double fractionOfStubsToKillInLayers_
Definition: StubKiller.h:68
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
std::map< DetId, float > deadModules_
Definition: StubKiller.h:72
std::vector< int > layersToKill_
Definition: StubKiller.h:61
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
double minRToKill_
Definition: StubKiller.h:66
Definition: DetId.h:17
const TrackerTopology * trackerTopology_
Definition: StubKiller.h:58
static constexpr auto TIB
unsigned int layer(const DetId &id) const
double minZToKill_
Definition: StubKiller.h:64
double maxZToKill_
Definition: StubKiller.h:65
double maxRToKill_
Definition: StubKiller.h:67
void tmtt::StubKiller::chooseModulesToKill ( )
private

Definition at line 227 of file StubKiller.cc.

References deadModules_, TrackerGeometry::detUnits(), fractionOfModulesToKillEverywhere_, TrackerTopology::isLower(), rndmEngine_, trackerGeometry_, and trackerTopology_.

Referenced by StubKiller().

227  {
228  for (const GeomDetUnit* gd : trackerGeometry_->detUnits()) {
229  if (!trackerTopology_->isLower(gd->geographicalId()))
230  continue;
232  deadModules_[gd->geographicalId()] = 1;
233  }
234  }
235  }
const TrackerGeometry * trackerGeometry_
Definition: StubKiller.h:59
double fractionOfModulesToKillEverywhere_
Definition: StubKiller.h:70
CLHEP::HepRandomEngine * rndmEngine_
Definition: StubKiller.h:75
const DetContainer & detUnits() const override
Returm a vector of all GeomDet.
bool isLower(const DetId &id) const
std::map< DetId, float > deadModules_
Definition: StubKiller.h:72
const TrackerTopology * trackerTopology_
Definition: StubKiller.h:58
bool tmtt::StubKiller::killStub ( const TTStub< Ref_Phase2TrackerDigi_ > *  stub) const

Definition at line 110 of file StubKiller.cc.

References fractionOfStubsToKillEverywhere_, fractionOfStubsToKillInLayers_, killScenario_, killStubInDeadModule(), layersToKill_, maxPhiToKill_, maxRToKill_, maxZToKill_, minPhiToKill_, minRToKill_, minZToKill_, and none.

Referenced by tmtt::Stub::setFrontend().

110  {
112  return false;
113  else {
114  // Check if stub is in dead region specified by *ToKill_
115  // And randomly kill stubs throughout tracker (not just thos in specific regions/modules)
116  bool killStubRandomly = killStub(stub,
120  minZToKill_,
121  maxZToKill_,
122  minRToKill_,
123  maxRToKill_,
126  // Kill modules in specifid modules
127  // Random modules throughout the tracker, and those modules in specific regions (so may already have been killed by killStub above)
128  bool killStubInDeadModules = killStubInDeadModule(stub);
129  return killStubRandomly || killStubInDeadModules;
130  }
131  }
bool killStub(const TTStub< Ref_Phase2TrackerDigi_ > *stub) const
Definition: StubKiller.cc:110
double minPhiToKill_
Definition: StubKiller.h:62
double maxPhiToKill_
Definition: StubKiller.h:63
KillOptions killScenario_
Definition: StubKiller.h:57
bool killStubInDeadModule(const TTStub< Ref_Phase2TrackerDigi_ > *stub) const
Definition: StubKiller.cc:206
double fractionOfStubsToKillInLayers_
Definition: StubKiller.h:68
std::vector< int > layersToKill_
Definition: StubKiller.h:61
double minRToKill_
Definition: StubKiller.h:66
double fractionOfStubsToKillEverywhere_
Definition: StubKiller.h:69
double minZToKill_
Definition: StubKiller.h:64
double maxZToKill_
Definition: StubKiller.h:65
double maxRToKill_
Definition: StubKiller.h:67
bool tmtt::StubKiller::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 
) const

Definition at line 141 of file StubKiller.cc.

References TTStub< T >::clusterRef(), deadModules_, reco::deltaPhi(), spr::find(), TTStub< T >::getDetId(), TrackerGeometry::idToDetUnit(), TrackerTopology::layer(), Topology::localPosition(), PV3DBase< T, PVType, FrameType >::perp(), PV3DBase< T, PVType, FrameType >::phi(), DetId::rawId(), rndmEngine_, TrackerTopology::side(), PixelGeomDetUnit::specificTopology(), DetId::subdetId(), GeomDet::surface(), StripSubdetector::TIB, TrackerTopology::tidWheel(), StripSubdetector::TOB, Surface::toGlobal(), trackerGeometry_, trackerTopology_, and PV3DBase< T, PVType, FrameType >::z().

150  {
151  // Only kill stubs in specified layers
152  if (not layersToKill.empty()) {
153  // Get the layer the stub is in, and check if it's in the layer you want to kill
154  DetId stackDetid = stub->getDetId();
155  DetId geoDetId(stackDetid.rawId() + 1);
156 
157  // If this module is in the deadModule list, don't also try to kill the stub here
158  if (deadModules_.empty() || deadModules_.find(geoDetId) == deadModules_.end()) {
159  bool isInBarrel = geoDetId.subdetId() == StripSubdetector::TOB || geoDetId.subdetId() == StripSubdetector::TIB;
160 
161  int layerID = 0;
162  if (isInBarrel) {
163  layerID = trackerTopology_->layer(geoDetId);
164  } else {
165  layerID = 10 * trackerTopology_->side(geoDetId) + trackerTopology_->tidWheel(geoDetId);
166  }
167 
168  if (find(layersToKill.begin(), layersToKill.end(), layerID) != layersToKill.end()) {
169  // Get the phi and z of stub, and check if it's in the region you want to kill
170  const GeomDetUnit* det0 = trackerGeometry_->idToDetUnit(geoDetId);
171  const PixelGeomDetUnit* theGeomDet = dynamic_cast<const PixelGeomDetUnit*>(det0);
172  const PixelTopology* topol = dynamic_cast<const PixelTopology*>(&(theGeomDet->specificTopology()));
173  MeasurementPoint measurementPoint = stub->clusterRef(0)->findAverageLocalCoordinatesCentered();
174  LocalPoint clustlp = topol->localPosition(measurementPoint);
175  GlobalPoint pos = theGeomDet->surface().toGlobal(clustlp);
176 
177  double stubPhi = reco::deltaPhi(pos.phi(), 0.);
178 
179  if (stubPhi > minPhiToKill && stubPhi < maxPhiToKill && pos.z() > minZToKill && pos.z() < maxZToKill &&
180  pos.perp() > minRToKill && pos.perp() < maxRToKill) {
181  // Kill fraction of stubs
182  if (fractionOfStubsToKillInLayers == 1) {
183  return true;
184  } else {
185  if (rndmEngine_->flat() < fractionOfStubsToKillInLayers) {
186  return true;
187  }
188  }
189  }
190  }
191  }
192  }
193 
194  // Kill fraction of stubs throughout tracker
195  if (fractionOfStubsToKillEverywhere > 0) {
196  if (rndmEngine_->flat() < fractionOfStubsToKillEverywhere) {
197  return true;
198  }
199  }
200 
201  return false;
202  }
constexpr double deltaPhi(double phi1, double phi2)
Definition: deltaPhi.h:26
GlobalPoint toGlobal(const Point2DBase< Scalar, LocalTag > lp) const
Definition: Surface.h:79
const TrackerGeometry * trackerGeometry_
Definition: StubKiller.h:59
T perp() const
Definition: PV3DBase.h:69
const TrackerGeomDet * idToDetUnit(DetId) const override
Return the pointer to the GeomDetUnit corresponding to a given DetId.
virtual LocalPoint localPosition(const MeasurementPoint &) const =0
CLHEP::HepRandomEngine * rndmEngine_
Definition: StubKiller.h:75
Geom::Phi< T > phi() const
Definition: PV3DBase.h:66
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
unsigned int side(const DetId &id) const
unsigned int tidWheel(const DetId &id) const
DetId getDetId() const
Detector element.
Definition: TTStub.h:44
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
const Plane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:37
std::map< DetId, float > deadModules_
Definition: StubKiller.h:72
T z() const
Definition: PV3DBase.h:61
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
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/SLHCTrackerTriggerSWTo...
Definition: TTStub.h:150
Definition: DetId.h:17
const TrackerTopology * trackerTopology_
Definition: StubKiller.h:58
static constexpr auto TIB
virtual const PixelTopology & specificTopology() const
Returns a reference to the pixel proxy topology.
unsigned int layer(const DetId &id) const
bool tmtt::StubKiller::killStubInDeadModule ( const TTStub< Ref_Phase2TrackerDigi_ > *  stub) const

Definition at line 206 of file StubKiller.cc.

References deadModules_, TTStub< T >::getDetId(), DetId::rawId(), and rndmEngine_.

Referenced by killStub().

206  {
207  if (not deadModules_.empty()) {
208  DetId stackDetid = stub->getDetId();
209  DetId geoDetId(stackDetid.rawId() + 1);
210  auto deadModule = deadModules_.find(geoDetId);
211  if (deadModule != deadModules_.end()) {
212  if (deadModule->second == 1) {
213  return true;
214  } else {
215  if (rndmEngine_->flat() < deadModule->second) {
216  return true;
217  }
218  }
219  }
220  }
221 
222  return false;
223  }
CLHEP::HepRandomEngine * rndmEngine_
Definition: StubKiller.h:75
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
DetId getDetId() const
Detector element.
Definition: TTStub.h:44
std::map< DetId, float > deadModules_
Definition: StubKiller.h:72
Definition: DetId.h:17
const std::map<DetId, float>& tmtt::StubKiller::listOfDeadModules ( ) const
inline

Definition at line 49 of file StubKiller.h.

References deadModules_.

49 { return deadModules_; }
std::map< DetId, float > deadModules_
Definition: StubKiller.h:72

Member Data Documentation

std::map<DetId, float> tmtt::StubKiller::deadModules_
private
double tmtt::StubKiller::fractionOfModulesToKillEverywhere_
private

Definition at line 70 of file StubKiller.h.

Referenced by chooseModulesToKill(), and StubKiller().

double tmtt::StubKiller::fractionOfStubsToKillEverywhere_
private

Definition at line 69 of file StubKiller.h.

Referenced by killStub(), and StubKiller().

double tmtt::StubKiller::fractionOfStubsToKillInLayers_
private

Definition at line 68 of file StubKiller.h.

Referenced by addDeadLayerModulesToDeadModuleList(), killStub(), and StubKiller().

KillOptions tmtt::StubKiller::killScenario_
private

Definition at line 57 of file StubKiller.h.

Referenced by killStub(), and StubKiller().

std::vector<int> tmtt::StubKiller::layersToKill_
private

Definition at line 61 of file StubKiller.h.

Referenced by addDeadLayerModulesToDeadModuleList(), killStub(), and StubKiller().

double tmtt::StubKiller::maxPhiToKill_
private

Definition at line 63 of file StubKiller.h.

Referenced by killStub(), and StubKiller().

double tmtt::StubKiller::maxRToKill_
private

Definition at line 67 of file StubKiller.h.

Referenced by addDeadLayerModulesToDeadModuleList(), killStub(), and StubKiller().

double tmtt::StubKiller::maxZToKill_
private

Definition at line 65 of file StubKiller.h.

Referenced by addDeadLayerModulesToDeadModuleList(), killStub(), and StubKiller().

double tmtt::StubKiller::minPhiToKill_
private

Definition at line 62 of file StubKiller.h.

Referenced by addDeadLayerModulesToDeadModuleList(), killStub(), and StubKiller().

double tmtt::StubKiller::minRToKill_
private

Definition at line 66 of file StubKiller.h.

Referenced by addDeadLayerModulesToDeadModuleList(), killStub(), and StubKiller().

double tmtt::StubKiller::minZToKill_
private

Definition at line 64 of file StubKiller.h.

Referenced by addDeadLayerModulesToDeadModuleList(), killStub(), and StubKiller().

CLHEP::HepRandomEngine* tmtt::StubKiller::rndmEngine_
private

Definition at line 75 of file StubKiller.h.

Referenced by chooseModulesToKill(), killStub(), killStubInDeadModule(), and StubKiller().

edm::Service<edm::RandomNumberGenerator> tmtt::StubKiller::rndmService_
private

Definition at line 74 of file StubKiller.h.

Referenced by StubKiller().

const TrackerGeometry* tmtt::StubKiller::trackerGeometry_
private

Definition at line 59 of file StubKiller.h.

Referenced by addDeadLayerModulesToDeadModuleList(), chooseModulesToKill(), and killStub().

const TrackerTopology* tmtt::StubKiller::trackerTopology_
private

Definition at line 58 of file StubKiller.h.

Referenced by addDeadLayerModulesToDeadModuleList(), chooseModulesToKill(), and killStub().