CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_2_SLHC2/src/SimDataFormats/SLHC/interface/L1TkStub.h

Go to the documentation of this file.
00001 
00002 
00003 
00004 
00005 
00006 
00007 
00008 
00009 
00010 
00011 
00012 
00013 
00014 
00015 #ifndef STACKED_TRACKER_L1TK_STUB_FORMAT_H
00016 #define STACKED_TRACKER_L1TK_STUB_FORMAT_H
00017 
00018 #include "DataFormats/SiPixelDetId/interface/StackedTrackerDetId.h"
00019 #include "CLHEP/Units/PhysicalConstants.h"
00020 
00021 #include "DataFormats/DetId/interface/DetId.h"
00022 #include "DataFormats/Common/interface/Ref.h"
00023 #include "DataFormats/GeometryCommonDetAlgo/interface/MeasurementPoint.h"
00024 #include "DataFormats/GeometryVector/interface/GlobalPoint.h"
00025 #include "DataFormats/GeometryVector/interface/GlobalVector.h"
00026 #include "DataFormats/Common/interface/DetSetVector.h"
00027 #include "SimDataFormats/TrackerDigiSimLink/interface/PixelDigiSimLink.h"
00028 #include "SimDataFormats/TrackingHit/interface/PSimHitContainer.h"
00029 #include "SimDataFormats/SLHC/interface/L1TkCluster.h"
00030 
00037   template< typename T >
00038   class L1TkStub
00039   {
00040     public:
00042       L1TkStub();
00043       L1TkStub( DetId aDetId );
00044 
00046       ~L1TkStub();
00047 
00050 
00052       std::vector< edm::Ptr< L1TkCluster< T > > > getClusterPtrs() const;
00053       const edm::Ptr< L1TkCluster< T > >&         getClusterPtr( unsigned int hitIdentifier ) const;
00054       void                                        addClusterPtr( edm::Ptr< L1TkCluster< T > > aL1TkCluster );
00055 
00059 
00061       DetId getDetId() const;
00062       void  setDetId( DetId aDetId );
00063 
00065       double getTriggerDisplacement() const; 
00066       void   setTriggerDisplacement( int aDisplacement ); 
00067       double getTriggerOffset() const; 
00068       void   setTriggerOffset( int anOffset ); 
00069 
00071       edm::Ptr< SimTrack > getSimTrackPtr() const;
00072       bool                 isGenuine() const;
00073       bool                 isCombinatoric() const;
00074       bool                 isUnknown() const;
00075       int                  findType() const;
00076       unsigned int         findSimTrackId() const;
00077 
00079       void checkSimTrack();
00080 
00082       std::string print( unsigned int i=0 ) const;
00083 
00084     private:
00086       DetId                                       theDetId;
00087       std::vector< edm::Ptr< L1TkCluster< T > > > theClusters;
00088       edm::Ptr< SimTrack >                        theSimTrack;
00089       int                                         theDisplacement;
00090       int                                         theOffset;
00091 
00092   }; 
00093 
00100 
00101   template< typename T >
00102   L1TkStub< T >::L1TkStub()
00103   {
00105     theDetId = 0;
00106     theClusters.clear();
00107     theDisplacement = 999999;
00108     theOffset = 0;
00109   }
00110 
00112   template< typename T >
00113   L1TkStub< T >::L1TkStub( DetId aDetId )
00114   {
00116     theDetId = aDetId;
00117     theClusters.clear();
00118     theDisplacement = 999999;
00119     theOffset = 0;
00120   }
00121 
00123   template< typename T >
00124   L1TkStub< T >::~L1TkStub(){}
00125 
00127   template< typename T >
00128   std::vector< edm::Ptr< L1TkCluster< T > > > L1TkStub< T >::getClusterPtrs() const { return theClusters; }
00129 
00131   template< typename T >
00132   const edm::Ptr< L1TkCluster< T > >& L1TkStub< T >::getClusterPtr( unsigned int hitIdentifier ) const
00133   {
00134     typename std::vector< edm::Ptr< L1TkCluster< T > > >::const_iterator clusIter;
00135     for ( clusIter = theClusters.begin();
00136           clusIter != theClusters.end();
00137           ++clusIter )
00138     {
00139       if ( (*clusIter)->getStackMember() == hitIdentifier )
00140         return *clusIter;
00141     }
00142     return edm::Ptr< L1TkCluster< T > >();
00143   }
00144 
00146   template< typename T >
00147   void L1TkStub< T >::addClusterPtr( edm::Ptr< L1TkCluster< T > > aL1TkCluster )
00148   {
00152     theClusters.push_back( aL1TkCluster );
00153   }
00154 
00156   template< typename T >
00157   DetId L1TkStub< T >::getDetId() const { return theDetId; }
00158 
00159   template< typename T >
00160   void L1TkStub< T >::setDetId( DetId aDetId ) { theDetId = aDetId; }
00161 
00163   template< typename T >
00164   double L1TkStub< T >::getTriggerDisplacement() const { return 0.5*theDisplacement; }
00165 
00166   template< typename T >
00167   void L1TkStub< T >::setTriggerDisplacement( int aDisplacement ) { theDisplacement = aDisplacement; }
00168 
00169   template< typename T >
00170   double L1TkStub< T >::getTriggerOffset() const { return 0.5*theOffset; }
00171 
00172   template< typename T >
00173   void L1TkStub< T >::setTriggerOffset( int anOffset ) { theOffset = anOffset; }
00174 
00176   template< typename T >
00177   edm::Ptr< SimTrack > L1TkStub< T >::getSimTrackPtr() const { return theSimTrack; }
00178 
00179   template< typename T >
00180   bool L1TkStub< T >::isGenuine() const
00181   {
00182     /*
00186     if ( theClusters.at(0)->isUnknown() || theClusters.at(1)->isUnknown() )
00189       return false;
00190 
00191     else
00192     {
00196       if ( theClusters.at(0)->isGenuine() && theClusters.at(1)->isGenuine() )
00197       {
00198         if ( theClusters.at(0)->findSimTrackId() == theClusters.at(1)->findSimTrackId() )
00200           return true;
00201         else
00202           return false;
00203       }
00204       else
00205       {
00207         int prevTrack = -99999; // SimTrackId storage
00208         std::vector< edm::Ptr< SimTrack > > innerSimTracks = theClusters.at(0)->getSimTrackPtrs();
00209         std::vector< edm::Ptr< SimTrack > > outerSimTracks = theClusters.at(1)->getSimTrackPtrs();
00210         for ( unsigned int i = 0; i < innerSimTracks.size(); i++ )
00211         {
00213           if ( innerSimTracks.at(i).isNull() );
00214             continue;
00215           for ( unsigned int j = 0; j < outerSimTracks.size(); j++ )
00216           {
00218             if ( outerSimTracks.at(j).isNull() );
00219               continue;
00220 
00221             if ( innerSimTracks.at(i)->trackId() == outerSimTracks.at(j)->trackId() )
00222             {
00224               if ( prevTrack < 0 )
00225                 prevTrack = outerSimTracks.at(j)->trackId();
00226 
00227               if ( prevTrack != (int)outerSimTracks.at(j)->trackId() )
00230                 return false;
00231             }
00232           }
00233         }
00234         if ( prevTrack < 0 )
00236           return false;
00237         else
00242           return true;
00243       }
00244     }
00247     std::cerr << "W A R N I N G! L1TkStub::isGenuine() \t we should never get here" << std::endl;
00248     return true;
00249     */
00250 
00251     if ( theSimTrack.isNull() )
00252       return false;
00253 
00254     return true;
00255   }
00256 
00257   template< typename T >
00258   bool L1TkStub< T >::isCombinatoric() const
00259   {
00260     if ( this->isGenuine() )
00261       return false;
00262 
00266     if ( theClusters.at(0)->isUnknown() && theClusters.at(1)->isUnknown() )
00268       return false;
00269 
00270     else if ( theClusters.at(0)->isUnknown() || theClusters.at(1)->isUnknown() )
00273       return true;
00274 
00275     else
00276     {
00280       if ( theClusters.at(0)->isGenuine() && theClusters.at(1)->isGenuine() )
00281       {
00282         if ( theClusters.at(0)->findSimTrackId() == theClusters.at(1)->findSimTrackId() )
00284           return false;
00285         else
00286           return true;
00287       }
00288       else
00289       {
00291         int prevTrack = -99999; // SimTrackId storage
00292         std::vector< edm::Ptr< SimTrack > > innerSimTracks = theClusters.at(0)->getSimTrackPtrs();
00293         std::vector< edm::Ptr< SimTrack > > outerSimTracks = theClusters.at(1)->getSimTrackPtrs();
00294         for ( unsigned int i = 0; i < innerSimTracks.size(); i++ )
00295         {
00297           if ( innerSimTracks.at(i).isNull() );
00298             continue;
00299           for ( unsigned int j = 0; j < outerSimTracks.size(); j++ )
00300           {
00302             if ( outerSimTracks.at(j).isNull() );
00303               continue;
00304 
00305             if ( innerSimTracks.at(i)->trackId() == outerSimTracks.at(j)->trackId() )
00306             {
00308               if ( prevTrack < 0 )
00309                 prevTrack = outerSimTracks.at(j)->trackId();
00310 
00311               if ( prevTrack != (int)outerSimTracks.at(j)->trackId() )
00314                 return true;
00315             }
00316           }
00317         }
00318         if ( prevTrack < 0 )
00320           return true;
00321         else
00326           return false;
00327       }
00328     }
00331     std::cerr << "W A R N I N G! L1TkStub::isCombinatoric() \t we should never get here" << std::endl;
00332     return false; 
00333   }
00334 
00335   template< typename T >
00336   bool L1TkStub< T >::isUnknown() const
00337   {
00339     return ( theClusters.at(0)->isUnknown() && theClusters.at(1)->isUnknown() );
00340   }
00341 
00342   template< typename T >
00343   int L1TkStub< T >::findType() const
00344   {
00345     if ( theSimTrack.isNull() )
00346       return 999999999;
00347     return theSimTrack->type();
00348   }
00349 
00350   template< typename T >
00351   unsigned int L1TkStub< T >::findSimTrackId() const
00352   {
00353     if ( theSimTrack.isNull() )
00354       return 0;
00355     return theSimTrack->trackId();
00356   }
00357 
00359   template< typename T >
00360   void L1TkStub< T >::checkSimTrack()
00361   {
00364 
00368     if ( theClusters.at(0)->isUnknown() || theClusters.at(1)->isUnknown() )
00372       return;
00373 
00374     else
00375     {
00379       if ( theClusters.at(0)->isGenuine() && theClusters.at(1)->isGenuine() )
00380       {
00381         if ( theClusters.at(0)->findSimTrackId() == theClusters.at(1)->findSimTrackId() )
00382         {
00384           std::vector< edm::Ptr< SimTrack > > curSimTracks = theClusters.at(0)->getSimTrackPtrs();
00385           for ( unsigned int k = 0; k < curSimTracks.size(); k++ )
00386           {
00387             if ( curSimTracks.at(k).isNull() == false )
00388             {
00389               theSimTrack = curSimTracks.at(k);
00390               return;
00391             }
00392           }
00393         }
00394         else
00395           return;
00396       }
00397       else
00398       {
00400         int prevTrack = -99999; // SimTrackId storage
00401         unsigned int whichSimTrack = 0;
00402         std::vector< edm::Ptr< SimTrack > > innerSimTracks = theClusters.at(0)->getSimTrackPtrs();
00403         std::vector< edm::Ptr< SimTrack > > outerSimTracks = theClusters.at(1)->getSimTrackPtrs();
00404         for ( unsigned int i = 0; i < innerSimTracks.size(); i++ )
00405         {
00407           if ( innerSimTracks.at(i).isNull() );
00408             continue;
00409           for ( unsigned int j = 0; j < outerSimTracks.size(); j++ )
00410           {
00412             if ( outerSimTracks.at(j).isNull() );
00413               continue;
00414 
00415             if ( innerSimTracks.at(i)->trackId() == outerSimTracks.at(j)->trackId() )
00416             {
00418               if ( prevTrack < 0 )
00419               {
00420                 prevTrack = outerSimTracks.at(j)->trackId();
00421                 whichSimTrack = j;
00422               }
00423 
00424               if ( prevTrack != (int)outerSimTracks.at(j)->trackId() )
00427                 return;
00428             }
00429           }
00430         }
00431         if ( prevTrack < 0 )
00433           return;
00434         else
00439           theSimTrack = outerSimTracks.at(whichSimTrack);
00440       }
00441     }
00442   }
00443 
00444 
00446   template< typename T >
00447   std::string L1TkStub< T >::print( unsigned int i ) const
00448   {
00449     std::string padding("");
00450     for ( unsigned int j=0; j!=i; ++j )
00451       padding+="\t";
00452     std::stringstream output;
00453     output<<padding<<"L1TkStub:\n";
00454     padding+='\t';
00455     output << padding << "DetId: " << theDetId.rawId() << '\n';
00456     unsigned int iClu = 0;
00457     typename std::vector< edm::Ptr< L1TkCluster< T > > >::const_iterator clusIter;
00458     for ( clusIter = theClusters.begin(); clusIter!= theClusters.end(); ++clusIter )
00459       output << padding << "cluster: " << iClu++ << ", member: " << (*clusIter)->getStackMember() << ", cluster size: " << (*clusIter)->getHits().size() << '\n';
00460     return output.str();
00461   }
00462 
00463   template< typename T >
00464   std::ostream& operator << (std::ostream& os, const L1TkStub< T >& aL1TkStub) { return ( os<<aL1TkStub.print() ); }
00465 
00466 
00467 
00468 #endif
00469