CMS 3D CMS Logo

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