Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2023-03-17 11:22:56

0001 #include "SimpleNavigationSchool.h"
0002 #include "RecoTracker/TkDetLayers/interface/GeometricSearchTracker.h"
0003 #include "FWCore/Utilities/interface/Visibility.h"
0004 
0005 #include <vector>
0006 
0007 /** Concrete navigation school for the Tracker, connecting disks only for traversing tracks : moslty beam halo muon 
0008  */
0009 class dso_hidden BeamHaloNavigationSchool final : public SimpleNavigationSchool {
0010 public:
0011   BeamHaloNavigationSchool(const GeometricSearchTracker* theTracker, const MagneticField* field);
0012   ~BeamHaloNavigationSchool() override { cleanMemory(); }
0013 
0014 protected:
0015   //addon to SimpleNavigationSchool
0016   void linkOtherEndLayers(SymmetricLayerFinder& symFinder);
0017   void addInward(const DetLayer* det, const FDLC& news);
0018   void addInward(const DetLayer* det, const ForwardDetLayer* newF);
0019   void establishInverseRelations() override;
0020   FDLC reachableFromHorizontal();
0021 };
0022 
0023 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0024 
0025 #include "TrackingTools/DetLayers/interface/DetLessZ.h"
0026 #include "TrackingTools/DetLayers/interface/BarrelDetLayer.h"
0027 #include "TrackingTools/DetLayers/interface/ForwardDetLayer.h"
0028 
0029 #include "SymmetricLayerFinder.h"
0030 #include "SimpleBarrelNavigableLayer.h"
0031 #include "SimpleForwardNavigableLayer.h"
0032 #include "SimpleNavigableLayer.h"
0033 
0034 using namespace std;
0035 
0036 BeamHaloNavigationSchool::BeamHaloNavigationSchool(const GeometricSearchTracker* theInputTracker,
0037                                                    const MagneticField* field) {
0038   edm::LogInfo("BeamHaloNavigationSchool") << "*********Running BeamHaloNavigationSchool *********";
0039   theBarrelLength = 0;
0040   theField = field;
0041   theTracker = theInputTracker;
0042   theAllDetLayersInSystem = &theInputTracker->allLayers();
0043   theAllNavigableLayer.resize(theInputTracker->allLayers().size(), nullptr);
0044 
0045   // Get barrel layers
0046   /*sideways does not need barrels*/
0047   /*  vector<BarrelDetLayer*> blc = theTracker->barrelLayers(); 
0048       for ( vector<BarrelDetLayer*>::iterator i = blc.begin(); i != blc.end(); i++) {
0049       theBarrelLayers.push_back( (*i) );
0050       }*/
0051 
0052   // get forward layers
0053   for (auto const& l : theTracker->forwardLayers()) {
0054     theForwardLayers.push_back(l);
0055   }
0056 
0057   FDLI middle =
0058       find_if(theForwardLayers.begin(), theForwardLayers.end(), [](auto const* a) { return a->position().z() >= 0.0; });
0059   theLeftLayers = FDLC(theForwardLayers.begin(), middle);
0060   theRightLayers = FDLC(middle, theForwardLayers.end());
0061 
0062   SymmetricLayerFinder symFinder(theForwardLayers);
0063 
0064   // only work on positive Z side; negative by mirror symmetry later
0065   /*sideways does not need barrels*/
0066   //  linkBarrelLayers( symFinder);
0067 
0068   linkForwardLayers(symFinder);
0069 
0070   setState(navigableLayers());
0071 
0072   LogDebug("BeamHaloNavigationSchool") << "inverse relation";
0073   establishInverseRelations();
0074 
0075   //add the necessary inward links to end caps
0076   LogDebug("BeamHaloNavigationSchool") << "linkOtherEndLayer";
0077 
0078   linkOtherEndLayers(symFinder);
0079 
0080   //set checkCrossing = false to all layers
0081   SimpleNavigationSchool::StateType allLayers = navigableLayers();
0082   SimpleNavigationSchool::StateType::iterator layerIt = allLayers.begin();
0083   SimpleNavigationSchool::StateType::iterator layerIt_end = allLayers.end();
0084   for (; layerIt != layerIt_end; ++layerIt) {
0085     //convert to SimpleNavigableLayer
0086     SimpleNavigableLayer* snl = dynamic_cast<SimpleNavigableLayer*>(*layerIt);
0087     if (!snl) {
0088       edm::LogError("BeamHaloNavigationSchool") << "navigable layer not casting to simplenavigablelayer.";
0089       continue;
0090     }
0091     snl->setCheckCrossingSide(false);
0092   }
0093 }
0094 
0095 void BeamHaloNavigationSchool::establishInverseRelations() {
0096   // find for each layer which are the barrel and forward
0097   // layers that point to it
0098   typedef map<const DetLayer*, vector<BarrelDetLayer const*>, less<const DetLayer*> > BarrelMapType;
0099   typedef map<const DetLayer*, vector<ForwardDetLayer const*>, less<const DetLayer*> > ForwardMapType;
0100 
0101   BarrelMapType reachedBarrelLayersMap;
0102   ForwardMapType reachedForwardLayersMap;
0103 
0104   for (BDLI bli = theBarrelLayers.begin(); bli != theBarrelLayers.end(); bli++) {
0105     DLC reachedLC = nextLayers(**bli, insideOut);
0106     for (DLI i = reachedLC.begin(); i != reachedLC.end(); i++) {
0107       reachedBarrelLayersMap[*i].push_back(*bli);
0108     }
0109   }
0110 
0111   for (FDLI fli = theForwardLayers.begin(); fli != theForwardLayers.end(); fli++) {
0112     DLC reachedLC = nextLayers(**fli, insideOut);
0113     for (DLI i = reachedLC.begin(); i != reachedLC.end(); i++) {
0114       reachedForwardLayersMap[*i].push_back(*fli);
0115     }
0116   }
0117 
0118   for (auto const i : theTracker->allLayers()) {
0119     SimpleNavigableLayer* navigableLayer = dynamic_cast<SimpleNavigableLayer*>(theAllNavigableLayer[i->seqNum()]);
0120     if (!navigableLayer) {
0121       edm::LogInfo("BeamHaloNavigationSchool")
0122           << "a detlayer does not have a navigable layer, which is normal in beam halo navigation.";
0123     }
0124     if (navigableLayer) {
0125       navigableLayer->setInwardLinks(reachedBarrelLayersMap[i], reachedForwardLayersMap[i], TkLayerLess(outsideIn, i));
0126     }
0127   }
0128 }
0129 
0130 void BeamHaloNavigationSchool::linkOtherEndLayers(SymmetricLayerFinder& symFinder) {
0131   LogDebug("BeamHaloNavigationSchool") << "reachable from horizontal";
0132   //generally, on the right side, what are the forward layers reachable from the horizontal
0133   FDLC reachableFL = reachableFromHorizontal();
0134 
0135   //even simpler navigation from end to end.
0136   //for each of them
0137   for (FDLI fl = reachableFL.begin(); fl != reachableFL.end(); fl++) {
0138     LogDebug("BeamHaloNavigationSchool") << "adding inward from right";
0139     //link it inward to the mirror reachable from horizontal
0140     addInward(static_cast<DetLayer const*>(*fl), symFinder.mirror(*fl));
0141 
0142     LogDebug("BeamHaloNavigationSchool") << "adding inward from mirror of right (left?)";
0143     addInward(static_cast<DetLayer const*>(symFinder.mirror(*fl)), *fl);
0144   }
0145 }
0146 
0147 void BeamHaloNavigationSchool::addInward(const DetLayer* det, const ForwardDetLayer* newF) {
0148   //get the navigable layer for this DetLayer
0149   SimpleNavigableLayer* navigableLayer = dynamic_cast<SimpleNavigableLayer*>(theAllNavigableLayer[(det)->seqNum()]);
0150 
0151   LogDebug("BeamHaloNavigationSchool") << "retreive the nextlayer outsidein";
0152   //get the inward reachable layers.
0153   DLC inwardsLayers(navigableLayer->nextLayers(outsideIn));
0154 
0155   LogDebug("BeamHaloNavigationSchool") << "split them barrel/forward";
0156   // split barrel and forward layers
0157   BDLC inwardsBarrel;
0158   FDLC inwardsForward;
0159   for (DLC::iterator dli = inwardsLayers.begin(); dli != inwardsLayers.end(); dli++) {
0160     if ((**dli).location() == GeomDetEnumerators::barrel)
0161       inwardsBarrel.push_back(static_cast<const BarrelDetLayer*>(*dli));
0162     else
0163       inwardsForward.push_back(static_cast<const ForwardDetLayer*>(*dli));
0164   }
0165   LogDebug("BeamHaloNavigationSchool") << "add the new ones";
0166   //add the other forward layers provided
0167   inwardsForward.push_back(newF);
0168 
0169   LogDebug("BeamHaloNavigationSchool") << "no duplicate please";
0170   sort(inwardsForward.begin(), inwardsForward.end());  //if you don't sort, unique will not work
0171   //  FDLI read = inwardsForward.begin();
0172   //  std::stringstream showMe;
0173   //  for (;read !=inwardsForward.end();++read)  showMe<<" layer p: "<<*read<<"\n";
0174   //  LogDebug("BeamHaloNavigationSchool")<<"list of layer pointers: \n"<<showMe.str();
0175 
0176   FDLI new_end = unique(inwardsForward.begin(), inwardsForward.end());
0177   //  if (new_end!=inwardsForward.end()) LogDebug("BeamHaloNavigationSchool")<<"removing duplicates here";
0178   inwardsForward.erase(new_end, inwardsForward.end());
0179 
0180   LogDebug("BeamHaloNavigationSchool") << "set back the inward links (no duplicate)";
0181   //  set them back to the navigable layer
0182   navigableLayer->setInwardLinks(inwardsBarrel, inwardsForward, TkLayerLess(outsideIn, det));
0183 }
0184 
0185 void BeamHaloNavigationSchool::addInward(const DetLayer* det, const FDLC& news) {
0186   //get the navigable layer for this DetLayer
0187   SimpleNavigableLayer* navigableLayer = dynamic_cast<SimpleNavigableLayer*>(theAllNavigableLayer[(det)->seqNum()]);
0188 
0189   LogDebug("BeamHaloNavigationSchool") << "retreive the nextlayer outsidein";
0190   //get the inward reachable layers.
0191   DLC inwardsLayers(navigableLayer->nextLayers(outsideIn));
0192 
0193   LogDebug("BeamHaloNavigationSchool") << "split them barrel/forward";
0194   // split barrel and forward layers
0195   BDLC inwardsBarrel;
0196   FDLC inwardsForward;
0197   for (DLC::iterator dli = inwardsLayers.begin(); dli != inwardsLayers.end(); dli++) {
0198     if ((**dli).location() == GeomDetEnumerators::barrel)
0199       inwardsBarrel.push_back(static_cast<const BarrelDetLayer*>(*dli));
0200     else
0201       inwardsForward.push_back(static_cast<const ForwardDetLayer*>(*dli));
0202   }
0203 
0204   LogDebug("BeamHaloNavigationSchool") << "add the new ones";
0205   //add the other forward layers provided
0206   inwardsForward.insert(inwardsForward.end(), news.begin(), news.end());
0207 
0208   LogDebug("BeamHaloNavigationSchool") << "no duplicate please";
0209   FDLI new_end = unique(inwardsForward.begin(), inwardsForward.end());
0210   inwardsForward.erase(new_end, inwardsForward.end());
0211 
0212   LogDebug("BeamHaloNavigationSchool") << "set back the inward links (no duplicate)";
0213   //  set them back to the navigable layer
0214   navigableLayer->setInwardLinks(inwardsBarrel, inwardsForward, TkLayerLess(outsideIn, det));
0215 }
0216 
0217 BeamHaloNavigationSchool::FDLC BeamHaloNavigationSchool::reachableFromHorizontal() {
0218   //determine which is the list of forward layers that can be reached from inside-out
0219   //at horizontal direction
0220 
0221   FDLC myRightLayers(theRightLayers);
0222   FDLI begin = myRightLayers.begin();
0223   FDLI end = myRightLayers.end();
0224 
0225   //sort along Z to be sure
0226   sort(begin, end, isDetLessZ);
0227 
0228   FDLC reachableFL;
0229 
0230   begin = myRightLayers.begin();
0231   end = myRightLayers.end();
0232 
0233   //the first one is always reachable
0234   reachableFL.push_back(*begin);
0235   FDLI current = begin;
0236   for (FDLI i = begin + 1; i != end; i++) {
0237     //is the previous layer NOT masking this one
0238     //inner radius smaller OR outer radius bigger
0239     if ((**i).specificSurface().innerRadius() < (**current).specificSurface().innerRadius() ||
0240         (**i).specificSurface().outerRadius() > (**current).specificSurface().outerRadius()) {  //not masked
0241       reachableFL.push_back(*i);
0242       current = i;
0243     }
0244   }
0245   return reachableFL;
0246 }
0247 
0248 #include "FWCore/PluginManager/interface/ModuleDef.h"
0249 #include "FWCore/Framework/interface/MakerMacros.h"
0250 
0251 #include "NavigationSchoolFactory.h"
0252 #include "TrackingTools/DetLayers/interface/NavigationSchool.h"
0253 DEFINE_EDM_PLUGIN(NavigationSchoolFactory, BeamHaloNavigationSchool, "BeamHaloNavigationSchool");