Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 11:57:20

0001 #ifndef Alignment_ReferenceTrajectories_ReferenceTrajectoryBase_H
0002 #define Alignment_ReferenceTrajectories_ReferenceTrajectoryBase_H
0003 
0004 /**
0005  * Author     : Gero Flucke (based on code for ORCA by Edmund Widl)
0006  * date       : 2006/09/17
0007  * last update: $Date: 2010/09/10 12:11:39 $
0008  * by         : $Author: mussgill $
0009  *
0010  * Base class for reference 'trajectories' of single- or multiparticles
0011  * stated.
0012  * Inheriting classes have to calculate all relevant quantities accessed
0013  * through member functions of this base class:
0014  *
0015  * The local measured x/y coordinates on all crossed detectors as vector:
0016  *
0017  * m = (x1, y1, x2, y2, ..., xN, yN) [transposed vector shown]
0018  *
0019  * their covariance matrix (possibly containing correlations between hits
0020  * due to material effects which are not taken into account by the ideal
0021  * trajectory parametrisation, cf. below),
0022  *
0023  * similarly the local x/y cordinates of the reference trajectory with covariance,
0024  *
0025  * the parameters of the (ideal) 'trajectory' 
0026  * (e.g. 5 values for a single track or 9 for a two-track-state with vertex constraint),
0027  *
0028  * the derivatives of the local coordinates of the reference trajectory
0029  * with respect to the initial 'trajectory'-parameters, 
0030  * e.g. for n parameters 'p' and N hits with 'x/y' coordinates:
0031  * 
0032  *  D = ( dx1/dp1, dx1/dp2, ..., dx1/dpn,
0033  *
0034  *        dy1/dp1, dy1/dp2, ..., dy1/dpn,
0035  *
0036  *        dx2/dp1, dx2/dp2, ..., dx2/dpn,
0037  *
0038  *        dy2/dp1, dy2/dp2, ..., dy2/dpn,
0039  *
0040  *           .        .             .
0041  *
0042  *           .        .             .
0043  *
0044  *        dxN/dp1, dxN/dp2, ..., dxN/dpn,
0045  *
0046  *        dyN/dp1, dyN/dp2, ..., dyN/dpn )
0047  *
0048  * and finally the TrajectoryStateOnSurface's of the reference trajectory.
0049  * 
0050  * Take care to check validity of the calculation (isValid()) before using the results.
0051  *
0052  *.............................................................................
0053  *
0054  * 090730 C. Kleinwort: 'Break Points' introduced for better description of multiple 
0055  *                      scattering (correlations)
0056  *
0057  * For each detector layer (in the hit collection) two ortogonal scattering angles are
0058  * introduced as new local track parameters ("break points"). To constrain those to the
0059  * expected mean (=0.) and variance (~1/p^2 X/X0) corresponding measurements are added.
0060  * Invalid hits may be used to produce break points too (UseInvalidHits = True).
0061  *
0062  * Break points have been implemented for: ReferenceTrajectory, BzeroReferenceTrajectory,
0063  *    DualReferenceTrajectory and DualBzeroReferenceTrajectory
0064  * For 'Dual' trajectories they make no sense as only the correlations between the hits  
0065  * in a track half are accounted for and not those between the halves! 
0066  *
0067  * Break Points are selected by TrajectoryFactory.MaterialEffects = "BreakPoints"
0068  *
0069  * 090909 C. Kleinwort: 'Broken Lines' introduced for description of multiple scattering
0070  *      (V. Blobel, Nuclear Instruments and Methods A, 566 (2006), pp. 14-17)
0071  * Fine Broken Lines are selected by TrajectoryFactory.MaterialEffects = "BrokenLinesFine"
0072  *      (exact derivatives)
0073  * Coarse Broken Lines are selected by TrajectoryFactory.MaterialEffects = "BrokenLines[Coarse]"
0074  *      (approximate derivatives, closeby hits (ds<1cm) combined)
0075  *
0076  * 091127 C. Kleinwort: 
0077  *   1) 'Broken Lines' extended to PCA,
0078  *      required for beamspot constraint and TwoBodyDecayTrajectory,
0079  *      selected with "BrokenLines[Coarse]Pca" or "BrokenLinesFinePca"
0080  *   2) For coarse Broken Lines linear interpolation is used for combined hits
0081  *   3) TwoBodyDecayTrajectory implemented for break points and Broken Lines
0082  *
0083  * 141103 C. Kleinwort: 'General Broken Lines' introduced for description of multiple scattering
0084  *       (C. Kleinwort, Nuclear Instruments and Methods A, 673 (2012), pp. 107-110)
0085  *        needs GBL version >= V01-13-00 (from svnsrv.desy.de)
0086  *        Selected by TrajectoryFactory.MaterialEffects = "LocalGBL" or = "CurvlinGBL"
0087  *        (for trajectory constructed in local or curvilinear system)
0088  */
0089 
0090 #include "DataFormats/GeometrySurface/interface/ReferenceCounted.h"
0091 #include "DataFormats/TrajectorySeed/interface/PropagationDirection.h"
0092 
0093 // for AlgebraicVector, -Matrix and -SymMatrix:
0094 #include "DataFormats/CLHEP/interface/AlgebraicObjects.h"
0095 #include "TrackingTools/TransientTrackingRecHit/interface/TransientTrackingRecHit.h"
0096 #include "TrackingTools/TrajectoryState/interface/TrajectoryStateOnSurface.h"
0097 #include "TrackingTools/TrajectoryParametrization/interface/LocalTrajectoryError.h"
0098 
0099 #include <vector>
0100 
0101 #include <Eigen/Dense>
0102 
0103 #include "GblTrajectory.h"
0104 
0105 class ReferenceTrajectoryBase : public ReferenceCounted {
0106 public:
0107   typedef ReferenceCountingPointer<ReferenceTrajectoryBase> ReferenceTrajectoryPtr;
0108 
0109   enum MaterialEffects {
0110     none,
0111     multipleScattering,
0112     energyLoss,
0113     combined,
0114     breakPoints,
0115     brokenLinesCoarse,
0116     brokenLinesFine,
0117     localGBL,
0118     curvlinGBL
0119   };
0120 
0121   struct Config {
0122     Config(MaterialEffects matEff,
0123            PropagationDirection direction,
0124            double m = -std::numeric_limits<double>::infinity(),
0125            double est = -std::numeric_limits<double>::infinity())
0126         : materialEffects(matEff), propDir(direction), mass(m), momentumEstimate(est) {}
0127 
0128     MaterialEffects materialEffects;
0129     PropagationDirection propDir;
0130     double mass;
0131     double momentumEstimate;
0132     bool useBeamSpot{false};
0133     bool hitsAreReverse{false};
0134     bool useRefittedState{false};
0135     bool constructTsosWithErrors{false};
0136     bool includeAPEs{false};
0137     bool allowZeroMaterial{false};
0138   };
0139 
0140   ~ReferenceTrajectoryBase() override {}
0141 
0142   bool isValid() { return theValidityFlag; }
0143 
0144   /** Returns the measurements in local coordinates.
0145    */
0146   const AlgebraicVector& measurements() const { return theMeasurements; }
0147 
0148   /** Returns the full covariance matrix of the measurements. ORCA-equivalent: covariance()
0149    */
0150   const AlgebraicSymMatrix& measurementErrors() const { return theMeasurementsCov; }
0151 
0152   /** Returns the local coordinates of the reference trajectory.
0153       ORCA-equivalent: referenceTrack()
0154    */
0155   const AlgebraicVector& trajectoryPositions() const { return theTrajectoryPositions; }
0156 
0157   /** Returns the covariance matrix of the reference trajectory.
0158    */
0159   const AlgebraicSymMatrix& trajectoryPositionErrors() const { return theTrajectoryPositionCov; }
0160 
0161   /** Returns the derivatives of the local coordinates of the reference
0162    *  trajectory (i.e. trajectoryPositions) w.r.t. the initial 'trajectory'-parameters.
0163    */
0164   const AlgebraicMatrix& derivatives() const { return theDerivatives; }
0165 
0166   /** Returns the transformation of tracjectory to curvilinear parameters
0167    */
0168   const AlgebraicMatrix& trajectoryToCurv() const { return theInnerTrajectoryToCurvilinear; }
0169   /** Returns the transformation of local to tracjectory parameters
0170    */
0171   const AlgebraicMatrix& localToTrajectory() const { return theInnerLocalToTrajectory; }
0172 
0173   /** Returns the GBL input
0174    */
0175   std::vector<std::pair<std::vector<gbl::GblPoint>, Eigen::MatrixXd> >& gblInput() { return theGblInput; }
0176 
0177   /** Returns the GBL external derivatives.
0178    */
0179   const Eigen::MatrixXd& gblExtDerivatives() const { return theGblExtDerivatives; }
0180 
0181   /** Returns the GBL external derivatives.
0182    */
0183   const Eigen::VectorXd& gblExtMeasurements() const { return theGblExtMeasurements; }
0184 
0185   /** Returns the GBL external derivatives.
0186    */
0187   const Eigen::VectorXd& gblExtPrecisions() const { return theGblExtPrecisions; }
0188 
0189   /** Returns the set of 'track'-parameters.
0190    */
0191   const AlgebraicVector& parameters() const { return theParameters; }
0192 
0193   /** Returns true if the covariance matrix of the 'track'-parameters is set.
0194    */
0195   inline bool parameterErrorsAvailable() const { return theParamCovFlag; }
0196 
0197   /** Set the covariance matrix of the 'track'-parameters.
0198    */
0199   inline void setParameterErrors(const AlgebraicSymMatrix& error) {
0200     theParameterCov = error;
0201     theParamCovFlag = true;
0202   }
0203 
0204   /** Returns the covariance matrix of the 'track'-parameters.
0205    */
0206   inline const AlgebraicSymMatrix& parameterErrors() const { return theParameterCov; }
0207 
0208   /** Returns the Tsos at the surfaces of the hits, parallel to recHits()
0209    */
0210   const std::vector<TrajectoryStateOnSurface>& trajectoryStates() const { return theTsosVec; }
0211 
0212   /** Returns the TransientTrackingRecHits (as ConstRecHitPointer's), 
0213       order might depend on concrete implementation of inheriting class
0214    */
0215   const TransientTrackingRecHit::ConstRecHitContainer& recHits() const { return theRecHits; }
0216 
0217   inline unsigned int numberOfHits() const { return theNumberOfHits; }
0218   inline unsigned int numberOfPar() const { return theNumberOfPars; }
0219   inline unsigned int numberOfVirtualMeas() const { return theNumberOfVirtualMeas; }
0220   inline unsigned int numberOfVirtualPar() const { return theNumberOfVirtualPars; }
0221   inline unsigned int numberOfHitMeas() const { return theNumberOfHits * nMeasPerHit; }
0222   inline int nominalField() const { return theNomField; }
0223 
0224   virtual ReferenceTrajectoryBase* clone() const = 0;
0225 
0226 protected:
0227   explicit ReferenceTrajectoryBase(unsigned int nPar,
0228                                    unsigned int nHits,
0229                                    unsigned int nVirtualPar,
0230                                    unsigned int nVirtualMeas);
0231 
0232   unsigned int numberOfUsedRecHits(const TransientTrackingRecHit::ConstRecHitContainer& recHits) const;
0233   bool useRecHit(const TransientTrackingRecHit::ConstRecHitPointer& hitPtr) const;
0234 
0235   bool theValidityFlag;
0236   bool theParamCovFlag;
0237 
0238   unsigned int theNumberOfHits;         // number of (measurements from) hits
0239   unsigned int theNumberOfPars;         // number of (track) parameters
0240   unsigned int theNumberOfVirtualMeas;  // number of virtual measurements
0241   unsigned int theNumberOfVirtualPars;  // number of parameters for virtual measurements
0242 
0243   std::vector<TrajectoryStateOnSurface> theTsosVec;
0244   TransientTrackingRecHit::ConstRecHitContainer theRecHits;
0245 
0246   AlgebraicVector theMeasurements;
0247   AlgebraicSymMatrix theMeasurementsCov;
0248 
0249   AlgebraicVector theTrajectoryPositions;
0250   AlgebraicSymMatrix theTrajectoryPositionCov;
0251 
0252   AlgebraicVector theParameters;
0253   AlgebraicSymMatrix theParameterCov;
0254 
0255   AlgebraicMatrix theDerivatives;
0256 
0257   // CHK for beamspot   transformation trajectory parameter to curvilinear at refTSos
0258   AlgebraicMatrix theInnerTrajectoryToCurvilinear;
0259   // CHK for TwoBodyD.  transformation local to trajectory parameter at refTsos
0260   AlgebraicMatrix theInnerLocalToTrajectory;
0261   // CHK GBL input:     list of (list of points on trajectory and transformation at inner (first) point)
0262   std::vector<std::pair<std::vector<gbl::GblPoint>, Eigen::MatrixXd> > theGblInput;
0263   int theNomField;
0264   // CHK GBL TBD:       virtual (mass) measurement
0265   Eigen::MatrixXd theGblExtDerivatives;
0266   Eigen::VectorXd theGblExtMeasurements;
0267   Eigen::VectorXd theGblExtPrecisions;
0268 
0269   static constexpr unsigned int nMeasPerHit{2};
0270 };
0271 
0272 #endif  // REFERENCE_TRAJECTORY_BASE_H