CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_1_8_patch9/src/Alignment/ReferenceTrajectories/src/DualKalmanTrajectory.cc

Go to the documentation of this file.
00001 
00002 
00003 
00004 
00005 
00006 
00007 
00008 
00009  
00010 #include "Alignment/ReferenceTrajectories/interface/DualKalmanTrajectory.h"
00011 
00012 #include "TrackingTools/TrajectoryState/interface/TrajectoryStateOnSurface.h"
00013 #include "TrackingTools/PatternTools/interface/Trajectory.h"
00014 #include "TrackingTools/PatternTools/interface/TrajectoryMeasurement.h"
00015 
00016 #include "TrackingTools/TrackFitters/interface/TrajectoryStateCombiner.h"
00017 #include "TrackingTools/TransientTrackingRecHit/interface/HelpertRecHit2DLocalPos.h"
00018 
00019 #include "DataFormats/CLHEP/interface/AlgebraicObjects.h" 
00020 #include "DataFormats/TrajectoryState/interface/LocalTrajectoryParameters.h"
00021 #include "DataFormats/GeometrySurface/interface/LocalError.h"
00022 #include "DataFormats/GeometryVector/interface/LocalPoint.h"
00023 
00024 #include "Alignment/ReferenceTrajectories/interface/ReferenceTrajectory.h"
00025 
00026 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00027 #include "FWCore/Utilities/interface/Exception.h"
00028 
00029 typedef TransientTrackingRecHit::ConstRecHitContainer ConstRecHitContainer;
00030 typedef TransientTrackingRecHit::ConstRecHitPointer ConstRecHitPointer;
00031 
00032 //-----------------------------------------------------------------------------------------------
00033 DualKalmanTrajectory::DualKalmanTrajectory(const Trajectory::DataContainer &trajMeasurements,
00034                                            const TrajectoryStateOnSurface &referenceTsos,
00035                                            const std::vector<unsigned int> &forwardRecHitNums,
00036                                            const std::vector<unsigned int> &backwardRecHitNums,
00037                                            const MagneticField *magField,
00038                                            MaterialEffects materialEffects,
00039                                            PropagationDirection propDir,
00040                                            double mass, 
00041                                            bool useBeamSpot, const reco::BeamSpot &beamSpot,
00042                                            int residualMethod)
00043   : ReferenceTrajectoryBase(referenceTsos.localParameters().mixedFormatVector().kSize,
00044                             forwardRecHitNums.size() + backwardRecHitNums.size() - 1, 0, 0)
00045 {
00046   theValidityFlag = this->construct(trajMeasurements, referenceTsos,
00047                                     forwardRecHitNums, backwardRecHitNums,
00048                                     mass, materialEffects, propDir, magField,
00049                                     useBeamSpot, beamSpot, residualMethod);
00050 }
00051 
00052 
00053 
00054 //-----------------------------------------------------------------------------------------------
00055 DualKalmanTrajectory::DualKalmanTrajectory(unsigned int nPar, unsigned int nHits)
00056   : ReferenceTrajectoryBase(nPar, nHits, 0, 0)
00057 {}
00058 
00059 
00060 //-----------------------------------------------------------------------------------------------
00061 bool DualKalmanTrajectory::construct(const Trajectory::DataContainer &trajMeasurements,
00062                                      const TrajectoryStateOnSurface &refTsos, 
00063                                      const std::vector<unsigned int> &forwardRecHitNums,
00064                                      const std::vector<unsigned int> &backwardRecHitNums,
00065                                      double mass, MaterialEffects materialEffects,
00066                                      const PropagationDirection propDir,
00067                                      const MagneticField *magField, 
00068                                      bool useBeamSpot, const reco::BeamSpot &beamSpot,
00069                                      int residualMethod)
00070 {
00071 
00072   ReferenceTrajectoryPtr fwdTraj = this->construct(trajMeasurements, refTsos, forwardRecHitNums,
00073                                                    mass, materialEffects,
00074                                                    propDir, magField, useBeamSpot, beamSpot);
00075   
00076   ReferenceTrajectoryPtr bwdTraj = this->construct(trajMeasurements, refTsos, backwardRecHitNums,
00077                                                    mass, materialEffects,
00078                                                    this->oppositeDirection(propDir), magField,
00079                                                    useBeamSpot, beamSpot);
00080 
00081   //   edm::LogError("Alignment") << "@SUB=DualKalmanTrajectory::construct"
00082   //                         << "valid fwd/bwd " << fwdTraj->isValid() 
00083   //                         << "/" << bwdTraj->isValid();
00084   
00085   if (!fwdTraj->isValid() || !bwdTraj->isValid()) {
00086     return false;
00087   }
00088 
00089   //
00090   // Combine both reference trajactories to a dual reference trajectory
00091   //
00092 
00093 //   const std::vector<TrajectoryStateOnSurface>& fwdTsosVec = fwdTraj->trajectoryStates();
00094 //   const std::vector<TrajectoryStateOnSurface>& bwdTsosVec = bwdTraj->trajectoryStates();
00095 //   theTsosVec.insert( theTsosVec.end(), fwdTsosVec.begin(), fwdTsosVec.end() );
00096 //   theTsosVec.insert( theTsosVec.end(), ++bwdTsosVec.begin(), bwdTsosVec.end() );
00097 
00098   // Take hits as they come from the Kalman fit.
00099   const ConstRecHitContainer &fwdRecHits = fwdTraj->recHits(); 
00100   const ConstRecHitContainer &bwdRecHits = bwdTraj->recHits();
00101   theRecHits.insert(theRecHits.end(), fwdRecHits.begin(), fwdRecHits.end());
00102   theRecHits.insert(theRecHits.end(), ++bwdRecHits.begin(), bwdRecHits.end());
00103 
00104   theParameters = this->extractParameters(refTsos);
00105 
00106   unsigned int nParam   = theParameters.num_row();
00107   unsigned int nFwdMeas = nMeasPerHit * fwdTraj->numberOfHits();
00108   unsigned int nBwdMeas = nMeasPerHit * bwdTraj->numberOfHits();
00109 
00110 //   theMeasurements.sub( 1, fwdTraj->measurements() );
00111 //   theMeasurements.sub( nFwdMeas+1, bwdTraj->measurements().sub( nMeasPerHit+1, nBwdMeas ) );
00112 
00113 //   theMeasurementsCov.sub( 1, fwdTraj->measurementErrors() );
00114 //   theMeasurementsCov.sub( nFwdMeas+1, bwdTraj->measurementErrors().sub( nMeasPerHit+1, nBwdMeas ) );
00115 
00116 //   theTrajectoryPositions.sub( 1, fwdTraj->trajectoryPositions() );
00117 //   theTrajectoryPositions.sub( nFwdMeas+1, bwdTraj->trajectoryPositions().sub( nMeasPerHit+1, nBwdMeas ) );
00118 
00119 //   theTrajectoryPositionCov.sub( 1, fwdTraj->trajectoryPositionErrors() );
00120 //   theTrajectoryPositionCov.sub( nFwdMeas+1, bwdTraj->trajectoryPositionErrors().sub( nMeasPerHit+1, nBwdMeas ) );
00121 
00122   theDerivatives.sub(1, 1, fwdTraj->derivatives());
00123   theDerivatives.sub(nFwdMeas+1, 1,
00124                      bwdTraj->derivatives().sub(nMeasPerHit+1, nBwdMeas, 1, nParam));
00125 
00126   // FIXME: next lines stolen from ReferenceTrajectory, no idea whether OK or not...
00127   if (refTsos.hasError()) {
00128     AlgebraicSymMatrix parameterCov = asHepMatrix<5>(refTsos.localError().matrix());
00129     theTrajectoryPositionCov = parameterCov.similarity(theDerivatives);
00130   } else {
00131     theTrajectoryPositionCov = AlgebraicSymMatrix(theDerivatives.num_row(), 1);
00132   }
00133 
00134   // Fill Kalman part, first for forward, then for backward part.
00135   return this->fillKalmanPart(trajMeasurements, forwardRecHitNums, true, 0, residualMethod) 
00136     && this->fillKalmanPart(trajMeasurements, backwardRecHitNums, false,
00137                             forwardRecHitNums.size(), residualMethod);
00138 }
00139 
00140 
00141 //-----------------------------------------------------------------------------------------------
00142 ReferenceTrajectory*
00143 DualKalmanTrajectory::construct(const Trajectory::DataContainer &trajMeasurements,
00144                                 const TrajectoryStateOnSurface &referenceTsos, 
00145                                 const std::vector<unsigned int> &recHitNums,
00146                                 double mass, MaterialEffects materialEffects,
00147                                 const PropagationDirection propDir,
00148                                 const MagneticField *magField,
00149                                 bool useBeamSpot, const reco::BeamSpot &beamSpot) const
00150 {
00151 
00152   ConstRecHitContainer recHits;
00153   recHits.reserve(recHitNums.size());
00154   for (unsigned int i = 0; i < recHitNums.size(); ++i) {
00155     recHits.push_back(trajMeasurements[recHitNums[i]].recHit());
00156   }
00157 
00158   return new ReferenceTrajectory(referenceTsos, recHits, false, // hits are already ordered
00159                                  magField, materialEffects, propDir, mass,
00160                                  useBeamSpot, beamSpot);
00161 }
00162 
00163 //-----------------------------------------------------------------------------------------------
00164 bool DualKalmanTrajectory::fillKalmanPart(const Trajectory::DataContainer &trajMeasurements,
00165                                           const std::vector<unsigned int> &recHitNums,
00166                                           bool startFirst, unsigned int iNextHit,
00167                                           int residualMethod)
00168 {
00169   // startFirst==false: skip first hit of recHitNums
00170   // iNextHit: first hit number to fill data members with
00171   //
00172   // Three approaches, choosen by 'residualMethod':
00173   // 0: Just take hit error and the updated state as millepede needs it to redo the last fit step.
00174   // 1: Use the unbiased residuals as for residual monitoring.
00175   // 2: Use the _updated_ state and calculate the sigma that is part of
00176   //    the pull as sqrt(sigma_hit^2 - sigma_tsos^2).
00177   //    This should lead to the pull as defined on p. 236 of Blobel's book.
00178   //    Not sure whether this is 100% applicable/correct here, at least it might lead to 
00179   //    numerical instabilities, cf.:
00180   //    https://hypernews.cern.ch/HyperNews/CMS/get/recoTracking/517/1/1/1/1/1.html
00181 
00182   for (unsigned int iMeas = (startFirst ? 0 : 1); iMeas < recHitNums.size(); ++iMeas, ++iNextHit) {
00183     const TrajectoryMeasurement &trajMeasurement = trajMeasurements[recHitNums[iMeas]];
00184 
00185     TrajectoryStateOnSurface tsos; // depends on method
00186     switch (residualMethod) {
00187     case 0:
00188       tsos = this->fillMeasurementAndError2(trajMeasurement.recHit(), iNextHit, trajMeasurement,
00189                                             false); // plain hit error, not pull
00190       break;
00191     case 1:
00192       tsos = this->fillMeasurementAndError1(trajMeasurement.recHit(), iNextHit, trajMeasurement);
00193       break;
00194     case 2:
00195       tsos = this->fillMeasurementAndError2(trajMeasurement.recHit(), iNextHit, trajMeasurement,
00196                                             true); // error of pull
00197       break;
00198     default:
00199       throw cms::Exception("BadConfig")
00200         << "[DualKalmanTrajectory::fillKalmanPart] expect residualMethod == 0, 1 or 2, not "
00201         << residualMethod << ".";
00202     }
00203 
00204     if (!tsos.isValid()) return false;
00205     theTsosVec.push_back(tsos);
00206     
00207     this->fillTrajectoryPositions(trajMeasurement.recHit()->projectionMatrix(),
00208                                   tsos, iNextHit);
00209   }
00210 
00211   return true;
00212 } 
00213 
00214 //-----------------------------------------------------------------------------------------------
00215 TrajectoryStateOnSurface 
00216 DualKalmanTrajectory::fillMeasurementAndError1(const ConstRecHitPointer &hitPtr,
00217                                                unsigned int iHit,
00218                                                const TrajectoryMeasurement &trajMeasurement)
00219 {
00220   // Get the measurements and their errors.
00221   // We have to add error from hit and that tsos that is combination of fwd and bwd state.
00222 
00223   TrajectoryStateCombiner tsosComb;
00224   const TrajectoryStateOnSurface tsos = tsosComb(trajMeasurement.forwardPredictedState(),
00225                                                  trajMeasurement.backwardPredictedState());
00226   if (tsos.isValid()) {
00227     // No update of hit with tsos: it comes already from fwd+bwd tsos combination.
00228     // See also https://hypernews.cern.ch/HyperNews/CMS/get/recoTracking/517/1.html .
00229     // ConstRecHitPointer newHitPtr(hitPtr->canImproveWithTrack() ? hitPtr->clone(tsos) : hitPtr);
00230     
00231     const LocalPoint localMeasurement(hitPtr->localPosition());
00232     const LocalError hitErr(this->hitErrorWithAPE(hitPtr));
00233     // tsos prediction includes APE of other hits:
00234     const LocalError tsosErr(tsos.localError().positionError());
00235     const LocalError localMeasurementCov(hitErr.xx() + tsosErr.xx(),
00236                                          hitErr.xy() + tsosErr.xy(),
00237                                          hitErr.yy() + tsosErr.yy());
00238     
00239     theMeasurements[nMeasPerHit*iHit]   = localMeasurement.x();
00240     theMeasurements[nMeasPerHit*iHit+1] = localMeasurement.y();
00241     theMeasurementsCov[nMeasPerHit*iHit]  [nMeasPerHit*iHit]   = localMeasurementCov.xx();
00242     theMeasurementsCov[nMeasPerHit*iHit]  [nMeasPerHit*iHit+1] = localMeasurementCov.xy();
00243     theMeasurementsCov[nMeasPerHit*iHit+1][nMeasPerHit*iHit+1] = localMeasurementCov.yy();
00244   }
00245 
00246   return tsos;
00247 }
00248 
00249 //-----------------------------------------------------------------------------------------------
00250 TrajectoryStateOnSurface
00251 DualKalmanTrajectory::fillMeasurementAndError2(const ConstRecHitPointer &hitPtr,
00252                                                unsigned int iHit,
00253                                                const TrajectoryMeasurement &trajMeasurement,
00254                                                bool doPull)
00255 {
00256   // Get the measurements and their errors.
00257   // We have to subtract error of updated tsos from hit.
00258 
00259   const TrajectoryStateOnSurface tsos = trajMeasurement.updatedState();
00260   if (tsos.isValid()) {
00261     // No further update of hit: 
00262     // - The Kalman fit used hitPtr as it comes here (besides APE, see below).
00263     // - If the hit errors improve, we might get (rare) problems of negative diagonal elements, see below.
00264     // ConstRecHitPointer newHitPtr(hitPtr->canImproveWithTrack() ? hitPtr->clone(tsos) : hitPtr);
00265 
00266     const LocalPoint localMeasurement(hitPtr->localPosition());
00267     // Add APE manually to avoid that the hit error might be smaller than tsos error:
00268     // - hit local errors are always without APE,
00269     // - the tsos errors include APE since they come from track fit.
00270     const LocalError hitErr(this->hitErrorWithAPE(hitPtr));
00271 
00272     LocalError localMeasurementCov(hitErr.xx(), hitErr.xy(), hitErr.yy());
00273     if (doPull) {
00274       const LocalError tsosErr(tsos.localError().positionError());
00275       // Should not be possible to become negative if all is correct - see above.
00276       if (hitErr.xx() < tsosErr.xx() || hitErr.yy() < tsosErr.yy()) {
00277         edm::LogError("Alignment") << "@SUB=DualKalmanTrajectory::fillMeasurementAndError2"
00278                                    << "not OK in subdet " << hitPtr->geographicalId().subdetId()
00279                                    << "\ns_x " << sqrt(hitErr.xx()) << " " << sqrt(tsosErr.xx())
00280                                    << "\ns_xy " << hitErr.xy() << " " << tsosErr.xy()
00281                                    << "\ns_y " << sqrt(hitErr.yy()) << " " << sqrt(tsosErr.yy());
00282         return TrajectoryStateOnSurface(); // is invalid state
00283       }
00284       // cf. Blobel/Lohrmann, p. 236:
00285       // But numerical stability? Should we return false if difference very small compared to values?
00286       localMeasurementCov = LocalError(hitErr.xx() - tsosErr.xx(), // tsos puts correlation in,   
00287                                        hitErr.xy() - tsosErr.xy(), // even for 1D strip!            
00288                                        hitErr.yy() - tsosErr.yy());
00289     }
00290 
00291     theMeasurements[nMeasPerHit*iHit]   = localMeasurement.x();
00292     theMeasurements[nMeasPerHit*iHit+1] = localMeasurement.y();
00293     theMeasurementsCov[nMeasPerHit*iHit]  [nMeasPerHit*iHit]   = localMeasurementCov.xx();
00294     theMeasurementsCov[nMeasPerHit*iHit]  [nMeasPerHit*iHit+1] = localMeasurementCov.xy();
00295     theMeasurementsCov[nMeasPerHit*iHit+1][nMeasPerHit*iHit+1] = localMeasurementCov.yy();
00296   }
00297 
00298   return tsos;
00299 }
00300 
00301 //-----------------------------------------------------------------------------------------------
00302 void DualKalmanTrajectory::fillTrajectoryPositions(const AlgebraicMatrix &projection, 
00303                                                    const TrajectoryStateOnSurface &tsos, 
00304                                                    unsigned int iHit)
00305 {
00306   // get the local coordinates of the reference trajectory
00307   // (~copied from ReferenceTrajectory::fillTrajectoryPositions)
00308   AlgebraicVector mixedLocalParams = asHepVector<5>(tsos.localParameters().mixedFormatVector());
00309   const AlgebraicVector localPosition(projection * mixedLocalParams);
00310 
00311   theTrajectoryPositions[nMeasPerHit*iHit]   = localPosition[0];
00312   theTrajectoryPositions[nMeasPerHit*iHit+1] = localPosition[1];
00313 }
00314 
00315 //-----------------------------------------------------------------------------------------------
00316 LocalError DualKalmanTrajectory::hitErrorWithAPE(const ConstRecHitPointer &hitPtr) const
00317 {
00318   if (hitPtr->det() && hitPtr->det()->alignmentPositionError()) {
00319     HelpertRecHit2DLocalPos help;
00320     const AlgebraicSymMatrix errMat(help.parError(hitPtr->localPositionError(), *(hitPtr->det())));
00321     return LocalError(errMat[0][0], errMat[0][1], errMat[1][1]);
00322   } else {
00323     return hitPtr->localPositionError();
00324   }
00325 }
00326 
00327 //-----------------------------------------------------------------------------------------------
00328 AlgebraicVector
00329 DualKalmanTrajectory::extractParameters(const TrajectoryStateOnSurface &referenceTsos) const
00330 {
00331   return asHepVector<5>( referenceTsos.localParameters().mixedFormatVector() );
00332 }