New issue
Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.
By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.
Already on GitHub? Sign in to your account
Migration of L1 tracker primitives from SLHCDEV #13696
Merged
davidlange6
merged 19 commits into
cms-sw:CMSSW_8_1_X
from
delaere:L1primitivesMigration_step1
May 3, 2016
Merged
Changes from all commits
Commits
Show all changes
19 commits
Select commit
Hold shift + click to select a range
aa77a65
step 1: Dataformat and L1Trigger code
delaere 12965ad
Fix BuildFile and class_def. Temporary checksums.
delaere d0861be
Merged L1primitivesMigration_step1 from repository delaere
delaere a6dc7ef
some cleaning of classes.h and classes_def.xml
delaere 59a5889
Merged L1primitivesMigration_step1 from repository delaere
delaere 48152ef
adapt to the change of case in the TrackerTopology methods
delaere cc7b342
fix the classes_def.xml
delaere 4f1ccc3
Merged L1primitivesMigration_step1 from repository delaere
delaere 0e569f9
Cleaning and proper configuration of TTClusterBuilder
delaere 2b25d6d
Fix the TTCluster and TTStub builders
delaere af8fff5
First step in adapting the Stub reco algos to new numbering
delaere 3adfbb2
Reenable tab2013 stub algo.
delaere 80e7e7f
Small adaptation to port TrackTrigger devs to 81X
delaere bedff85
Merged L1primitivesMigration_step1 from repository delaere
delaere 8a14745
Review the implementation of TTCluster and TTStub
delaere 5e3b101
Work on Cluster and Stub builders. More to come.
delaere eadb300
replace boost::shared_ptr by the std version
delaere 7d078bf
Comments from David...
delaere e067815
Cleaning of unused algorithms
delaere File filter
Filter by extension
Conversations
Failed to load comments.
Jump to
Jump to file
Failed to load files.
Diff view
Diff view
There are no files selected for viewing
This file contains bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
Original file line number | Diff line number | Diff line change |
---|---|---|
@@ -0,0 +1,9 @@ | ||
<use name="DataFormats/Common"/> | ||
<use name="DataFormats/DetId"/> | ||
<use name="DataFormats/GeometryCommonDetAlgo"/> | ||
<use name="DataFormats/GeometryVector"/> | ||
<use name="DataFormats/Phase2TrackerDigi"/> | ||
<export> | ||
<lib name="1"/> | ||
</export> | ||
|
This file contains bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
Original file line number | Diff line number | Diff line change |
---|---|---|
@@ -0,0 +1,185 @@ | ||
/*! \class TTCluster | ||
* \brief Class to store the L1 Track Trigger clusters | ||
* \details After moving from SimDataFormats to DataFormats, | ||
* the template structure of the class was maintained | ||
* in order to accomodate any types other than Phase2TrackerDigis | ||
* in case there is such a need in the future. | ||
* | ||
* \author Nicola Pozzobon | ||
* \author Emmanuele Salvati | ||
* \date 2013, Jul 12 | ||
* | ||
*/ | ||
|
||
#ifndef L1_TRACK_TRIGGER_CLUSTER_FORMAT_H | ||
#define L1_TRACK_TRIGGER_CLUSTER_FORMAT_H | ||
|
||
#include "DataFormats/Common/interface/Ref.h" | ||
#include "DataFormats/Common/interface/Ptr.h" | ||
#include "DataFormats/Common/interface/DetSet.h" | ||
#include "DataFormats/Common/interface/DetSetVector.h" | ||
#include "DataFormats/DetId/interface/DetId.h" | ||
#include "DataFormats/Phase2TrackerDigi/interface/Phase2TrackerDigi.h" | ||
#include "DataFormats/GeometryCommonDetAlgo/interface/MeasurementPoint.h" | ||
#include "DataFormats/GeometryVector/interface/GlobalPoint.h" /// NOTE: this is needed even if it seems not | ||
|
||
template< typename T > | ||
class TTCluster | ||
{ | ||
public: | ||
/// Constructors | ||
TTCluster(); | ||
TTCluster( std::vector< T > aHits, | ||
DetId aDetId, | ||
unsigned int aStackMember, | ||
bool storeLocal ); | ||
|
||
/// Destructor | ||
~TTCluster(); | ||
|
||
/// Data members: getABC( ... ) | ||
/// Helper methods: findABC( ... ) | ||
|
||
/// Hits in the Cluster | ||
std::vector< T > getHits() const { return theHits; } | ||
void setHits( std::vector< T > aHits ) { theHits = aHits; } | ||
|
||
/// Detector element | ||
DetId getDetId() const { return theDetId; } | ||
void setDetId( DetId aDetId ) { theDetId = aDetId; } | ||
unsigned int getStackMember() const { return theStackMember; } | ||
void setStackMember( unsigned int aStackMember ) { theStackMember = aStackMember; } | ||
|
||
/// Rows and columns to get rid of Digi collection | ||
std::vector< int > findRows() const; | ||
std::vector< int > findCols() const; | ||
void setCoordinates( std::vector< int > a, std::vector< int > b ) { theRows = a; theCols = b; } | ||
std::vector< int > getRows() const { return theRows; } | ||
std::vector< int > getCols() const { return theCols; } | ||
|
||
/// Cluster width | ||
unsigned int findWidth() const; | ||
|
||
/// Single hit coordinates | ||
/// Average cluster coordinates | ||
MeasurementPoint findHitLocalCoordinates( unsigned int hitIdx ) const; | ||
MeasurementPoint findAverageLocalCoordinates() const; | ||
|
||
/// Information | ||
std::string print( unsigned int i = 0 ) const; | ||
|
||
private: | ||
/// Data members | ||
std::vector< T > theHits; | ||
DetId theDetId; | ||
unsigned int theStackMember; | ||
|
||
std::vector< int > theRows; | ||
std::vector< int > theCols; | ||
|
||
}; /// Close class | ||
|
||
/*! \brief Implementation of methods | ||
* \details Here, in the header file, the methods which do not depend | ||
* on the specific type <T> that can fit the template. | ||
* Other methods, with type-specific features, are implemented | ||
* in the source file. | ||
*/ | ||
|
||
/// Default Constructor | ||
/// NOTE: to be used with setSomething(...) methods | ||
template< typename T > | ||
TTCluster< T >::TTCluster() | ||
{ | ||
/// Set default data members | ||
theHits.clear(); | ||
theDetId = 0; | ||
theStackMember = 0; | ||
|
||
theRows.clear(); | ||
theCols.clear(); | ||
} | ||
|
||
/// Another Constructor | ||
template< typename T > | ||
TTCluster< T >::TTCluster( std::vector< T > aHits, | ||
DetId aDetId, | ||
unsigned int aStackMember, | ||
bool storeLocal ) | ||
{ | ||
/// Set data members | ||
this->setHits( aHits ); | ||
this->setDetId( aDetId ); | ||
this->setStackMember( aStackMember ); | ||
|
||
theRows.clear(); | ||
theCols.clear(); | ||
if ( storeLocal ) | ||
{ | ||
this->setCoordinates( this->findRows(), this->findCols() ); | ||
} | ||
} | ||
|
||
/// Destructor | ||
template< typename T > | ||
TTCluster< T >::~TTCluster(){} | ||
|
||
/// Cluster width | ||
template< > | ||
unsigned int TTCluster< edm::Ref< edm::DetSetVector< Phase2TrackerDigi >, Phase2TrackerDigi > >::findWidth() const; | ||
|
||
/// Single hit coordinates | ||
/// Average cluster coordinates | ||
template< > | ||
MeasurementPoint TTCluster< edm::Ref< edm::DetSetVector< Phase2TrackerDigi >, Phase2TrackerDigi > >::findHitLocalCoordinates( unsigned int hitIdx ) const; | ||
|
||
template< > | ||
MeasurementPoint TTCluster< edm::Ref< edm::DetSetVector< Phase2TrackerDigi >, Phase2TrackerDigi > >::findAverageLocalCoordinates() const; | ||
|
||
/// Operations with coordinates stored locally | ||
template< typename T > | ||
std::vector< int > TTCluster< T >::findRows() const | ||
{ | ||
std::vector< int > temp; | ||
return temp; | ||
} | ||
|
||
template< typename T > | ||
std::vector< int > TTCluster< T >::findCols() const | ||
{ | ||
std::vector< int > temp; | ||
return temp; | ||
} | ||
|
||
template< > | ||
std::vector< int > TTCluster< edm::Ref< edm::DetSetVector< Phase2TrackerDigi >, Phase2TrackerDigi > >::findRows() const; | ||
|
||
template< > | ||
std::vector< int > TTCluster< edm::Ref< edm::DetSetVector< Phase2TrackerDigi >, Phase2TrackerDigi > >::findCols() const; | ||
|
||
/// Information | ||
template< typename T > | ||
std::string TTCluster< T >::print( unsigned int i ) const | ||
{ | ||
std::string padding(""); | ||
for ( unsigned int j = 0; j != i; ++j ) | ||
{ | ||
padding+="\t"; | ||
} | ||
|
||
std::stringstream output; | ||
output<<padding<<"TTCluster:\n"; | ||
padding+='\t'; | ||
output << padding << "DetId: " << theDetId.rawId() << '\n'; | ||
output << padding << "member: " << theStackMember << ", cluster size: " << theHits.size() << '\n'; | ||
return output.str(); | ||
} | ||
|
||
template< typename T > | ||
std::ostream& operator << ( std::ostream& os, const TTCluster< T >& aTTCluster ) | ||
{ | ||
return ( os << aTTCluster.print() ); | ||
} | ||
|
||
#endif | ||
|
This file contains bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
Original file line number | Diff line number | Diff line change |
---|---|---|
@@ -0,0 +1,171 @@ | ||
/*! \class TTStub | ||
* \brief Class to store the L1 Track Trigger stubs | ||
* \details After moving from SimDataFormats to DataFormats, | ||
* the template structure of the class was maintained | ||
* in order to accomodate any types other than Phase2TrackerDigis | ||
* in case there is such a need in the future. | ||
* | ||
* \author Andrew W. Rose | ||
* \author Nicola Pozzobon | ||
* \date 2013, Jul 12 | ||
* | ||
*/ | ||
|
||
#ifndef L1_TRACK_TRIGGER_STUB_FORMAT_H | ||
#define L1_TRACK_TRIGGER_STUB_FORMAT_H | ||
|
||
#include "DataFormats/GeometryVector/interface/GlobalVector.h" | ||
#include "DataFormats/L1TrackTrigger/interface/TTCluster.h" | ||
#include "DataFormats/Common/interface/DetSetVectorNew.h" | ||
|
||
template< typename T > | ||
class TTStub | ||
{ | ||
public: | ||
/// Constructors | ||
TTStub(); | ||
TTStub( DetId aDetId ); | ||
|
||
/// Destructor | ||
~TTStub(); | ||
|
||
/// Data members: getABC( ... ) | ||
/// Helper methods: findABC( ... ) | ||
|
||
/// Clusters composing the Stub | ||
const edm::Ref< edmNew::DetSetVector< TTCluster< T > >, TTCluster< T > >& getClusterRef( unsigned int hitIdentifier ) const; | ||
void addClusterRef( edm::Ref< edmNew::DetSetVector< TTCluster< T > >, TTCluster< T > > aTTCluster ); | ||
|
||
/// Detector element | ||
DetId getDetId() const { return theDetId; } | ||
void setDetId( DetId aDetId ) { theDetId = aDetId; } | ||
|
||
/// Trigger information | ||
double getTriggerDisplacement() const; /// In FULL-STRIP units! (hence, not implemented herein) | ||
void setTriggerDisplacement( int aDisplacement ); /// In HALF-STRIP units! | ||
double getTriggerOffset() const; /// In FULL-STRIP units! (hence, not implemented herein) | ||
void setTriggerOffset( int anOffset ); /// In HALF-STRIP units! | ||
|
||
/// CBC3-style trigger information | ||
/// for sake of simplicity, these methods are | ||
/// slightly out of the getABC(...)/findABC(...) rule | ||
double getTriggerPosition() const; /// In FULL-STRIP units! | ||
double getTriggerBend() const; /// In FULL-STRIP units! | ||
|
||
/// Information | ||
std::string print( unsigned int i = 0 ) const; | ||
|
||
private: | ||
/// Data members | ||
DetId theDetId; | ||
edm::Ref< edmNew::DetSetVector< TTCluster< T > >, TTCluster< T > > theClusterRef0; | ||
edm::Ref< edmNew::DetSetVector< TTCluster< T > >, TTCluster< T > > theClusterRef1; | ||
int theDisplacement; | ||
int theOffset; | ||
|
||
}; /// Close class | ||
|
||
/*! \brief Implementation of methods | ||
* \details Here, in the header file, the methods which do not depend | ||
* on the specific type <T> that can fit the template. | ||
* Other methods, with type-specific features, are implemented | ||
* in the source file. | ||
*/ | ||
|
||
/// Default Constructor | ||
template< typename T > | ||
TTStub< T >::TTStub() | ||
{ | ||
/// Set default data members | ||
theDetId = 0; | ||
theDisplacement = 999999; | ||
theOffset = 0; | ||
} | ||
|
||
/// Another Constructor | ||
template< typename T > | ||
TTStub< T >::TTStub( DetId aDetId ) | ||
{ | ||
/// Set data members | ||
this->setDetId( aDetId ); | ||
|
||
/// Set default data members | ||
theDisplacement = 999999; | ||
theOffset = 0; | ||
} | ||
|
||
/// Destructor | ||
template< typename T > | ||
TTStub< T >::~TTStub(){} | ||
|
||
/// Get the Reference to a Cluster | ||
template< typename T > | ||
const edm::Ref< edmNew::DetSetVector< TTCluster< T > >, TTCluster< T > >& TTStub< T >::getClusterRef( unsigned int hitIdentifier ) const | ||
{ | ||
return hitIdentifier==0 ? theClusterRef0 : theClusterRef1; | ||
} | ||
|
||
template< typename T > | ||
void TTStub< T >::addClusterRef( edm::Ref< edmNew::DetSetVector< TTCluster< T > >, TTCluster< T > > aTTCluster ) | ||
{ | ||
if(aTTCluster->getStackMember() == 0) theClusterRef0 = aTTCluster; | ||
else if (aTTCluster->getStackMember() == 1) theClusterRef1 = aTTCluster; | ||
} | ||
|
||
/// Trigger info | ||
template< typename T > | ||
double TTStub< T >::getTriggerDisplacement() const { return 0.5*theDisplacement; } | ||
|
||
template< typename T > | ||
void TTStub< T >::setTriggerDisplacement( int aDisplacement ) { theDisplacement = aDisplacement; } | ||
|
||
template< typename T > | ||
double TTStub< T >::getTriggerOffset() const { return 0.5*theOffset; } | ||
|
||
template< typename T > | ||
void TTStub< T >::setTriggerOffset( int anOffset ) { theOffset = anOffset; } | ||
|
||
/// CBC3-style trigger info | ||
template< typename T > | ||
double TTStub< T >::getTriggerPosition() const | ||
{ | ||
return this->getClusterRef(0)->findAverageLocalCoordinates().x(); | ||
} | ||
|
||
template< typename T > | ||
double TTStub< T >::getTriggerBend() const | ||
{ | ||
if ( theDisplacement == 999999 ) | ||
return theDisplacement; | ||
|
||
return 0.5*( theDisplacement - theOffset ); | ||
} | ||
|
||
/// Information | ||
template< typename T > | ||
std::string TTStub< T >::print( unsigned int i ) const | ||
{ | ||
std::string padding(""); | ||
for ( unsigned int j = 0; j != i; ++j ) | ||
{ | ||
padding+="\t"; | ||
} | ||
|
||
std::stringstream output; | ||
output<<padding<<"TTStub:\n"; | ||
padding+='\t'; | ||
output << padding << "DetId: " << theDetId.rawId() << ", position: " << this->getTriggerPosition(); | ||
output << ", bend: " << this->getTriggerBend() << '\n'; | ||
unsigned int iClu = 0; | ||
output << padding << "cluster 0: address: " << theClusterRef0.get(); | ||
output << ", cluster size: " << theClusterRef0->getHits().size() << '\n'; | ||
output << padding << "cluster 1: address: " << theClusterRef1.get(); | ||
output << ", cluster size: " << theClusterRef1->getHits().size() << '\n'; | ||
return output.str(); | ||
} | ||
|
||
template< typename T > | ||
std::ostream& operator << ( std::ostream& os, const TTStub< T >& aTTStub ) { return ( os << aTTStub.print() ); } | ||
|
||
#endif | ||
|
This file contains bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
Original file line number | Diff line number | Diff line change |
---|---|---|
@@ -0,0 +1,24 @@ | ||
/*! \brief Definition of all the relevant data types | ||
* | ||
* \author Andrew W. Rose | ||
* \author Nicola Pozzobon | ||
* \date 2013, Jul 12 | ||
* | ||
*/ | ||
|
||
#ifndef L1_TRACK_TRIGGER_TYPES_H | ||
#define L1_TRACK_TRIGGER_TYPES_H | ||
|
||
/// Standard CMS Formats | ||
#include "DataFormats/Common/interface/Ref.h" | ||
#include "DataFormats/Common/interface/Ptr.h" | ||
#include "DataFormats/Common/interface/DetSetVector.h" | ||
#include "DataFormats/Phase2TrackerDigi/interface/Phase2TrackerDigi.h" | ||
#include "DataFormats/L1TrackTrigger/interface/TTCluster.h" | ||
#include "DataFormats/L1TrackTrigger/interface/TTStub.h" | ||
|
||
/// The reference types | ||
typedef edm::Ref< edm::DetSetVector< Phase2TrackerDigi >, Phase2TrackerDigi > Ref_Phase2TrackerDigi_; | ||
|
||
#endif | ||
|
Oops, something went wrong.
Add this suggestion to a batch that can be applied as a single commit.
This suggestion is invalid because no changes were made to the code.
Suggestions cannot be applied while the pull request is closed.
Suggestions cannot be applied while viewing a subset of changes.
Only one suggestion per line can be applied in a batch.
Add this suggestion to a batch that can be applied as a single commit.
Applying suggestions on deleted lines is not supported.
You must change the existing code in this line in order to create a valid suggestion.
Outdated suggestions cannot be applied.
This suggestion has been applied or marked resolved.
Suggestions cannot be applied from pending reviews.
Suggestions cannot be applied on multi-line comments.
Suggestions cannot be applied while the pull request is queued to merge.
Suggestion cannot be applied right now. Please check back later.
There was a problem hiding this comment.
Choose a reason for hiding this comment
The reason will be displayed to describe this comment to others. Learn more.
@delaere - is the maximum size of these vectors small and known?
[given the size of the PR, expect things little by little]
There was a problem hiding this comment.
Choose a reason for hiding this comment
The reason will be displayed to describe this comment to others. Learn more.
Unfortunately not. In general the size is small but it can get very large in some cases. The only formal maximum is the number of strips in the module.