Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:28:51

0001 #include "SimpleNavigationSchool.h"
0002 
0003 #include "RecoTracker/TkDetLayers/interface/GeometricSearchTracker.h"
0004 
0005 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0006 #include "FWCore/Utilities/interface/Visibility.h"
0007 
0008 #include <vector>
0009 
0010 //class FakeDetLayer;
0011 
0012 /** Concrete navigation school for cosmics in the Tracker
0013  */
0014 
0015 class dso_hidden CosmicNavigationSchool : public SimpleNavigationSchool {
0016 public:
0017   CosmicNavigationSchool(const GeometricSearchTracker* theTracker, const MagneticField* field);
0018   ~CosmicNavigationSchool() override { cleanMemory(); }
0019 
0020   class CosmicNavigationSchoolConfiguration {
0021   public:
0022     CosmicNavigationSchoolConfiguration()
0023         : noPXB(false),
0024           noPXF(false),
0025           noTOB(false),
0026           noTIB(false),
0027           noTEC(false),
0028           noTID(false),
0029           self(false),
0030           allSelf(false) {}
0031     CosmicNavigationSchoolConfiguration(const edm::ParameterSet& conf);
0032     bool noPXB;
0033     bool noPXF;
0034     bool noTOB;
0035     bool noTIB;
0036     bool noTEC;
0037     bool noTID;
0038 
0039     bool self;
0040     bool allSelf;
0041   };
0042 
0043   void build(const GeometricSearchTracker* theTracker,
0044              const MagneticField* field,
0045              const CosmicNavigationSchoolConfiguration conf);
0046 
0047 protected:
0048   CosmicNavigationSchool() {}
0049 
0050 private:
0051   //FakeDetLayer* theFakeDetLayer;
0052   void linkBarrelLayers(SymmetricLayerFinder& symFinder) override;
0053   //void linkForwardLayers( SymmetricLayerFinder& symFinder);
0054   using SimpleNavigationSchool::establishInverseRelations;
0055   void establishInverseRelations(SymmetricLayerFinder& symFinder);
0056   void buildAdditionalBarrelLinks();
0057   void buildAdditionalForwardLinks(SymmetricLayerFinder& symFinder);
0058 };
0059 
0060 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0061 
0062 #include "SimpleBarrelNavigableLayer.h"
0063 #include "SimpleForwardNavigableLayer.h"
0064 #include "SimpleNavigableLayer.h"
0065 #include "SymmetricLayerFinder.h"
0066 
0067 #include "TrackingTools/DetLayers/interface/BarrelDetLayer.h"
0068 #include "TrackingTools/DetLayers/interface/ForwardDetLayer.h"
0069 #include "TrackingTools/DetLayers/interface/DetLessZ.h"
0070 
0071 #include <functional>
0072 #include <algorithm>
0073 #include <map>
0074 #include <cmath>
0075 
0076 using namespace std;
0077 
0078 CosmicNavigationSchool::CosmicNavigationSchoolConfiguration::CosmicNavigationSchoolConfiguration(
0079     const edm::ParameterSet& conf) {
0080   noPXB = conf.getParameter<bool>("noPXB");
0081   noPXF = conf.getParameter<bool>("noPXF");
0082   noTIB = conf.getParameter<bool>("noTIB");
0083   noTID = conf.getParameter<bool>("noTID");
0084   noTOB = conf.getParameter<bool>("noTOB");
0085   noTEC = conf.getParameter<bool>("noTEC");
0086   self = conf.getParameter<bool>("selfSearch");
0087   allSelf = conf.getParameter<bool>("allSelf");
0088 }
0089 
0090 CosmicNavigationSchool::CosmicNavigationSchool(const GeometricSearchTracker* theInputTracker,
0091                                                const MagneticField* field) {
0092   build(theInputTracker, field, CosmicNavigationSchoolConfiguration());
0093 }
0094 
0095 void CosmicNavigationSchool::build(const GeometricSearchTracker* theInputTracker,
0096                                    const MagneticField* field,
0097                                    const CosmicNavigationSchoolConfiguration conf) {
0098   LogTrace("CosmicNavigationSchool") << "*********Running CosmicNavigationSchool***********";
0099   theBarrelLength = 0;
0100   theField = field;
0101   theTracker = theInputTracker;
0102 
0103   theAllDetLayersInSystem = &theInputTracker->allLayers();
0104   theAllNavigableLayer.resize(theInputTracker->allLayers().size(), nullptr);
0105 
0106   // Get barrel layers
0107   vector<BarrelDetLayer const*> const& blc = theTracker->barrelLayers();
0108   for (auto i = blc.begin(); i != blc.end(); i++) {
0109     if (conf.noPXB && GeomDetEnumerators::isTrackerPixel((*i)->subDetector()))
0110       continue;
0111     if (conf.noTOB && (*i)->subDetector() == GeomDetEnumerators::TOB)
0112       continue;
0113     if (conf.noTIB && (*i)->subDetector() == GeomDetEnumerators::TIB)
0114       continue;
0115     theBarrelLayers.push_back((*i));
0116   }
0117 
0118   // get forward layers
0119   vector<ForwardDetLayer const*> const& flc = theTracker->forwardLayers();
0120   for (auto i = flc.begin(); i != flc.end(); i++) {
0121     if (conf.noPXF && GeomDetEnumerators::isTrackerPixel((*i)->subDetector()))
0122       continue;
0123     if (conf.noTEC && (*i)->subDetector() == GeomDetEnumerators::TEC)
0124       continue;
0125     if (conf.noTID && (*i)->subDetector() == GeomDetEnumerators::TID)
0126       continue;
0127     theForwardLayers.push_back((*i));
0128   }
0129 
0130   FDLI middle =
0131       find_if(theForwardLayers.begin(), theForwardLayers.end(), [](auto const* a) { return a->position().z() >= 0.0; });
0132   theLeftLayers = FDLC(theForwardLayers.begin(), middle);
0133   theRightLayers = FDLC(middle, theForwardLayers.end());
0134 
0135   SymmetricLayerFinder symFinder(theForwardLayers);
0136 
0137   // only work on positive Z side; negative by mirror symmetry later
0138   linkBarrelLayers(symFinder);
0139   linkForwardLayers(symFinder);
0140   establishInverseRelations(symFinder);
0141 
0142   if (conf.self) {
0143     // set the self search by hand
0144     //   NavigationSetter setter(*this);
0145 
0146     //add TOB1->TOB1 inward link
0147     const std::vector<const BarrelDetLayer*>& tobL = theInputTracker->tobLayers();
0148     if (!tobL.empty()) {
0149       if (conf.allSelf) {
0150         LogDebug("CosmicNavigationSchool") << " adding all TOB self search.";
0151         for (auto lIt = tobL.begin(); lIt != tobL.end(); ++lIt)
0152           dynamic_cast<SimpleNavigableLayer*>(theAllNavigableLayer[(*lIt)->seqNum()])->theSelfSearch = true;
0153       } else {
0154         SimpleNavigableLayer* navigableLayer =
0155             dynamic_cast<SimpleNavigableLayer*>(theAllNavigableLayer[tobL.front()->seqNum()]);
0156         LogDebug("CosmicNavigationSchool") << " adding TOB1 to TOB1.";
0157         navigableLayer->theSelfSearch = true;
0158       }
0159     }
0160     const std::vector<const BarrelDetLayer*>& tibL = theInputTracker->tibLayers();
0161     if (!tibL.empty()) {
0162       if (conf.allSelf) {
0163         LogDebug("CosmicNavigationSchool") << " adding all TIB self search.";
0164         for (auto lIt = tibL.begin(); lIt != tibL.end(); ++lIt)
0165           dynamic_cast<SimpleNavigableLayer*>(theAllNavigableLayer[(*lIt)->seqNum()])->theSelfSearch = true;
0166       } else {
0167         SimpleNavigableLayer* navigableLayer =
0168             dynamic_cast<SimpleNavigableLayer*>(theAllNavigableLayer[tibL.front()->seqNum()]);
0169         LogDebug("CosmicNavigationSchool") << " adding tib1 to tib1.";
0170         navigableLayer->theSelfSearch = true;
0171       }
0172     }
0173     const std::vector<const BarrelDetLayer*>& pxbL = theInputTracker->pixelBarrelLayers();
0174     if (!pxbL.empty()) {
0175       if (conf.allSelf) {
0176         LogDebug("CosmicNavigationSchool") << " adding all PXB self search.";
0177         for (auto lIt = pxbL.begin(); lIt != pxbL.end(); ++lIt)
0178           dynamic_cast<SimpleNavigableLayer*>(theAllNavigableLayer[(*lIt)->seqNum()])->theSelfSearch = true;
0179       } else {
0180         SimpleNavigableLayer* navigableLayer =
0181             dynamic_cast<SimpleNavigableLayer*>(theAllNavigableLayer[pxbL.front()->seqNum()]);
0182         LogDebug("CosmicNavigationSchool") << " adding pxb1 to pxb1.";
0183         navigableLayer->theSelfSearch = true;
0184       }
0185     }
0186   }
0187 }
0188 
0189 void CosmicNavigationSchool::linkBarrelLayers(SymmetricLayerFinder& symFinder) {
0190   //identical to the SimpleNavigationSchool one, but it allows crossing over the tracker
0191   //is some non-standard link is needed, it should probably be added here
0192 
0193   // Link barrel layers outwards
0194   for (BDLI i = theBarrelLayers.begin(); i != theBarrelLayers.end(); i++) {
0195     BDLC reachableBL;
0196     FDLC leftFL;
0197     FDLC rightFL;
0198 
0199     // always add next barrel layer first
0200     if (i + 1 != theBarrelLayers.end())
0201       reachableBL.push_back(*(i + 1));
0202 
0203     // Add closest reachable forward layer (except for last BarrelLayer)
0204     if (i != theBarrelLayers.end() - 1) {
0205       linkNextForwardLayer(*i, rightFL);
0206     }
0207 
0208     // Add next BarrelLayer with length larger than the current BL
0209     if (i + 2 < theBarrelLayers.end()) {
0210       linkNextLargerLayer(i, theBarrelLayers.end(), reachableBL);
0211     }
0212 
0213     theBarrelNLC.push_back(
0214         new SimpleBarrelNavigableLayer(*i, reachableBL, symFinder.mirror(rightFL), rightFL, theField, 5., false));
0215   }
0216 }
0217 
0218 // identical to  SimpleNavigationSchool but for the last additional stuff
0219 void CosmicNavigationSchool::establishInverseRelations(SymmetricLayerFinder& symFinder) {
0220   //again: standard part is identical to SimpleNavigationSchool one.
0221   //After the standard link, special outsideIn links are added
0222 
0223   // NavigationSetter setter(*this);
0224 
0225   setState(navigableLayers());
0226 
0227   // find for each layer which are the barrel and forward
0228   // layers that point to it
0229   typedef map<const DetLayer*, vector<const BarrelDetLayer*>, less<const DetLayer*> > BarrelMapType;
0230   typedef map<const DetLayer*, vector<const ForwardDetLayer*>, less<const DetLayer*> > ForwardMapType;
0231 
0232   BarrelMapType reachedBarrelLayersMap;
0233   ForwardMapType reachedForwardLayersMap;
0234 
0235   for (auto bli : theBarrelLayers) {
0236     auto reachedLC = nextLayers(*bli, insideOut);
0237     for (auto i : reachedLC) {
0238       reachedBarrelLayersMap[i].push_back(bli);
0239     }
0240   }
0241 
0242   for (auto fli : theForwardLayers) {
0243     auto reachedLC = nextLayers(*fli, insideOut);
0244     for (auto i : reachedLC) {
0245       reachedForwardLayersMap[i].push_back(fli);
0246     }
0247   }
0248 
0249   for (auto nl : theAllNavigableLayer) {
0250     if (!nl)
0251       continue;
0252     auto navigableLayer = static_cast<SimpleNavigableLayer*>(nl);
0253     auto dl = nl->detLayer();
0254     navigableLayer->setInwardLinks(reachedBarrelLayersMap[dl], reachedForwardLayersMap[dl]);
0255   }
0256 
0257   //buildAdditionalBarrelLinks();
0258   buildAdditionalForwardLinks(symFinder);
0259 }
0260 
0261 void CosmicNavigationSchool::buildAdditionalBarrelLinks() {
0262   for (auto i = theBarrelLayers.begin(); i != theBarrelLayers.end(); i++) {
0263     SimpleNavigableLayer* navigableLayer = dynamic_cast<SimpleNavigableLayer*>(theAllNavigableLayer[(*i)->seqNum()]);
0264     if (i + 1 != theBarrelLayers.end())
0265       navigableLayer->setAdditionalLink(*(i + 1), outsideIn);
0266   }
0267 }
0268 
0269 void CosmicNavigationSchool::buildAdditionalForwardLinks(SymmetricLayerFinder& symFinder) {
0270   //the first layer of FPIX should not check the crossing side (since there are no inner layers to be tryed first)
0271   SimpleNavigableLayer* firstR =
0272       dynamic_cast<SimpleNavigableLayer*>(theAllNavigableLayer[theRightLayers.front()->seqNum()]);
0273   SimpleNavigableLayer* firstL =
0274       dynamic_cast<SimpleNavigableLayer*>(theAllNavigableLayer[theLeftLayers.front()->seqNum()]);
0275   firstR->setCheckCrossingSide(false);
0276   firstL->setCheckCrossingSide(false);
0277 
0278   for (auto i : theRightLayers) {
0279     //look for first bigger barrel layer and link to it outsideIn
0280     SimpleForwardNavigableLayer* nfl = dynamic_cast<SimpleForwardNavigableLayer*>(theAllNavigableLayer[(i)->seqNum()]);
0281     SimpleForwardNavigableLayer* mnfl =
0282         dynamic_cast<SimpleForwardNavigableLayer*>(theAllNavigableLayer[symFinder.mirror(i)->seqNum()]);
0283     for (auto j : theBarrelLayers) {
0284       if ((i)->specificSurface().outerRadius() < (j)->specificSurface().radius() &&
0285           fabs((i)->specificSurface().position().z()) < (j)->surface().bounds().length() / 2.) {
0286         nfl->setAdditionalLink(j, outsideIn);
0287         mnfl->setAdditionalLink(j, outsideIn);
0288         break;
0289       }
0290     }
0291   }
0292 }
0293 
0294 #include "FWCore/PluginManager/interface/ModuleDef.h"
0295 #include "FWCore/Framework/interface/MakerMacros.h"
0296 
0297 #include "NavigationSchoolFactory.h"
0298 #include "TrackingTools/DetLayers/interface/NavigationSchool.h"
0299 DEFINE_EDM_PLUGIN(NavigationSchoolFactory, CosmicNavigationSchool, "CosmicNavigationSchool");
0300 
0301 #include "RecoTracker/TkDetLayers/interface/GeometricSearchTracker.h"
0302 
0303 #include <vector>
0304 
0305 //class FakeDetLayer;
0306 
0307 /** Concrete navigation school for cosmics in the Tracker
0308  */
0309 
0310 class SkippingLayerCosmicNavigationSchool : public CosmicNavigationSchool {
0311 public:
0312   SkippingLayerCosmicNavigationSchool(const GeometricSearchTracker* theTracker,
0313                                       const MagneticField* field,
0314                                       const CosmicNavigationSchoolConfiguration conf);
0315 
0316   ~SkippingLayerCosmicNavigationSchool() override { cleanMemory(); };
0317 };
0318 
0319 SkippingLayerCosmicNavigationSchool::SkippingLayerCosmicNavigationSchool(
0320     const GeometricSearchTracker* theInputTracker,
0321     const MagneticField* field,
0322     const CosmicNavigationSchoolConfiguration conf) {
0323   build(theInputTracker, field, conf);
0324 }
0325 
0326 #include <FWCore/Utilities/interface/ESInputTag.h>
0327 
0328 #include <memory>
0329 
0330 // user include files
0331 #include "FWCore/Framework/interface/ModuleFactory.h"
0332 #include "FWCore/Framework/interface/ESProducer.h"
0333 
0334 #include "FWCore/Framework/interface/ESHandle.h"
0335 
0336 #include "NavigationSchoolFactory.h"
0337 #include "RecoTracker/Record/interface/NavigationSchoolRecord.h"
0338 
0339 //
0340 // class decleration
0341 //
0342 class dso_hidden SkippingLayerCosmicNavigationSchoolESProducer final : public edm::ESProducer {
0343 public:
0344   SkippingLayerCosmicNavigationSchoolESProducer(const edm::ParameterSet& iConfig) : config_(iConfig) {
0345     //the following line is needed to tell the framework what
0346     // data is being produced
0347     auto cc = setWhatProduced(this, iConfig.getParameter<std::string>("ComponentName"));
0348     magFieldToken_ = cc.consumes();
0349     geometricSearchTrackerToken_ = cc.consumes();
0350   }
0351 
0352   using ReturnType = std::unique_ptr<NavigationSchool>;
0353 
0354   ReturnType produce(const NavigationSchoolRecord&);
0355 
0356   static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0357 
0358 private:
0359   // ----------member data ---------------------------
0360   CosmicNavigationSchool::CosmicNavigationSchoolConfiguration const config_;
0361   edm::ESGetToken<MagneticField, IdealMagneticFieldRecord> magFieldToken_;
0362   edm::ESGetToken<GeometricSearchTracker, TrackerRecoGeometryRecord> geometricSearchTrackerToken_;
0363 };
0364 
0365 SkippingLayerCosmicNavigationSchoolESProducer::ReturnType SkippingLayerCosmicNavigationSchoolESProducer::produce(
0366     const NavigationSchoolRecord& iRecord) {
0367   return std::unique_ptr<NavigationSchool>(new SkippingLayerCosmicNavigationSchool(
0368       &iRecord.get(geometricSearchTrackerToken_), &iRecord.get(magFieldToken_), config_));
0369 }
0370 
0371 void SkippingLayerCosmicNavigationSchoolESProducer::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0372   edm::ParameterSetDescription desc;
0373   desc.add<std::string>("ComponentName");
0374   desc.add<bool>("noPXB");
0375   desc.add<bool>("noPXF");
0376   desc.add<bool>("noTIB");
0377   desc.add<bool>("noTID");
0378   desc.add<bool>("noTOB");
0379   desc.add<bool>("noTEC");
0380   desc.add<bool>("selfSearch");
0381   desc.add<bool>("allSelf");
0382 
0383   descriptions.addDefault(desc);
0384 }
0385 
0386 #include "FWCore/PluginManager/interface/ModuleDef.h"
0387 #include "FWCore/Framework/interface/MakerMacros.h"
0388 
0389 DEFINE_FWK_EVENTSETUP_MODULE(SkippingLayerCosmicNavigationSchoolESProducer);