CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_2_9_HLT1_bphpatch4/src/RecoTracker/NuclearSeedGenerator/src/SeedFromNuclearInteraction.cc

Go to the documentation of this file.
00001 #include "RecoTracker/NuclearSeedGenerator/interface/SeedFromNuclearInteraction.h"
00002 
00003 #include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
00004 #include "MagneticField/Engine/interface/MagneticField.h"
00005 
00006 #include "TrackingTools/TrajectoryState/interface/TrajectoryStateTransform.h"
00007 #include "TrackingTools/Records/interface/TrackingComponentsRecord.h"
00008 #include "TrackingTools/KalmanUpdators/interface/KFUpdator.h"
00009 
00010 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00011 
00012 
00013 SeedFromNuclearInteraction::SeedFromNuclearInteraction(const Propagator* prop, const TrackerGeometry* geom, const edm::ParameterSet& iConfig):
00014     ptMin(iConfig.getParameter<double>("ptMin")),
00015     thePropagator(prop), theTrackerGeom(geom)
00016     {
00017          isValid_=true;
00018          initialTSOS_ = boost::shared_ptr<TrajectoryStateOnSurface>(new TrajectoryStateOnSurface());
00019          updatedTSOS_ = boost::shared_ptr<TrajectoryStateOnSurface>(new TrajectoryStateOnSurface());
00020          freeTS_ = boost::shared_ptr<FreeTrajectoryState>(new FreeTrajectoryState());
00021          pTraj = boost::shared_ptr<PTrajectoryStateOnDet>(new PTrajectoryStateOnDet());
00022     }
00023 
00024 //----------------------------------------------------------------------
00025 void SeedFromNuclearInteraction::setMeasurements(const TSOS& inner_TSOS, ConstRecHitPointer ihit, ConstRecHitPointer ohit) {
00026 
00027        // delete pointer to TrackingRecHits
00028        theHits.clear();
00029 
00030        // get the inner and outer transient TrackingRecHits
00031        innerHit_ = ihit;
00032        outerHit_ = ohit;
00033 
00034        //theHits.push_back(  inner_TM.recHit() ); // put temporarily - TODO: remove this line
00035        theHits.push_back(  outerHit_  );
00036 
00037        initialTSOS_.reset( new TrajectoryStateOnSurface(inner_TSOS) );
00038 
00039        // calculate the initial FreeTrajectoryState.
00040        freeTS_.reset(stateWithError());
00041 
00042        // check transverse momentum
00043        if(freeTS_->momentum().perp() < ptMin) { isValid_ = false; }
00044        else {
00045           // convert freeTS_ into a persistent TSOS on the outer surface
00046           isValid_ = construct(); }
00047 }
00048 //----------------------------------------------------------------------
00049 void SeedFromNuclearInteraction::setMeasurements(TangentHelix& thePrimaryHelix, const TSOS& inner_TSOS, ConstRecHitPointer ihit, ConstRecHitPointer ohit) {
00050 
00051        // delete pointer to TrackingRecHits
00052        theHits.clear();
00053 
00054        // get the inner and outer transient TrackingRecHits
00055        innerHit_ = ihit;
00056        outerHit_ = ohit;
00057 
00058        GlobalPoint innerPos = theTrackerGeom->idToDet(innerHit_->geographicalId())->surface().toGlobal(innerHit_->localPosition());
00059        GlobalPoint outerPos = theTrackerGeom->idToDet(outerHit_->geographicalId())->surface().toGlobal(outerHit_->localPosition());
00060 
00061        TangentHelix helix(thePrimaryHelix, outerPos, innerPos);
00062 
00063        theHits.push_back( innerHit_ );
00064        theHits.push_back( outerHit_ );
00065 
00066        initialTSOS_.reset( new TrajectoryStateOnSurface(inner_TSOS) );
00067 
00068        // calculate the initial FreeTrajectoryState from the inner and outer TM assuming that the helix equation is already known.
00069        freeTS_.reset(stateWithError(helix));
00070 
00071        if(freeTS_->momentum().perp() < ptMin) { isValid_ = false; }
00072        else {
00073           // convert freeTS_ into a persistent TSOS on the outer surface
00074           isValid_ = construct(); }
00075 }
00076 //----------------------------------------------------------------------
00077 FreeTrajectoryState* SeedFromNuclearInteraction::stateWithError() const {
00078 
00079    // Calculation of the helix assuming that the secondary track has the same direction
00080    // than the primary track and pass through the inner and outer hits.
00081    GlobalVector direction = initialTSOS_->globalDirection();
00082    GlobalPoint inner = initialTSOS_->globalPosition();
00083    TangentHelix helix(direction, inner, outerHitPosition());
00084 
00085    return stateWithError(helix);
00086 }
00087 //----------------------------------------------------------------------
00088 FreeTrajectoryState* SeedFromNuclearInteraction::stateWithError(TangentHelix& helix) const {
00089 
00090 //   typedef TkRotation<float> Rotation;
00091 
00092    GlobalVector dirAtVtx = helix.directionAtVertex();
00093    const MagneticField& mag = initialTSOS_->globalParameters().magneticField();
00094 
00095    // Get the global parameters of the trajectory
00096    // we assume that the magnetic field at the vertex is equal to the magnetic field at the inner TM.
00097    GlobalTrajectoryParameters gtp(helix.vertexPoint(), dirAtVtx , helix.charge(mag.inTesla(helix.vertexPoint()).z())/helix.rho(), 0, &mag);
00098 
00099    // Error matrix in a frame where z is in the direction of the track at the vertex
00100    AlgebraicSymMatrix66 primaryError( initialTSOS_->cartesianError().matrix() );
00101    double p_max = initialTSOS_->globalParameters().momentum().mag();
00102    AlgebraicMatrix33 rot = this->rotationMatrix( dirAtVtx );
00103 
00104    AlgebraicMatrix66 globalRotation;
00105    globalRotation.Place_at(rot,0,0);
00106    globalRotation.Place_at(rot,3,3);
00107    AlgebraicSymMatrix66 primaryErrorInNewFrame = ROOT::Math::Similarity(globalRotation, primaryError);
00108 
00109    AlgebraicSymMatrix66 secondaryErrorInNewFrame = AlgebraicMatrixID();
00110    double p_perp_max = 2; // energy max of a secondary track emited perpendicularly to the 
00111                           // primary track is +/- 2 GeV
00112    secondaryErrorInNewFrame(0,0) = primaryErrorInNewFrame(0,0)+helix.vertexError()*p_perp_max/p_max;
00113    secondaryErrorInNewFrame(1,1) = primaryErrorInNewFrame(1,1)+helix.vertexError()*p_perp_max/p_max;
00114    secondaryErrorInNewFrame(2,2) = helix.vertexError() * helix.vertexError();
00115    secondaryErrorInNewFrame(3,3) = p_perp_max*p_perp_max;  
00116    secondaryErrorInNewFrame(4,4) = p_perp_max*p_perp_max; 
00117    secondaryErrorInNewFrame(5,5) = p_max*p_max;
00118 
00119    AlgebraicSymMatrix66 secondaryError = ROOT::Math::SimilarityT(globalRotation, secondaryErrorInNewFrame);
00120 
00121    return new FreeTrajectoryState( gtp, CartesianTrajectoryError(secondaryError) );
00122 }
00123 
00124 //----------------------------------------------------------------------
00125 bool SeedFromNuclearInteraction::construct() {
00126 
00127    // loop on all hits in theHits
00128    KFUpdator                 theUpdator;
00129 
00130    const TrackingRecHit* hit = 0;
00131 
00132    LogDebug("NuclearSeedGenerator") << "Seed ** initial state " << freeTS_->cartesianError().matrix();
00133 
00134    for ( unsigned int iHit = 0; iHit < theHits.size(); iHit++) {
00135      hit = theHits[iHit]->hit();
00136      TrajectoryStateOnSurface state = (iHit==0) ?
00137         thePropagator->propagate( *freeTS_, theTrackerGeom->idToDet(hit->geographicalId())->surface())
00138        : thePropagator->propagate( *updatedTSOS_, theTrackerGeom->idToDet(hit->geographicalId())->surface());
00139 
00140      if (!state.isValid()) return false;
00141 
00142      const TransientTrackingRecHit::ConstRecHitPointer& tth = theHits[iHit];
00143      updatedTSOS_.reset( new TrajectoryStateOnSurface(theUpdator.update(state, *tth)) );
00144 
00145    }
00146 
00147    TrajectoryStateTransform transformer;
00148 
00149    LogDebug("NuclearSeedGenerator") << "Seed ** updated state " << updatedTSOS_->cartesianError().matrix();
00150 
00151    pTraj = boost::shared_ptr<PTrajectoryStateOnDet>( transformer.persistentState(*updatedTSOS_, outerHitDetId().rawId()) );
00152    return true;
00153 }
00154 
00155 //----------------------------------------------------------------------
00156 edm::OwnVector<TrackingRecHit>  SeedFromNuclearInteraction::hits() const {
00157     recHitContainer      _hits;
00158     for( ConstRecHitContainer::const_iterator it = theHits.begin(); it!=theHits.end(); it++ ){
00159            _hits.push_back( it->get()->hit()->clone() );
00160     }
00161     return _hits;
00162 }
00163 //----------------------------------------------------------------------
00164 AlgebraicMatrix33 SeedFromNuclearInteraction::rotationMatrix(const GlobalVector& perp) const {
00165 
00166    AlgebraicMatrix33 result;
00167 
00168    // z axis coincides with perp
00169    GlobalVector zAxis = perp.unit();
00170 
00171    // x axis has no global Z component
00172    GlobalVector xAxis;
00173    if ( zAxis.x() != 0 || zAxis.y() != 0) {
00174      // precision is not an issue here, just protect against divizion by zero
00175      xAxis = GlobalVector( -zAxis.y(), zAxis.x(), 0).unit();
00176    }
00177    else { // perp coincides with global Z
00178      xAxis = GlobalVector( 1, 0, 0);
00179    }
00180 
00181    // y axis obtained by cross product
00182    GlobalVector yAxis( zAxis.cross( xAxis));
00183 
00184   result(0,0) = xAxis.x();
00185   result(0,1) = xAxis.y();
00186   result(0,2) = xAxis.z();
00187   result(1,0) = yAxis.x();
00188   result(1,1) = yAxis.y();
00189   result(1,2) = yAxis.z();
00190   result(2,0) = zAxis.x();
00191   result(2,1) = zAxis.y();
00192   result(2,2) = zAxis.z();
00193   return result;
00194 }
00195