Skip to content
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

CTPPS: tracking for the diamond detector #17817

Merged
merged 23 commits into from Mar 30, 2017
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
Show all changes
23 commits
Select commit Hold shift + click to select a range
5d39c91
Added the CTPPSDiamondLocalTrack object
forthommel Jan 31, 2017
17dfa01
New cfi file for the local tracks fitter
forthommel Jan 31, 2017
26dd3ac
first attempt of reco
nminafra Jan 16, 2017
d8d0af1
Added setters to CTPPS*LocalTrack and RecHit
nminafra Jan 16, 2017
c6926bb
New constructor for the track recognition algorithm object ; plugging…
forthommel Jan 31, 2017
7c3bd1f
Added the fitting parameters to the configuration file
forthommel Jan 31, 2017
7caf7bf
Combined both the totemRP and ctppsDiamond local reconstructions into…
forthommel Feb 8, 2017
b5120ae
Making use of the already present recoCTPPS object (instead of defini…
forthommel Feb 8, 2017
ffe12e6
Added the local tracks to the output definition
forthommel Feb 8, 2017
b2a13e3
Reco good to go
nminafra Feb 16, 2017
eb19b84
Changed variable names
nminafra Feb 16, 2017
10098d9
Storing local tracks information as floats instead of doubles ; match…
forthommel Feb 16, 2017
bbff525
Code documentation and cleanup
forthommel Feb 18, 2017
e6c6b05
Setting the default value to 0 for timing information in the local tr…
forthommel Feb 18, 2017
163c98c
Reco out of time bug fixed
nminafra Feb 20, 2017
921564a
RecHit and LocalTrack MH included
nminafra Feb 27, 2017
46800ae
Amended the diamonds local track recognition ; New local tracks const…
forthommel Mar 16, 2017
51530aa
Improved memory footprint of the local tracks recognition
forthommel Mar 18, 2017
3fc6e15
Use the precise pixel efficiency function instead of the fast one
forthommel Mar 20, 2017
fd4dd25
Passing by references ; range for loops
forthommel Mar 22, 2017
f45af97
Tracks sorting by temporal, then by spatial position
forthommel Mar 27, 2017
edf1427
Track recognition algorithms not defined as pointers
forthommel Mar 28, 2017
1744996
Unordered maps for the parameters and hitmp in local tracks fitter ; …
forthommel Mar 30, 2017
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Jump to
Jump to file
Failed to load files.
Diff view
Diff view
101 changes: 101 additions & 0 deletions DataFormats/CTPPSReco/interface/CTPPSDiamondLocalTrack.h
@@ -0,0 +1,101 @@
/****************************************************************************
*
* This is a part of CTPPS offline software.
* Authors:
* Laurent Forthomme (laurent.forthomme@cern.ch)
* Nicola Minafra nicola.minafra@cern.ch)
*
****************************************************************************/

#ifndef DataFormats_CTPPSReco_CTPPSDiamondLocalTrack
#define DataFormats_CTPPSReco_CTPPSDiamondLocalTrack

#include "DataFormats/Math/interface/Point3D.h"
#include "DataFormats/CTPPSReco/interface/CTPPSDiamondRecHit.h"

//----------------------------------------------------------------------------------------------------

class CTPPSDiamondLocalTrack
{
public:
CTPPSDiamondLocalTrack() :
chi_squared_( 0. ), valid_( true ), t_( 0. ), t_sigma_( 0. ), ts_index_( 0 ), mh_( 0 ) {}
CTPPSDiamondLocalTrack( const math::XYZPoint& pos0, const math::XYZPoint& pos0_sigma, float chisq, float t, float t_sigma, int oot_idx, int mult_hits ) :
pos0_( pos0 ), pos0_sigma_( pos0_sigma ),
chi_squared_( chisq ), valid_( false ),
t_( t ), t_sigma_( t_sigma ), ts_index_( oot_idx ), mh_( mult_hits ) {}
virtual ~CTPPSDiamondLocalTrack() {}

//--- spatial get'ters

inline float getX0() const { return pos0_.x(); }
inline float getX0Sigma() const { return pos0_sigma_.x(); }

inline float getY0() const { return pos0_.y(); }
inline float getY0Sigma() const { return pos0_sigma_.y(); }

inline float getZ0() const { return pos0_.z(); }

inline float getChiSquared() const { return chi_squared_; }

//--- spatial set'ters

inline void setPosition( const math::XYZPoint& pos0 ) { pos0_ = pos0; }
inline void setPositionSigma( const math::XYZPoint& pos0_sigma ) { pos0_sigma_ = pos0_sigma; }

inline void setChiSquared( const float chisq ) { chi_squared_ = chisq; }

inline bool isValid() const { return valid_; }
inline void setValid( bool valid ) { valid_ = valid; }

//--- temporal get'ters

inline float getT() const { return t_; }
inline float getTSigma() const { return t_sigma_; }

//--- temporal set'ters

inline void setT( const float t ) { t_ = t; }
inline void setTSigma( const float t_sigma ) { t_sigma_ = t_sigma; }

inline void setOOTIndex( const int i ) { ts_index_ = i; }
inline int getOOTIndex() const { return ts_index_; }

inline void setMultipleHits( const int i ) { mh_ = i; }
inline int getMultipleHits() const { return mh_; }

private:
//--- spatial information

/// initial track position
math::XYZPoint pos0_;
/// error on the initial track position
math::XYZPoint pos0_sigma_;
Copy link
Contributor

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

it's a bit odd to see diagonal-only elements of uncertainty.
Don't you need to have a more complete covariance, e.g.
typedef math::ErrorF<3>::type CovMatrix; ?

Copy link
Contributor

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

looking at the implementation, it's unclear to me if z is actually needed.
Is there going to be another reconstruction algo which sets it?

Copy link
Contributor

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Yes, it is needed for the proton reconstruction - this detector samples the proton trajectory at the given z position from the IP.

Copy link
Contributor

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

For the covariance matrix - Laurent and Nicola correct me if I'm wrong - as the sensor has rectangular pads with edges parallel to x and y, I can't see how one could get non-trivial off-diagonal elements. But I'm not expert.

Copy link
Contributor

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

What Jan says is true, moreover the size of the pads is ~mm, while the detector is parallel to x and y with an accuracy of <mRad. Hence the non-diagonal terms will always be negligible (~um) with respect to the precision of the tracks (~100 um).


/// fit chi^2
float chi_squared_;

/// fit valid?
bool valid_;

//--- timing information

float t_;
float t_sigma_;
/// Time slice index
int ts_index_;
/// Multiple hits counter
int mh_;

};

inline bool operator<( const CTPPSDiamondLocalTrack& lhs, const CTPPSDiamondLocalTrack& rhs )
{
// start to sort by temporal coordinate
if ( lhs.getT() < rhs.getT() ) return true;
if ( lhs.getT() > rhs.getT() ) return false;
// then sort by x-position
return ( lhs.getX0() < rhs.getX0() );
}

#endif
10 changes: 7 additions & 3 deletions DataFormats/CTPPSReco/interface/CTPPSDiamondRecHit.h
Expand Up @@ -19,12 +19,12 @@ class CTPPSDiamondRecHit
CTPPSDiamondRecHit() :
x_( 0. ), x_width_( 0. ), y_( 0. ), y_width_( 0. ),
t_( 0. ), tot_( 0. ),
ts_index_( 0 ), hptdc_err_( 0 )
ts_index_( 0 ), hptdc_err_( 0 ), mh_( false )
{}
CTPPSDiamondRecHit( float x, float x_width, float y, float y_width, float t, float tot, int oot_idx, const HPTDCErrorFlags& hptdc_err ) :
CTPPSDiamondRecHit( float x, float x_width, float y, float y_width, float t, float tot, int oot_idx, const HPTDCErrorFlags& hptdc_err, const bool mh ) :
x_( x ), x_width_( x_width ), y_( y ), y_width_( y_width ),
t_( t ), tot_( tot ),
ts_index_( oot_idx ), hptdc_err_( hptdc_err )
ts_index_( oot_idx ), hptdc_err_( hptdc_err ), mh_( mh )
{}

inline void setX( const float& x ) { x_ = x; }
Expand All @@ -47,6 +47,9 @@ class CTPPSDiamondRecHit

inline void setOOTIndex( const int& i ) { ts_index_ = i; }
inline int getOOTIndex() const { return ts_index_; }

inline void setMultipleHits( const bool mh ) { mh_ = mh; }
inline bool getMultipleHits() const { return mh_; }

inline void setHPTDCErrorFlags( const HPTDCErrorFlags& err ) { hptdc_err_ = err; }
inline HPTDCErrorFlags getHPTDCErrorFlags() const { return hptdc_err_; }
Expand All @@ -58,6 +61,7 @@ class CTPPSDiamondRecHit
/// Time slice index
int ts_index_;
HPTDCErrorFlags hptdc_err_;
bool mh_;
};

//----------------------------------------------------------------------------------------------------
Expand Down
18 changes: 15 additions & 3 deletions DataFormats/CTPPSReco/src/classes.h
Expand Up @@ -10,6 +10,7 @@
#include "DataFormats/CTPPSReco/interface/TotemRPLocalTrack.h"

#include "DataFormats/CTPPSReco/interface/CTPPSDiamondRecHit.h"
#include "DataFormats/CTPPSReco/interface/CTPPSDiamondLocalTrack.h"

#include "DataFormats/CTPPSReco/interface/CTPPSLocalTrackLite.h"

Expand Down Expand Up @@ -52,11 +53,22 @@ namespace DataFormats_CTPPSReco {
CTPPSDiamondRecHit ctd_rh;
edm::Ptr<CTPPSDiamondRecHit> ptr_ctd_rh;
edm::Wrapper<CTPPSDiamondRecHit> wrp_ctd_rh;
std::vector<CTPPSDiamondRecHit> vec_rh;
std::vector< edm::DetSet<CTPPSDiamondRecHit> > vec_ds_rh;
edm::DetSet<CTPPSDiamondRecHit> ds_rh;
std::vector<CTPPSDiamondRecHit> vec_ctd_rh;
std::vector< edm::DetSet<CTPPSDiamondRecHit> > vec_ds_ctd_rh;
edm::DetSet<CTPPSDiamondRecHit> ds_ctd_rh;
edm::DetSetVector<CTPPSDiamondRecHit> dsv_ctd_rh;
edm::Wrapper< edm::DetSetVector<CTPPSDiamondRecHit> > wrp_dsv_ctd_rh;
edm::Wrapper< std::vector<CTPPSDiamondRecHit> > wrp_vec_ctd_rh;

CTPPSDiamondLocalTrack ctd_lt;
edm::Ptr<CTPPSDiamondLocalTrack> ptr_ctd_lt;
edm::Wrapper<CTPPSDiamondLocalTrack> wrp_ctd_lt;
std::vector<CTPPSDiamondLocalTrack> vec_ctd_lt;
edm::DetSet<CTPPSDiamondLocalTrack> ds_ctd_lt;
std::vector< edm::DetSet<CTPPSDiamondLocalTrack> > vec_ds_ctd_lt;
edm::Wrapper< std::vector<CTPPSDiamondLocalTrack> > wrp_vec_ctd_lt;
edm::DetSetVector<CTPPSDiamondLocalTrack> dsv_ctd_lt;
edm::Wrapper<edm::DetSetVector<CTPPSDiamondLocalTrack> > wrp_dsv_ctd_lt;

//--- common objects

Expand Down
16 changes: 15 additions & 1 deletion DataFormats/CTPPSReco/src/classes_def.xml
Expand Up @@ -45,8 +45,9 @@

<!-- diamonds objects -->

<class name="CTPPSDiamondRecHit" ClassVersion="3">
<class name="CTPPSDiamondRecHit" ClassVersion="4">
<version ClassVersion="3" checksum="2377257106"/>
<version ClassVersion="4" checksum="2150979176"/>
</class>
<class name="edm::Ptr<CTPPSDiamondRecHit>"/>
<class name="edm::Wrapper<CTPPSDiamondRecHit>"/>
Expand All @@ -55,6 +56,19 @@
<class name="edm::DetSet<CTPPSDiamondRecHit>"/>
<class name="edm::DetSetVector<CTPPSDiamondRecHit>"/>
<class name="edm::Wrapper< edm::DetSetVector<CTPPSDiamondRecHit> >"/>
<class name="edm::Wrapper< std::vector<CTPPSDiamondRecHit> >"/>

<class name="CTPPSDiamondLocalTrack" ClassVersion="3">
<version ClassVersion="3" checksum="2739193553"/>
</class>
<class name="edm::Ptr<CTPPSDiamondLocalTrack>"/>
<class name="edm::Wrapper<CTPPSDiamondLocalTrack>"/>
<class name="std::vector<CTPPSDiamondLocalTrack>"/>
<class name="edm::DetSet<CTPPSDiamondLocalTrack>"/>
<class name="std::vector< edm::DetSet<CTPPSDiamondLocalTrack> >"/>
<class name="edm::Wrapper< std::vector<CTPPSDiamondLocalTrack> >"/>
<class name="edm::DetSetVector<CTPPSDiamondLocalTrack>"/>
<class name="edm::Wrapper< edm::DetSetVector<CTPPSDiamondLocalTrack> >"/>

<!-- common objects -->

Expand Down
3 changes: 3 additions & 0 deletions RecoCTPPS/Configuration/python/RecoCTPPS_EventContent_cff.py
Expand Up @@ -19,6 +19,7 @@
'keep CTPPSDiamondDigiedmDetSetVector_ctppsDiamondRawToDigi_*_*',
'keep TotemVFATStatusedmDetSetVector_ctppsDiamondRawToDigi_*_*',
'keep CTPPSDiamondRecHitedmDetSetVector_ctppsDiamondRecHits_*_*',
'keep CTPPSDiamondLocalTrackedmDetSetVector_ctppsDiamondLocalTracks_*_*',

# CTPPS common
'keep CTPPSLocalTrackLites_ctppsLocalTrackLiteProducer_*_*'
Expand All @@ -45,6 +46,7 @@
'keep CTPPSDiamondDigiedmDetSetVector_ctppsDiamondRawToDigi_*_*',
'keep TotemVFATStatusedmDetSetVector_ctppsDiamondRawToDigi_*_*',
'keep CTPPSDiamondRecHitedmDetSetVector_ctppsDiamondRecHits_*_*',
'keep CTPPSDiamondLocalTrackedmDetSetVector_ctppsDiamondLocalTracks_*_*',

# CTPPS common
'keep CTPPSLocalTrackLites_ctppsLocalTrackLiteProducer_*_*'
Expand All @@ -71,6 +73,7 @@
'keep CTPPSDiamondDigiedmDetSetVector_ctppsDiamondRawToDigi_*_*',
'keep TotemVFATStatusedmDetSetVector_ctppsDiamondRawToDigi_*_*',
'keep CTPPSDiamondRecHitedmDetSetVector_ctppsDiamondRecHits_*_*',
'keep CTPPSDiamondLocalTrackedmDetSetVector_ctppsDiamondLocalTracks_*_*',

# CTPPS common
'keep CTPPSLocalTrackLites_ctppsLocalTrackLiteProducer_*_*'
Expand Down
73 changes: 73 additions & 0 deletions RecoCTPPS/TotemRPLocal/interface/CTPPSDiamondTrackRecognition.h
@@ -0,0 +1,73 @@
/****************************************************************************
*
* This is a part of CTPPS offline software.
* Authors:
* Laurent Forthomme (laurent.forthomme@cern.ch)
* Nicola Minafra (nicola.minafra@cern.ch)
*
****************************************************************************/

#ifndef RecoCTPPS_TotemRPLocal_CTPPSDiamondTrackRecognition
#define RecoCTPPS_TotemRPLocal_CTPPSDiamondTrackRecognition

#include "FWCore/ParameterSet/interface/ParameterSet.h"

#include "DataFormats/Common/interface/DetSet.h"
#include "DataFormats/Common/interface/DetSetVector.h"

#include "DataFormats/CTPPSReco/interface/CTPPSDiamondRecHit.h"
#include "DataFormats/CTPPSReco/interface/CTPPSDiamondLocalTrack.h"

#include <vector>
#include <unordered_map>
#include "TF1.h"

/**
* \brief Class performing smart reconstruction for CTPPS Diamond Detectors.
* \date Jan 2017
**/
class CTPPSDiamondTrackRecognition
{
public:
CTPPSDiamondTrackRecognition( const edm::ParameterSet& );
~CTPPSDiamondTrackRecognition();

/// Reset the list of hits
void clear();

/// Feed a new hit to the tracks recognition algorithm
void addHit( const CTPPSDiamondRecHit& recHit );

/// Produce a collection of tracks for the current station, given its hits collection
int produceTracks( edm::DetSet<CTPPSDiamondLocalTrack>& tracks );

private:
struct HitParameters {
HitParameters( const float center, const float width ) :
center( center ), width( width ) {}
float center;
float width;
};
typedef std::vector<HitParameters> HitParametersVector;
typedef std::unordered_map<int,HitParametersVector> HitParametersVectorMap;

/// Default hit function accounting for the pad spatial efficiency
static const std::string pixelEfficiencyDefaultFunction_;

const float threshold_;
const float thresholdFromMaximum_;
const float resolution_;
const float sigma_;
const float startFromX_;
const float stopAtX_;

float yPosition_;
float yWidth_;

/// Function for pad efficiency
TF1 hit_f_;
HitParametersVectorMap hitParametersVectorMap_;
std::unordered_map<int,int> mhMap_;
};

#endif