Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 #ifndef TrackingRegion_H
0002 #define TrackingRegion_H
0003 
0004 /** \class TrackingRegion
0005  * kinematic data common to 
0006  * some concreate implementations of TrackingRegion.
0007  */
0008 
0009 #include "DataFormats/GeometryVector/interface/GlobalVector.h"
0010 #include "DataFormats/GeometryVector/interface/GlobalPoint.h"
0011 
0012 #include "RecoTracker/TkMSParametrization/interface/PixelRecoRange.h"
0013 #include "TrackingTools/DetLayers/interface/DetLayer.h"
0014 #include "RecoTracker/TkTrackingRegions/interface/HitRZCompatibility.h"
0015 #include "RecoTracker/TkTrackingRegions/interface/HitEtaCheck.h"
0016 #include "RecoTracker/TkTrackingRegions/interface/HitRCheck.h"
0017 #include "RecoTracker/TkTrackingRegions/interface/HitZCheck.h"
0018 
0019 #include "TrackingTools/TransientTrackingRecHit/interface/SeedingLayerSetsHits.h"
0020 
0021 #include "DataFormats/TrackReco/interface/Track.h"
0022 #include "DataFormats/TrackReco/interface/TrackFwd.h"
0023 
0024 #include <utility>
0025 
0026 #include <sstream>
0027 #include <memory>
0028 #include <vector>
0029 #include <string>
0030 
0031 class DetLayer;
0032 class HitRZCompatibility;
0033 class BarrelDetLayer;
0034 class ForwardDetLayer;
0035 
0036 namespace edm {
0037   class Event;
0038 }
0039 
0040 class TrackingRegion {
0041 public:
0042 public:
0043   virtual ~TrackingRegion() {}
0044   typedef PixelRecoRange<float> Range;
0045   typedef SeedingLayerSetsHits::ConstRecHitPointer Hit;
0046   typedef SeedingLayerSetsHits::Hits Hits;
0047 
0048 public:
0049   TrackingRegion() {}
0050 
0051   TrackingRegion(const GlobalVector& direction,
0052                  const GlobalPoint& originPos,
0053                  const Range& invPtRange,
0054                  const float& originRBound,
0055                  const float& originZBound)
0056       : theDirection(direction),
0057         theUnitDirection(direction.unit()),
0058         theVertexPos(originPos),
0059         theInvPtRange(invPtRange),
0060         thePhi(direction.barePhi()),
0061         thePtMin(1.f / std::max(std::abs(invPtRange.max()), std::abs(invPtRange.min()))),
0062         theVertexRBound(originRBound),
0063         theVertexZBound(originZBound) {}
0064 
0065   /// the direction around which region is constructed
0066   GlobalVector const& direction() const { return theDirection; }
0067   GlobalVector const& unitDirection() const { return theUnitDirection; }
0068 
0069   float phiDirection() const { return thePhi; }
0070 
0071   /** The origin (centre,vertex) of the region. <BR> 
0072   *  The origin with bounds is ment to constraint point of the <BR>
0073   *  closest approach of the track to the beam line
0074   */
0075   GlobalPoint const& origin() const { return theVertexPos; }
0076 
0077   /// bounds the particle vertex in the transverse plane
0078   float originRBound() const { return theVertexRBound; }
0079 
0080   /// bounds the particle vertex in the longitudinal plane
0081   float originZBound() const { return theVertexZBound; }
0082 
0083   /// minimal pt of interest
0084   float ptMin() const { return thePtMin; }
0085 
0086   /// inverse pt range
0087   Range invPtRange() const { return theInvPtRange; }
0088 
0089   /// utility to check eta/theta hit compatibility with region constraints
0090   /// and outer hit constraint
0091   virtual std::unique_ptr<HitRZCompatibility> checkRZ(const DetLayer* layer,
0092                                                       const Hit& outerHit,
0093                                                       const DetLayer* outerlayer = nullptr,
0094                                                       float lr = 0,
0095                                                       float gz = 0,
0096                                                       float dr = 0,
0097                                                       float dz = 0) const = 0;
0098 
0099   /// get hits from layer compatible with region constraints
0100   virtual Hits hits(const SeedingLayerSetsHits::SeedingLayer& layer) const = 0;
0101 
0102   /// Set the elements of the mask corresponding to the tracks that are compatable with the region.
0103   /// Does not reset the elements corresponding to the tracks that are not compatible.
0104   virtual void checkTracks(reco::TrackCollection const& tracks, std::vector<bool>& mask) const = 0;
0105 
0106   /// return a boolean mask over the TrackCollection reflecting the compatibility of each track with the region constraints
0107   std::vector<bool> checkTracks(reco::TrackCollection const& tracks) const {
0108     std::vector<bool> region_mask(tracks.size(), false);
0109     checkTracks(tracks, region_mask);
0110     return region_mask;
0111   }
0112   /// clone region with new vertex position
0113   std::unique_ptr<TrackingRegion> restrictedRegion(const GlobalPoint& originPos,
0114                                                    const float& originRBound,
0115                                                    const float& originZBound) const {
0116     auto restr = clone();
0117     restr->theVertexPos = originPos;
0118     restr->theVertexRBound = originRBound;
0119     restr->theVertexZBound = originZBound;
0120     return restr;
0121   }
0122 
0123   virtual std::unique_ptr<TrackingRegion> clone() const = 0;
0124 
0125   virtual std::string name() const { return "TrackingRegion"; }
0126   virtual std::string print() const {
0127     std::ostringstream str;
0128     str << name() << " dir:" << theDirection << " vtx:" << theVertexPos << " dr:" << theVertexRBound
0129         << " dz:" << theVertexZBound << " pt:" << 1. / theInvPtRange.max();
0130     return str.str();
0131   }
0132 
0133   // void setDirection(const GlobalVector & dir ) { theDirection = dir; }
0134 
0135 private:
0136   GlobalVector theDirection;      // this is not direction is momentum...
0137   GlobalVector theUnitDirection;  // the real direction
0138   GlobalPoint theVertexPos;
0139   Range theInvPtRange;
0140   float thePhi;
0141   float thePtMin;
0142   float theVertexRBound;
0143   float theVertexZBound;
0144 };
0145 
0146 using TrackingRegionBase = TrackingRegion;
0147 
0148 #endif