Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:22:01

0001 #include "L1Trigger/TrackFindingTracklet/interface/StubKiller.h"
0002 
0003 using namespace std;
0004 
0005 StubKiller::StubKiller()
0006     : killScenario_(0),
0007       trackerTopology_(nullptr),
0008       trackerGeometry_(nullptr),
0009       layersToKill_(vector<int>()),
0010       minPhiToKill_(0),
0011       maxPhiToKill_(0),
0012       minZToKill_(0),
0013       maxZToKill_(0),
0014       minRToKill_(0),
0015       maxRToKill_(0),
0016       fractionOfStubsToKillInLayers_(0),
0017       fractionOfStubsToKillEverywhere_(0),
0018       fractionOfModulesToKillEverywhere_(0) {}
0019 
0020 void StubKiller::initialise(unsigned int killScenario,
0021                             const TrackerTopology* trackerTopology,
0022                             const TrackerGeometry* trackerGeometry) {
0023   killScenario_ = killScenario;
0024   trackerTopology_ = trackerTopology;
0025   trackerGeometry_ = trackerGeometry;
0026 
0027   switch (killScenario_) {
0028     // kill layer 5 in one quadrant + 5% random module loss to connect to what was done before
0029     case 1:
0030       layersToKill_ = {5};
0031       minPhiToKill_ = 0.0;
0032       maxPhiToKill_ = TMath::PiOver2();
0033       minZToKill_ = -1000.0;
0034       maxZToKill_ = 0.0;
0035       minRToKill_ = 0.0;
0036       maxRToKill_ = 1000.0;
0037       fractionOfStubsToKillInLayers_ = 1;
0038       fractionOfStubsToKillEverywhere_ = 0;
0039       fractionOfModulesToKillEverywhere_ = 0.05;
0040       break;
0041 
0042     // kill layer 1 in one quadrant + 5% random module loss
0043     case 2:
0044       layersToKill_ = {1};
0045       minPhiToKill_ = 0.0;
0046       maxPhiToKill_ = TMath::PiOver2();
0047       minZToKill_ = -1000.0;
0048       maxZToKill_ = 0.0;
0049       minRToKill_ = 0.0;
0050       maxRToKill_ = 1000.0;
0051       fractionOfStubsToKillInLayers_ = 1;
0052       fractionOfStubsToKillEverywhere_ = 0;
0053       fractionOfModulesToKillEverywhere_ = 0.05;
0054       break;
0055 
0056     // kill layer 1 + layer 2, both in same quadrant
0057     case 3:
0058       layersToKill_ = {1, 2};
0059       minPhiToKill_ = 0.0;
0060       maxPhiToKill_ = TMath::PiOver2();
0061       minZToKill_ = -1000.0;
0062       maxZToKill_ = 0.0;
0063       minRToKill_ = 0.0;
0064       maxRToKill_ = 1000.0;
0065       fractionOfStubsToKillInLayers_ = 1;
0066       fractionOfStubsToKillEverywhere_ = 0;
0067       fractionOfModulesToKillEverywhere_ = 0;
0068       break;
0069 
0070     // kill layer 1 and disk 1, both in same quadrant
0071     case 4:
0072       layersToKill_ = {1, 11};
0073       minPhiToKill_ = 0.0;
0074       maxPhiToKill_ = TMath::PiOver2();
0075       minZToKill_ = -1000.0;
0076       maxZToKill_ = 0.0;
0077       minRToKill_ = 0.0;
0078       maxRToKill_ = 66.5;
0079       fractionOfStubsToKillInLayers_ = 1;
0080       fractionOfStubsToKillEverywhere_ = 0;
0081       fractionOfModulesToKillEverywhere_ = 0;
0082       break;
0083 
0084     // 5% random module loss throughout tracker
0085     case 5:
0086       layersToKill_ = {};
0087       fractionOfStubsToKillInLayers_ = 0;
0088       fractionOfStubsToKillEverywhere_ = 0;
0089       fractionOfModulesToKillEverywhere_ = 0.05;
0090       break;
0091 
0092     // 1% random module loss throughout tracker
0093     case 6:
0094       layersToKill_ = {};
0095       fractionOfStubsToKillInLayers_ = 0;
0096       fractionOfStubsToKillEverywhere_ = 0;
0097       fractionOfModulesToKillEverywhere_ = 0.01;
0098       break;
0099 
0100     // kill layer 5 in one quadrant + 1% random module loss
0101     case 7:
0102       layersToKill_ = {5};
0103       minPhiToKill_ = 0.0;
0104       maxPhiToKill_ = TMath::PiOver2();
0105       minZToKill_ = -1000.0;
0106       maxZToKill_ = 0.0;
0107       minRToKill_ = 0.0;
0108       maxRToKill_ = 1000.0;
0109       fractionOfStubsToKillInLayers_ = 1;
0110       fractionOfStubsToKillEverywhere_ = 0;
0111       fractionOfModulesToKillEverywhere_ = 0.01;
0112       break;
0113 
0114     // kill layer 1 in one quadrant +1 % random module loss
0115     case 8:
0116       layersToKill_ = {1};
0117       minPhiToKill_ = 0.0;
0118       maxPhiToKill_ = TMath::PiOver2();
0119       minZToKill_ = -1000.0;
0120       maxZToKill_ = 0.0;
0121       minRToKill_ = 0.0;
0122       maxRToKill_ = 1000.0;
0123       fractionOfStubsToKillInLayers_ = 1;
0124       fractionOfStubsToKillEverywhere_ = 0;
0125       fractionOfModulesToKillEverywhere_ = 0.01;
0126       break;
0127 
0128     // 10% random module loss throughout tracker
0129     case 9:
0130       layersToKill_ = {};
0131       fractionOfStubsToKillInLayers_ = 0;
0132       fractionOfStubsToKillEverywhere_ = 0;
0133       fractionOfModulesToKillEverywhere_ = 0.10;
0134       break;
0135   }
0136 
0137   deadModules_.clear();
0138   if (fractionOfModulesToKillEverywhere_ > 0) {
0139     this->chooseModulesToKill();
0140   }
0141   this->addDeadLayerModulesToDeadModuleList();
0142 }
0143 
0144 void StubKiller::chooseModulesToKill() {
0145   TRandom3* randomGenerator = new TRandom3();
0146   randomGenerator->SetSeed(0);
0147 
0148   for (const GeomDetUnit* gd : trackerGeometry_->detUnits()) {
0149     if (!trackerTopology_->isLower(gd->geographicalId()))
0150       continue;
0151     if (randomGenerator->Uniform(0.0, 1.0) < fractionOfModulesToKillEverywhere_) {
0152       deadModules_[gd->geographicalId()] = 1;
0153     }
0154   }
0155 }
0156 
0157 void StubKiller::addDeadLayerModulesToDeadModuleList() {
0158   for (const GeomDetUnit* gd : trackerGeometry_->detUnits()) {
0159     float moduleR = gd->position().perp();
0160     float moduleZ = gd->position().z();
0161     float modulePhi = gd->position().phi();
0162     DetId geoDetId = gd->geographicalId();
0163     bool isInBarrel = geoDetId.subdetId() == StripSubdetector::TOB || geoDetId.subdetId() == StripSubdetector::TIB;
0164 
0165     int layerID = 0;
0166     if (isInBarrel) {
0167       layerID = trackerTopology_->layer(geoDetId);
0168     } else {
0169       layerID = 10 * trackerTopology_->side(geoDetId) + trackerTopology_->tidWheel(geoDetId);
0170     }
0171 
0172     if (find(layersToKill_.begin(), layersToKill_.end(), layerID) != layersToKill_.end()) {
0173       if (modulePhi < -1.0 * TMath::Pi())
0174         modulePhi += 2.0 * TMath::Pi();
0175       else if (modulePhi > TMath::Pi())
0176         modulePhi -= 2.0 * TMath::Pi();
0177 
0178       if (modulePhi > minPhiToKill_ && modulePhi < maxPhiToKill_ && moduleZ > minZToKill_ && moduleZ < maxZToKill_ &&
0179           moduleR > minRToKill_ && moduleR < maxRToKill_) {
0180         if (deadModules_.find(gd->geographicalId()) == deadModules_.end()) {
0181           deadModules_[gd->geographicalId()] = fractionOfStubsToKillInLayers_;
0182         }
0183       }
0184     }
0185   }
0186 }
0187 
0188 bool StubKiller::killStub(const TTStub<Ref_Phase2TrackerDigi_>* stub) {
0189   if (killScenario_ == 0)
0190     return false;
0191   else {
0192     bool killStubRandomly = killStub(stub,
0193                                      layersToKill_,
0194                                      minPhiToKill_,
0195                                      maxPhiToKill_,
0196                                      minZToKill_,
0197                                      maxZToKill_,
0198                                      minRToKill_,
0199                                      maxRToKill_,
0200                                      fractionOfStubsToKillInLayers_,
0201                                      fractionOfStubsToKillEverywhere_);
0202     bool killStubInDeadModules = killStubInDeadModule(stub);
0203     return killStubRandomly || killStubInDeadModules;
0204   }
0205 }
0206 
0207 // layersToKill - a vector stating the layers we are killing stubs in.  Can be an empty vector.
0208 // Barrel layers are encoded as 1-6. The endcap layers are encoded as 11-15 (-z) and 21-25 (+z)
0209 // min/max Phi/Z/R - stubs within the region specified by these boundaries and layersToKill are flagged for killing
0210 // fractionOfStubsToKillInLayers - The fraction of stubs to kill in the specified layers/region.
0211 // fractionOfStubsToKillEverywhere - The fraction of stubs to kill throughout the tracker
0212 bool StubKiller::killStub(const TTStub<Ref_Phase2TrackerDigi_>* stub,
0213                           const vector<int> layersToKill,
0214                           const double minPhiToKill,
0215                           const double maxPhiToKill,
0216                           const double minZToKill,
0217                           const double maxZToKill,
0218                           const double minRToKill,
0219                           const double maxRToKill,
0220                           const double fractionOfStubsToKillInLayers,
0221                           const double fractionOfStubsToKillEverywhere) {
0222   // Only kill stubs in specified layers
0223   if (!layersToKill.empty()) {
0224     // Get the layer the stub is in, and check if it's in the layer you want to kill
0225     DetId stackDetid = stub->getDetId();
0226     DetId geoDetId(stackDetid.rawId() + 1);
0227 
0228     bool isInBarrel = geoDetId.subdetId() == StripSubdetector::TOB || geoDetId.subdetId() == StripSubdetector::TIB;
0229 
0230     int layerID = 0;
0231     if (isInBarrel) {
0232       layerID = trackerTopology_->layer(geoDetId);
0233     } else {
0234       layerID = 10 * trackerTopology_->side(geoDetId) + trackerTopology_->tidWheel(geoDetId);
0235     }
0236 
0237     if (find(layersToKill.begin(), layersToKill.end(), layerID) != layersToKill.end()) {
0238       // Get the phi and z of stub, and check if it's in the region you want to kill
0239       const GeomDetUnit* det0 = trackerGeometry_->idToDetUnit(geoDetId);
0240       const PixelGeomDetUnit* theGeomDet = dynamic_cast<const PixelGeomDetUnit*>(det0);
0241       const PixelTopology* topol = dynamic_cast<const PixelTopology*>(&(theGeomDet->specificTopology()));
0242       MeasurementPoint measurementPoint = stub->clusterRef(0)->findAverageLocalCoordinatesCentered();
0243       LocalPoint clustlp = topol->localPosition(measurementPoint);
0244       GlobalPoint pos = theGeomDet->surface().toGlobal(clustlp);
0245 
0246       // Just in case phi is outside of -pi -> pi
0247       double stubPhi = pos.phi();
0248       if (stubPhi < -1.0 * TMath::Pi())
0249         stubPhi += 2.0 * TMath::Pi();
0250       else if (stubPhi > TMath::Pi())
0251         stubPhi -= 2.0 * TMath::Pi();
0252 
0253       if (stubPhi > minPhiToKill && stubPhi < maxPhiToKill && pos.z() > minZToKill && pos.z() < maxZToKill &&
0254           pos.perp() > minRToKill && pos.perp() < maxRToKill) {
0255         // Kill fraction of stubs
0256         if (fractionOfStubsToKillInLayers == 1) {
0257           return true;
0258         } else {
0259           static TRandom randomGenerator;
0260           if (randomGenerator.Rndm() < fractionOfStubsToKillInLayers) {
0261             return true;
0262           }
0263         }
0264       }
0265     }
0266   }
0267 
0268   // Kill fraction of stubs throughout tracker
0269   if (fractionOfStubsToKillEverywhere > 0) {
0270     static TRandom randomGenerator;
0271     if (randomGenerator.Rndm() < fractionOfStubsToKillEverywhere) {
0272       return true;
0273     }
0274   }
0275 
0276   return false;
0277 }
0278 
0279 bool StubKiller::killStubInDeadModule(const TTStub<Ref_Phase2TrackerDigi_>* stub) {
0280   if (!deadModules_.empty()) {
0281     DetId stackDetid = stub->getDetId();
0282     DetId geoDetId(stackDetid.rawId() + 1);
0283     if (deadModules_.find(geoDetId) != deadModules_.end())
0284       return true;
0285   }
0286 
0287   return false;
0288 }