CMS 3D CMS Logo

/data/doxygen/doxygen-1.7.3/gen/CMSSW_4_2_8/src/Alignment/ReferenceTrajectories/src/ReferenceTrajectory.cc

Go to the documentation of this file.
00001 //  Author     : Gero Flucke (based on code by Edmund Widl replacing ORCA's TkReferenceTrack)
00002 //  date       : 2006/09/17
00003 //  last update: $Date: 2010/11/26 07:54:53 $
00004 //  by         : $Author: mussgill $
00005 
00006 #include <memory>
00007 
00008 #include "Alignment/ReferenceTrajectories/interface/ReferenceTrajectory.h"
00009 
00010 #include "DataFormats/GeometrySurface/interface/Surface.h" 
00011 #include "DataFormats/GeometrySurface/interface/Plane.h"
00012 #include "DataFormats/GeometrySurface/interface/OpenBounds.h"
00013 
00014 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00015 
00016 #include "DataFormats/CLHEP/interface/AlgebraicObjects.h" 
00017 #include "DataFormats/GeometrySurface/interface/LocalError.h"
00018 #include "DataFormats/GeometryVector/interface/LocalPoint.h"
00019 #include "Geometry/CommonDetUnit/interface/GeomDet.h"
00020 
00021 #include "DataFormats/TrajectoryState/interface/LocalTrajectoryParameters.h"
00022 
00023 #include "TrackingTools/AnalyticalJacobians/interface/AnalyticalCurvilinearJacobian.h"
00024 #include "TrackingTools/AnalyticalJacobians/interface/JacobianLocalToCurvilinear.h"
00025 #include "TrackingTools/AnalyticalJacobians/interface/JacobianCurvilinearToLocal.h"
00026 
00027 #include "TrackingTools/GeomPropagators/interface/AnalyticalPropagator.h"
00028 #include "TrackPropagation/RungeKutta/interface/RKTestPropagator.h"
00029 
00030 #include "TrackingTools/TrajectoryState/interface/TrajectoryStateOnSurface.h"
00031 #include "TrackingTools/TrajectoryParametrization/interface/GlobalTrajectoryParameters.h"
00032 
00033 #include "TrackingTools/MaterialEffects/interface/MultipleScatteringUpdator.h"
00034 #include "TrackingTools/MaterialEffects/interface/EnergyLossUpdator.h"
00035 #include "TrackingTools/MaterialEffects/interface/CombinedMaterialEffectsUpdator.h"
00036 #include <TrackingTools/PatternTools/interface/TSCPBuilderNoMaterial.h>
00037 #include <TrackingTools/PatternTools/interface/TSCBLBuilderNoMaterial.h>
00038 #include "TrackingTools/TrajectoryState/interface/TrajectoryStateClosestToPoint.h"
00039 #include "TrackingTools/TrajectoryState/interface/TrajectoryStateClosestToBeamLine.h"
00040 
00041 #include "MagneticField/Engine/interface/MagneticField.h"
00042 
00043 #include "BeamSpotTransientTrackingRecHit.h"
00044 #include "BeamSpotGeomDet.h"
00045 
00046 //__________________________________________________________________________________
00047 
00048 ReferenceTrajectory::ReferenceTrajectory(const TrajectoryStateOnSurface &refTsos,
00049                                          const TransientTrackingRecHit::ConstRecHitContainer
00050                                          &recHits, bool hitsAreReverse,
00051                                          const MagneticField *magField, 
00052                                          MaterialEffects materialEffects,
00053                                          PropagationDirection propDir,
00054                                          double mass,
00055                                          bool useBeamSpot, const reco::BeamSpot &beamSpot) 
00056  : ReferenceTrajectoryBase( 
00057    (materialEffects >= brokenLinesCoarse) ? 1 : refTsos.localParameters().mixedFormatVector().kSize, 
00058    (useBeamSpot == true) ? recHits.size()+1 : recHits.size(),
00059    (materialEffects >= brokenLinesCoarse) ? 
00060        2*((useBeamSpot == true) ? recHits.size()+1 : recHits.size())   :
00061    ( (materialEffects == breakPoints) ? 2*((useBeamSpot == true) ? recHits.size()+1 : recHits.size())-2 : 0) , 
00062    (materialEffects >= brokenLinesCoarse) ? 
00063        2*((useBeamSpot == true) ? recHits.size()+1 : recHits.size())-4 : 
00064    ( (materialEffects == breakPoints) ? 2*((useBeamSpot == true) ? recHits.size()+1 : recHits.size())-2 : 0) )
00065 {
00066   // no check against magField == 0  
00067   theParameters = asHepVector<5>( refTsos.localParameters().mixedFormatVector() );
00068   
00069   if (hitsAreReverse) {
00070     TransientTrackingRecHit::ConstRecHitContainer fwdRecHits;
00071     fwdRecHits.reserve(recHits.size());
00072     for (TransientTrackingRecHit::ConstRecHitContainer::const_reverse_iterator it=recHits.rbegin();
00073          it != recHits.rend(); ++it) {
00074       fwdRecHits.push_back(*it);
00075     }
00076     theValidityFlag = this->construct(refTsos, fwdRecHits, mass, materialEffects,
00077                                       propDir, magField,
00078                                       useBeamSpot, beamSpot);
00079   } else {
00080     theValidityFlag = this->construct(refTsos, recHits, mass, materialEffects,
00081                                       propDir, magField,
00082                                       useBeamSpot, beamSpot);
00083   }
00084 }
00085 
00086 
00087 //__________________________________________________________________________________
00088 
00089 ReferenceTrajectory::ReferenceTrajectory( unsigned int nPar, unsigned int nHits,
00090                                           MaterialEffects materialEffects)
00091  : ReferenceTrajectoryBase( 
00092    (materialEffects >= brokenLinesCoarse) ? 1 : nPar, 
00093    nHits, 
00094    (materialEffects >= brokenLinesCoarse) ? 2*nHits   : ( (materialEffects == breakPoints) ? 2*nHits-2 : 0 ), 
00095    (materialEffects >= brokenLinesCoarse) ? 2*nHits-4 : ( (materialEffects == breakPoints) ? 2*nHits-2 : 0 ) )
00096 {}
00097 
00098 
00099 //__________________________________________________________________________________
00100 
00101 bool ReferenceTrajectory::construct(const TrajectoryStateOnSurface &refTsos, 
00102                                     const TransientTrackingRecHit::ConstRecHitContainer &recHits,
00103                                     double mass, MaterialEffects materialEffects,
00104                                     const PropagationDirection propDir,
00105                                     const MagneticField *magField,
00106                                     bool useBeamSpot,
00107                                     const reco::BeamSpot &beamSpot)
00108 {   
00109   TrajectoryStateOnSurface theRefTsos = refTsos;
00110 
00111   const SurfaceSide surfaceSide = this->surfaceSide(propDir);
00112   // auto_ptr to avoid memory leaks in case of not reaching delete at end of method:
00113   std::auto_ptr<MaterialEffectsUpdator> aMaterialEffectsUpdator
00114     (this->createUpdator(materialEffects, mass));
00115   if (!aMaterialEffectsUpdator.get()) return false; // empty auto_ptr
00116 
00117   AlgebraicMatrix                 fullJacobian(theParameters.num_row(), theParameters.num_row());
00118   std::vector<AlgebraicMatrix>    allJacobians; 
00119   allJacobians.reserve(theNumberOfHits);
00120 
00121   TransientTrackingRecHit::ConstRecHitPointer  previousHitPtr;
00122   TrajectoryStateOnSurface                     previousTsos;
00123   AlgebraicSymMatrix              previousChangeInCurvature(theParameters.num_row(), 1);
00124   std::vector<AlgebraicSymMatrix> allCurvatureChanges; 
00125   allCurvatureChanges.reserve(theNumberOfHits);
00126   
00127   const LocalTrajectoryError zeroErrors(0., 0., 0., 0., 0.);
00128 
00129   std::vector<AlgebraicMatrix> allProjections;
00130   allProjections.reserve(theNumberOfHits);
00131   std::vector<AlgebraicSymMatrix> allDeltaParameterCovs;
00132   allDeltaParameterCovs.reserve(theNumberOfHits);
00133 
00134   // CHK
00135   std::vector<AlgebraicMatrix> allLocalToCurv;
00136   allLocalToCurv.reserve(theNumberOfHits); 
00137   std::vector<double> allSteps;
00138   allSteps.reserve(theNumberOfHits); 
00139   std::vector<AlgebraicMatrix>    allCurvlinJacobians; 
00140   allCurvlinJacobians.reserve(theNumberOfHits);
00141   
00142   AlgebraicMatrix firstCurvlinJacobian(5, 5, 1);
00143   
00144   unsigned int iRow = 0;
00145 
00146   // local storage vector of all rechits (including rechit for beam spot in case it is used)
00147   TransientTrackingRecHit::ConstRecHitContainer allRecHits;
00148 
00149   if (useBeamSpot && propDir==alongMomentum) {
00150     
00151     GlobalPoint bs(beamSpot.x0(), beamSpot.y0(), beamSpot.z0());
00152     
00153     TrajectoryStateClosestToBeamLine tsctbl(TSCBLBuilderNoMaterial()(*(refTsos.freeState()), beamSpot));
00154     FreeTrajectoryState pcaFts = tsctbl.trackStateAtPCA();
00155     GlobalVector bd(beamSpot.dxdz(), beamSpot.dydz(), 1.0);
00156     
00157     //propagation FIXME: Should use same propagator everywhere...
00158     AnalyticalPropagator propagator(magField);
00159     std::pair< TrajectoryStateOnSurface, double > tsosWithPath =
00160       propagator.propagateWithPath(pcaFts, refTsos.surface());
00161     
00162     if (!tsosWithPath.first.isValid()) return false;
00163     
00164     GlobalVector momDir(pcaFts.momentum());
00165     GlobalVector perpDir(bd.cross(momDir));
00166     BoundPlane::RotationType rotation(perpDir, bd);
00167     
00168     BeamSpotGeomDet * bsGeom = new BeamSpotGeomDet(BoundPlane::build(bs, rotation, OpenBounds()));
00169 
00170     // There is also a constructor taking the magentic field. Use this one instead?
00171     theRefTsos = TrajectoryStateOnSurface(pcaFts, bsGeom->surface());
00172     
00173     TransientTrackingRecHit::ConstRecHitPointer bsRecHit = 
00174       new BeamSpotTransientTrackingRecHit(beamSpot,
00175                                           bsGeom,
00176                                           theRefTsos.freeState()->momentum().phi());
00177     allRecHits.push_back(bsRecHit);
00178 
00179   }
00180   
00181   // copy all rechits to the local storage vector
00182   TransientTrackingRecHit::ConstRecHitContainer::const_iterator itRecHit;
00183   for ( itRecHit = recHits.begin(); itRecHit != recHits.end(); ++itRecHit ) { 
00184     const TransientTrackingRecHit::ConstRecHitPointer &hitPtr = *itRecHit;
00185     allRecHits.push_back(hitPtr);
00186   }
00187 
00188   for ( itRecHit = allRecHits.begin(); itRecHit != allRecHits.end(); ++itRecHit ) { 
00189     
00190     const TransientTrackingRecHit::ConstRecHitPointer &hitPtr = *itRecHit;
00191     theRecHits.push_back(hitPtr);
00192 
00193     if (0 == iRow) {
00194 
00195       // compute the derivatives of the reference-track's parameters w.r.t. the initial ones
00196       // derivative of the initial reference-track parameters w.r.t. themselves is of course the identity 
00197       fullJacobian = AlgebraicMatrix(theParameters.num_row(), theParameters.num_row(), 1);
00198       allJacobians.push_back(fullJacobian);
00199       theTsosVec.push_back(theRefTsos);
00200       const JacobianLocalToCurvilinear startTrafo(hitPtr->det()->surface(), theRefTsos.localParameters(), *magField);
00201       const AlgebraicMatrix localToCurvilinear =  asHepMatrix<5>(startTrafo.jacobian());
00202       if (materialEffects <= breakPoints) {
00203          theInnerTrajectoryToCurvilinear = asHepMatrix<5>(startTrafo.jacobian());
00204          theInnerLocalToTrajectory = AlgebraicMatrix(5, 5, 1);
00205       }  
00206       allLocalToCurv.push_back(localToCurvilinear);
00207       allSteps.push_back(0.);
00208       allCurvlinJacobians.push_back(firstCurvlinJacobian);
00209 
00210     } else {
00211 
00212       AlgebraicMatrix nextJacobian;
00213       AlgebraicMatrix nextCurvlinJacobian;
00214       double nextStep = 0.;
00215       TrajectoryStateOnSurface nextTsos;
00216 
00217       if (!this->propagate(previousHitPtr->det()->surface(), previousTsos,
00218                            hitPtr->det()->surface(), nextTsos,
00219                            nextJacobian, nextCurvlinJacobian, nextStep, propDir, magField)) {
00220         return false; // stop if problem...// no delete aMaterialEffectsUpdator needed
00221       }
00222       
00223       allJacobians.push_back(nextJacobian);
00224       fullJacobian = nextJacobian * previousChangeInCurvature * fullJacobian;
00225       theTsosVec.push_back(nextTsos);
00226       
00227       const JacobianLocalToCurvilinear startTrafo(hitPtr->det()->surface(), nextTsos.localParameters(), *magField);
00228       const AlgebraicMatrix localToCurvilinear =  asHepMatrix<5>(startTrafo.jacobian());
00229       allLocalToCurv.push_back(localToCurvilinear);
00230       if (nextStep == 0.) { 
00231         edm::LogError("Alignment") << "@SUB=ReferenceTrajectory::construct"
00232                                    << "step 0. from id " << previousHitPtr->geographicalId()
00233                                    << " to " << hitPtr->det()->geographicalId() << ".";
00234         // brokenLinesFine will not work, brokenLinesCoarse combines close by layers
00235         if (materialEffects == brokenLinesFine) {
00236           edm::LogError("Alignment") << "@SUB=ReferenceTrajectory::construct" << "Skip track.";
00237           return false;
00238         }
00239       }
00240       allSteps.push_back(nextStep);
00241       allCurvlinJacobians.push_back(nextCurvlinJacobian);
00242 
00243     }
00244 
00245     // take material effects into account. since trajectory-state is constructed with errors equal zero,
00246     // the updated state contains only the uncertainties due to interactions in the current layer.
00247     const TrajectoryStateOnSurface tmpTsos(theTsosVec.back().localParameters(), zeroErrors,
00248                                            theTsosVec.back().surface(), magField, surfaceSide);
00249     const TrajectoryStateOnSurface updatedTsos = aMaterialEffectsUpdator->updateState(tmpTsos, propDir);
00250 
00251     if ( !updatedTsos.isValid() ) return false;// no delete aMaterialEffectsUpdator needed
00252     
00253     if ( theTsosVec.back().localParameters().charge() )
00254     {
00255       previousChangeInCurvature[0][0] = updatedTsos.localParameters().signedInverseMomentum() 
00256         / theTsosVec.back().localParameters().signedInverseMomentum();
00257     }
00258     
00259     // get multiple-scattering covariance-matrix
00260     allDeltaParameterCovs.push_back( asHepMatrix<5>(updatedTsos.localError().matrix()) );
00261     allCurvatureChanges.push_back(previousChangeInCurvature);
00262     
00263     // projection-matrix tsos-parameters -> measurement-coordinates
00264     if ( useRecHit( hitPtr ) ) { allProjections.push_back( hitPtr->projectionMatrix() ); }
00265     else                       { allProjections.push_back( AlgebraicMatrix( 2, 5, 0) ); } // get size from ???
00266     
00267     // set start-parameters for next propagation. trajectory-state without error
00268     //  - no error propagation needed here.
00269     previousHitPtr = hitPtr;
00270     previousTsos   = TrajectoryStateOnSurface(updatedTsos.globalParameters(),
00271                                               updatedTsos.surface(), surfaceSide);
00272     
00273     if (materialEffects < brokenLinesCoarse) {
00274       this->fillDerivatives(allProjections.back(), fullJacobian, iRow);
00275     }
00276 
00277     AlgebraicVector mixedLocalParams = asHepVector<5>(theTsosVec.back().localParameters().mixedFormatVector());
00278     this->fillTrajectoryPositions(allProjections.back(), mixedLocalParams, iRow);
00279     if ( useRecHit( hitPtr ) ) this->fillMeasurementAndError(hitPtr, iRow, updatedTsos);
00280 
00281     iRow += nMeasPerHit;
00282   } // end of loop on hits
00283 
00284   bool msOK = true;
00285   switch (materialEffects) {
00286   case none:
00287     break;
00288   case multipleScattering:
00289   case energyLoss:
00290   case combined:
00291     msOK = this->addMaterialEffectsCov(allJacobians, allProjections, allCurvatureChanges,
00292                                        allDeltaParameterCovs);
00293     break;
00294   case breakPoints:
00295     msOK = this->addMaterialEffectsBp(allJacobians, allProjections, allCurvatureChanges,
00296                                       allDeltaParameterCovs, allLocalToCurv);
00297     break;
00298   case brokenLinesCoarse:
00299     msOK = this->addMaterialEffectsBrl(allProjections, allDeltaParameterCovs, allLocalToCurv,
00300                                        allSteps, refTsos.globalParameters());
00301     break;
00302   case brokenLinesFine:
00303     msOK = this->addMaterialEffectsBrl(allCurvlinJacobians, allProjections, allCurvatureChanges,
00304                                        allDeltaParameterCovs, allLocalToCurv, refTsos.globalParameters());
00305   }
00306   if (!msOK) return false;
00307  
00308   if (refTsos.hasError()) {
00309     AlgebraicSymMatrix parameterCov = asHepMatrix<5>(refTsos.localError().matrix());
00310     AlgebraicMatrix  parDeriv;
00311     if (theNumberOfMsPars>0) { 
00312       parDeriv = theDerivatives.sub( 1, nMeasPerHit*allJacobians.size(), 1, theParameters.num_row() ); 
00313     } else {
00314       parDeriv = theDerivatives; 
00315     }
00316     theTrajectoryPositionCov = parameterCov.similarity(parDeriv);
00317   } else {
00318     theTrajectoryPositionCov = AlgebraicSymMatrix(theDerivatives.num_row(), 1);
00319   }
00320 
00321   return true;
00322 }
00323 
00324 //__________________________________________________________________________________
00325 
00326 MaterialEffectsUpdator*
00327 ReferenceTrajectory::createUpdator(MaterialEffects materialEffects, double mass) const
00328 {
00329   switch (materialEffects) {
00330     // MultipleScatteringUpdator doesn't change the trajectory-state
00331     // during update and can therefore be used if material effects should be ignored:
00332   case none:
00333   case multipleScattering: 
00334     return new MultipleScatteringUpdator(mass);
00335   case energyLoss:
00336     return new EnergyLossUpdator(mass);
00337   case combined:
00338     return new CombinedMaterialEffectsUpdator(mass);
00339   case breakPoints:
00340     return new CombinedMaterialEffectsUpdator(mass);
00341   case brokenLinesCoarse:
00342   case brokenLinesFine:
00343     return new CombinedMaterialEffectsUpdator(mass);
00344 }
00345 
00346   return 0;
00347 }
00348 
00349 //__________________________________________________________________________________
00350 
00351 bool ReferenceTrajectory::propagate(const BoundPlane &previousSurface, const TrajectoryStateOnSurface &previousTsos,
00352                                     const BoundPlane &newSurface, TrajectoryStateOnSurface &newTsos, AlgebraicMatrix &newJacobian, 
00353                                     AlgebraicMatrix &newCurvlinJacobian, double &nextStep,
00354                                     const PropagationDirection propDir, const MagneticField *magField) const
00355 {
00356   // propagate to next layer
00360   //AnalyticalPropagator aPropagator(magField, propDir);
00361   // Hard coded RungeKutta instead Analytical (avoid bias in TEC), but
00362   // work around TrackPropagation/RungeKutta/interface/RKTestPropagator.h and
00363   // https://www.parashift.com/c++-faq-lite/strange-inheritance.html#faq-23.9
00364   RKTestPropagator bPropagator(magField, propDir); //double tolerance = 5.e-5)
00365   Propagator &aPropagator = bPropagator;
00366   const std::pair<TrajectoryStateOnSurface, double> tsosWithPath =
00367     aPropagator.propagateWithPath(previousTsos, newSurface);
00368 
00369   // stop if propagation wasn't successful
00370   if (!tsosWithPath.first.isValid()) return false;
00371   
00372   nextStep = tsosWithPath.second;
00373   // calculate derivative of reference-track parameters on the actual layer w.r.t. the ones
00374   // on the previous layer (both in global coordinates)
00375   const AnalyticalCurvilinearJacobian aJacobian(previousTsos.globalParameters(), 
00376                                                 tsosWithPath.first.globalPosition(),
00377                                                 tsosWithPath.first.globalMomentum(),
00378                                                 tsosWithPath.second);
00379   const AlgebraicMatrix curvilinearJacobian = asHepMatrix<5,5>(aJacobian.jacobian());
00380 
00381   // jacobian of the track parameters on the previous layer for local->global transformation
00382   const JacobianLocalToCurvilinear startTrafo(previousSurface, previousTsos.localParameters(), *magField);
00383   const AlgebraicMatrix localToCurvilinear = asHepMatrix<5>(startTrafo.jacobian());
00384     
00385   // jacobian of the track parameters on the actual layer for global->local transformation
00386   const JacobianCurvilinearToLocal endTrafo(newSurface, tsosWithPath.first.localParameters(), *magField);
00387   const AlgebraicMatrix curvilinearToLocal = asHepMatrix<5>(endTrafo.jacobian());
00388   
00389   // compute derivative of reference-track parameters on the actual layer w.r.t. the ones on
00390   // the previous layer (both in their local representation)
00391   newCurvlinJacobian = curvilinearJacobian;
00392   newJacobian = curvilinearToLocal * curvilinearJacobian * localToCurvilinear;
00393   newTsos     = tsosWithPath.first;
00394 
00395   return true;
00396 }
00397 
00398 //__________________________________________________________________________________
00399 
00400 void ReferenceTrajectory::fillMeasurementAndError(const TransientTrackingRecHit::ConstRecHitPointer &hitPtr,
00401                                                   unsigned int iRow,
00402                                                   const TrajectoryStateOnSurface &updatedTsos)
00403 {
00404   // get the measurements and their errors, use information updated with tsos if improving
00405   // (GF: Also for measurements or only for errors or do the former not change?)
00406   // GF 10/2008: I doubt that it makes sense to update the hit with the tsos here:
00407   //             That is an analytical extrapolation and not the best guess of the real 
00408   //             track state on the module, but the latter should be better to get the best
00409   //             hit uncertainty estimate!
00410   TransientTrackingRecHit::ConstRecHitPointer newHitPtr(hitPtr->canImproveWithTrack() ?
00411                                                         hitPtr->clone(updatedTsos) : hitPtr);
00412 
00413   const LocalPoint localMeasurement    = newHitPtr->localPosition();
00414   const LocalError localMeasurementCov = newHitPtr->localPositionError();
00415   
00416   theMeasurements[iRow]   = localMeasurement.x();
00417   theMeasurements[iRow+1] = localMeasurement.y();
00418   theMeasurementsCov[iRow][iRow]     = localMeasurementCov.xx();
00419   theMeasurementsCov[iRow][iRow+1]   = localMeasurementCov.xy();
00420   theMeasurementsCov[iRow+1][iRow+1] = localMeasurementCov.yy();
00421   // GF: Should be a loop once the hit dimension is not hardcoded as nMeasPerHit (to be checked):
00422   // for (int i = 0; i < hitPtr->dimension(); ++i) {
00423   //   theMeasurements[iRow+i]   = hitPtr->parameters()[i]; // fixme: parameters() is by value!
00424   //   for (int j = i; j < hitPtr->dimension(); ++j) {
00425   //     theMeasurementsCov[iRow+i][iRow+j] = hitPtr->parametersError()[i][j];
00426   //   }
00427   // }
00428 }
00429 
00430 //__________________________________________________________________________________
00431 
00432 void ReferenceTrajectory::fillDerivatives(const AlgebraicMatrix &projection,
00433                                           const AlgebraicMatrix &fullJacobian,
00434                                           unsigned int iRow)
00435 {
00436   // derivatives of the local coordinates of the reference track w.r.t. to the inital track-parameters
00437   const AlgebraicMatrix projectedJacobian(projection * fullJacobian);
00438   for (int i = 0; i < parameters().num_row(); ++i) {
00439     theDerivatives[iRow  ][i] = projectedJacobian[0][i];
00440     theDerivatives[iRow+1][i] = projectedJacobian[1][i];
00441     // GF: Should be a loop once the hit dimension is not hardcoded as nMeasPerHit (to be checked):
00442     // for (int j = 0; j < projection.num_col(); ++j) {
00443     //   theDerivatives[iRow+j][i] = projectedJacobian[j][i];
00444     // }
00445   }
00446 }
00447 
00448 //__________________________________________________________________________________
00449 
00450 void ReferenceTrajectory::fillTrajectoryPositions(const AlgebraicMatrix &projection, 
00451                                                   const AlgebraicVector &mixedLocalParams, 
00452                                                   unsigned int iRow)
00453 {
00454   // get the local coordinates of the reference trajectory
00455   const AlgebraicVector localPosition(projection * mixedLocalParams);
00456   theTrajectoryPositions[iRow] = localPosition[0];
00457   theTrajectoryPositions[iRow+1] = localPosition[1];
00458   // GF: Should be a loop once the hit dimension is not hardcoded as nMeasPerHit (to be checked):
00459   // for (int j = 0; j < projection.num_col(); ++j) {
00460   //   theTrajectoryPositions[iRow+j] = localPosition[j];
00461   // }
00462 }
00463 
00464 //__________________________________________________________________________________
00465 
00466 bool ReferenceTrajectory::addMaterialEffectsCov(const std::vector<AlgebraicMatrix> &allJacobians, 
00467                                                 const std::vector<AlgebraicMatrix> &allProjections,
00468                                                 const std::vector<AlgebraicSymMatrix> &allCurvatureChanges,
00469                                                 const std::vector<AlgebraicSymMatrix> &allDeltaParameterCovs)
00470 {
00471   // the uncertainty due to multiple scattering is 'transferred' to the error matrix of the hits
00472 
00473   // GF: Needs update once hit dimension is not hardcoded as nMeasPerHit!
00474 
00475   AlgebraicSymMatrix materialEffectsCov(nMeasPerHit * allJacobians.size(), 0);
00476 
00477   // additional covariance-matrix of the measurements due to material-effects (single measurement)
00478   AlgebraicSymMatrix deltaMaterialEffectsCov;
00479 
00480   // additional covariance-matrix of the parameters due to material-effects
00481   AlgebraicSymMatrix paramMaterialEffectsCov(allDeltaParameterCovs[0]); //initialization
00482   // error-propagation to state after energy loss
00483   //GFback  paramMaterialEffectsCov = paramMaterialEffectsCov.similarity(allCurvatureChanges[0]);
00484 
00485   AlgebraicMatrix tempParameterCov;
00486   AlgebraicMatrix tempMeasurementCov;
00487 
00488   for (unsigned int k = 1; k < allJacobians.size(); ++k) {
00489     // error-propagation to next layer
00490     paramMaterialEffectsCov = paramMaterialEffectsCov.similarity(allJacobians[k]);
00491 
00492     // get dependences for the measurements
00493     deltaMaterialEffectsCov = paramMaterialEffectsCov.similarity(allProjections[k]);
00494     materialEffectsCov[nMeasPerHit*k  ][nMeasPerHit*k  ] = deltaMaterialEffectsCov[0][0];
00495     materialEffectsCov[nMeasPerHit*k  ][nMeasPerHit*k+1] = deltaMaterialEffectsCov[0][1];
00496     materialEffectsCov[nMeasPerHit*k+1][nMeasPerHit*k  ] = deltaMaterialEffectsCov[1][0];
00497     materialEffectsCov[nMeasPerHit*k+1][nMeasPerHit*k+1] = deltaMaterialEffectsCov[1][1];
00498 
00499     // GFback add uncertainties for the following layers due to scattering at this layer
00500     paramMaterialEffectsCov += allDeltaParameterCovs[k];
00501     // end GFback
00502     tempParameterCov = paramMaterialEffectsCov;
00503 
00504     // compute "inter-layer-dependencies"
00505     for (unsigned int l = k+1; l < allJacobians.size(); ++l) {
00506       tempParameterCov   = allJacobians[l]   * allCurvatureChanges[l] * tempParameterCov;
00507       tempMeasurementCov = allProjections[l] * tempParameterCov       * allProjections[k].T();
00508 
00509       materialEffectsCov[nMeasPerHit*l][nMeasPerHit*k] = tempMeasurementCov[0][0];
00510       materialEffectsCov[nMeasPerHit*k][nMeasPerHit*l] = tempMeasurementCov[0][0];
00511 
00512       materialEffectsCov[nMeasPerHit*l][nMeasPerHit*k+1] = tempMeasurementCov[0][1];
00513       materialEffectsCov[nMeasPerHit*k+1][nMeasPerHit*l] = tempMeasurementCov[0][1];
00514 
00515       materialEffectsCov[nMeasPerHit*l+1][nMeasPerHit*k] = tempMeasurementCov[1][0];
00516       materialEffectsCov[nMeasPerHit*k][nMeasPerHit*l+1] = tempMeasurementCov[1][0];
00517 
00518       materialEffectsCov[nMeasPerHit*l+1][nMeasPerHit*k+1] = tempMeasurementCov[1][1];
00519       materialEffectsCov[nMeasPerHit*k+1][nMeasPerHit*l+1] = tempMeasurementCov[1][1];
00520     } 
00521     // add uncertainties for the following layers due to scattering at this layer
00522     // GFback paramMaterialEffectsCov += allDeltaParameterCovs[k];    
00523     // error-propagation to state after energy loss
00524     paramMaterialEffectsCov = paramMaterialEffectsCov.similarity(allCurvatureChanges[k]);
00525      
00526   }
00527   theMeasurementsCov += materialEffectsCov;
00528 
00529   return true; // cannot fail
00530 }
00531 
00532 //__________________________________________________________________________________
00533 
00534 bool ReferenceTrajectory::addMaterialEffectsBp(const std::vector<AlgebraicMatrix> &allJacobians, 
00535                                                const std::vector<AlgebraicMatrix> &allProjections,
00536                                                const std::vector<AlgebraicSymMatrix> &allCurvatureChanges,
00537                                                const std::vector<AlgebraicSymMatrix> &allDeltaParameterCovs,
00538                                                const std::vector<AlgebraicMatrix> &allLocalToCurv)
00539 {
00540 //CHK: add material effects using break points
00541   int offsetPar = theNumberOfPars; 
00542   int offsetMeas = nMeasPerHit * allJacobians.size();
00543   int ierr = 0; 
00544 
00545   AlgebraicMatrix tempJacobian;
00546   AlgebraicMatrix MSprojection(2,5);
00547   MSprojection[0][1] = 1;
00548   MSprojection[1][2] = 1;
00549   AlgebraicSymMatrix tempMSCov; 
00550   AlgebraicSymMatrix tempMSCovProj;
00551   AlgebraicMatrix tempMSJacProj;   
00552 
00553   for (unsigned int k = 1; k < allJacobians.size(); ++k) {
00554 // CHK 
00555     int kbp = k-1;
00556     tempJacobian = allJacobians[k] * allCurvatureChanges[k];
00557     tempMSCov = allDeltaParameterCovs[k-1].similarity(allLocalToCurv[k-1]);
00558     tempMSCovProj = tempMSCov.similarity(MSprojection);
00559     theMeasurementsCov[offsetMeas+nMeasPerHit*kbp  ][offsetMeas+nMeasPerHit*kbp  ] = tempMSCovProj[0][0];
00560     theMeasurementsCov[offsetMeas+nMeasPerHit*kbp+1][offsetMeas+nMeasPerHit*kbp+1]=  tempMSCovProj[1][1];
00561     theDerivatives[offsetMeas+nMeasPerHit*kbp  ][offsetPar+2*kbp  ] = 1.0;
00562     theDerivatives[offsetMeas+nMeasPerHit*kbp+1][offsetPar+2*kbp+1] = 1.0 ; 
00563 
00564     tempMSJacProj = (allProjections[k] * ( tempJacobian * allLocalToCurv[k-1].inverse(ierr) )) * MSprojection.T();
00565     if (ierr) {
00566       edm::LogError("Alignment") << "@SUB=ReferenceTrajectory::addMaterialEffectsBp"
00567                                  << "Inversion 1 for break points failed: " << ierr;
00568       return false;
00569     }
00570     theDerivatives[nMeasPerHit*k  ][offsetPar+2*kbp  ] =  tempMSJacProj[0][0];  
00571     theDerivatives[nMeasPerHit*k  ][offsetPar+2*kbp+1] =  tempMSJacProj[0][1]; 
00572     theDerivatives[nMeasPerHit*k+1][offsetPar+2*kbp  ] =  tempMSJacProj[1][0];  
00573     theDerivatives[nMeasPerHit*k+1][offsetPar+2*kbp+1] =  tempMSJacProj[1][1];
00574               
00575     for (unsigned int l = k+1; l < allJacobians.size(); ++l) {
00576 // CHK    
00577       int kbp = k-1;
00578       tempJacobian = allJacobians[l] * allCurvatureChanges[l] * tempJacobian; 
00579       tempMSJacProj = (allProjections[l] * ( tempJacobian * allLocalToCurv[k-1].inverse(ierr) )) * MSprojection.T();
00580       if (ierr) {
00581         edm::LogError("Alignment") << "@SUB=ReferenceTrajectory::addMaterialEffectsBp"
00582                                    << "Inversion 2 for break points failed: " << ierr;
00583         return false;
00584       }
00585       theDerivatives[nMeasPerHit*l  ][offsetPar+2*kbp  ] =  tempMSJacProj[0][0];  
00586       theDerivatives[nMeasPerHit*l  ][offsetPar+2*kbp+1] =  tempMSJacProj[0][1]; 
00587       theDerivatives[nMeasPerHit*l+1][offsetPar+2*kbp  ] =  tempMSJacProj[1][0];  
00588       theDerivatives[nMeasPerHit*l+1][offsetPar+2*kbp+1] =  tempMSJacProj[1][1]; 
00589 
00590     }
00591 
00592   }
00593 
00594   return true;
00595 }
00596 
00597 //__________________________________________________________________________________
00598 
00599 bool ReferenceTrajectory::addMaterialEffectsBrl(const std::vector<AlgebraicMatrix> &allCurvlinJacobians, 
00600                                                 const std::vector<AlgebraicMatrix> &allProjections,
00601                                                 const std::vector<AlgebraicSymMatrix> &allCurvatureChanges,
00602                                                 const std::vector<AlgebraicSymMatrix> &allDeltaParameterCovs,
00603                                                 const std::vector<AlgebraicMatrix> &allLocalToCurv,
00604                                                 const GlobalTrajectoryParameters &gtp)
00605 {
00606 //CHK: add material effects using broken lines
00607 //fine: use exact Jacobians, all detectors  
00608 //broken lines: pair of offsets (u1,u2) = (xt,yt) (in curvilinear frame (q/p,lambda,phi,xt,yt)) at each layer
00609 //              scattering angles (alpha1,alpha2) = (cosLambda*dPhi, dLambda) (cosLambda cancels in Chi2)
00610 //              DU' = (dU'/dU)*DU + (dU'/dAlpha)*DAlpha + (dU'/dQbyp)*DQbyp (propagation of U)
00611 //                  = J*DU + S*DAlpha + d*DQbyp
00612 //           => DAlpha = S^-1 (DU' - J*DU - d*DQbyp)
00613 
00614   int offsetPar = theNumberOfPars;
00615   int offsetMeas = nMeasPerHit*allCurvlinJacobians.size();
00616   int ierr = 0;
00617 
00618   GlobalVector p = gtp.momentum();
00619   double cosLambda = sqrt((p.x()*p.x()+p.y()*p.y())/(p.x()*p.x()+p.y()*p.y()+p.z()*p.z()));
00620 
00621 // transformations Curvilinear <-> BrokenLines
00622   AlgebraicMatrix QbypToCurv(5,1);   // dCurv/dQbyp
00623   QbypToCurv[0][0] = 1.;             // dQbyp/dQbyp
00624   AlgebraicMatrix AngleToCurv(5,2);  // dCurv/dAlpha
00625   AngleToCurv[1][1] = 1.;            // dlambda/dalpha2
00626   AngleToCurv[2][0] = 1./cosLambda;  // dphi/dalpha1
00627   AlgebraicMatrix CurvToAngle(2,5);  // dAlpha/dCurv
00628   CurvToAngle[1][1] = 1.;            // dalpha2/dlambda
00629   CurvToAngle[0][2] = cosLambda;     // dalpha1/dphi
00630   AlgebraicMatrix OffsetToCurv(5,2); // dCurv/dU
00631   OffsetToCurv[3][0] = 1.;           // dxt/du1
00632   OffsetToCurv[4][1] = 1.;           // dyt/du2
00633   AlgebraicMatrix CurvToOffset(2,5); // dU/dCurv
00634   CurvToOffset[0][3] = 1.;           // du1/dxt
00635   CurvToOffset[1][4] = 1.;           // du2/dyt
00636 
00637 // transformations  trajectory to components (Qbyp, U1, U2)
00638   AlgebraicMatrix TrajToQbyp(1,5);
00639   TrajToQbyp[0][0] = 1.;
00640   AlgebraicMatrix TrajToOff1(2,5);
00641   TrajToOff1[0][1] = 1.;
00642   TrajToOff1[1][2] = 1.;
00643   AlgebraicMatrix TrajToOff2(2,5);
00644   TrajToOff2[0][3] = 1.;
00645   TrajToOff2[1][4] = 1.;
00646 
00647   AlgebraicMatrix JacOffsetToAngleC, JacQbypToAngleC;
00648   AlgebraicMatrix JacCurvToOffsetL, JacOffsetToOffsetL, JacAngleToOffsetL, JacQbypToOffsetL, JacOffsetToAngleL;
00649   AlgebraicMatrix JacCurvToOffsetN, JacOffsetToOffsetN, JacAngleToOffsetN, JacQbypToOffsetN, JacOffsetToAngleN;
00650 
00651 // transformation from trajectory to curvilinear parameters
00652 
00653   JacCurvToOffsetN = CurvToOffset * allCurvlinJacobians[1]; // (dU'/dCurv') * (dCurv'/dCurv) @ 2nd point
00654   JacOffsetToOffsetN = JacCurvToOffsetN * OffsetToCurv; // J: (dU'/dU)     = (dU'/dCurv') * (dCurv'/dCurv) * (dCurv/dU)
00655   JacAngleToOffsetN  = JacCurvToOffsetN * AngleToCurv;  // S: (dU'/dAlpha) = (dU'/dCurv') * (dCurv'/dCurv) * (dCurv/dAlpha)
00656   JacQbypToOffsetN   = JacCurvToOffsetN * QbypToCurv;   // d: (dU'/dQbyp)  = (dU'/dCurv') * (dCurv'/dCurv) * (dCurv/dQbyp)
00657   JacOffsetToAngleN  = JacAngleToOffsetN.inverse(ierr); // W
00658   if (ierr) {
00659      edm::LogError("Alignment") << "@SUB=ReferenceTrajectory::addMaterialEffectsBrl"
00660                                 << "Inversion 1 for fine broken lines failed: " << ierr;
00661      return false;
00662   }
00663   JacOffsetToAngleC = -(JacOffsetToAngleN * JacOffsetToOffsetN); // (dAlpha/dU)
00664   JacQbypToAngleC   = -(JacOffsetToAngleN * JacQbypToOffsetN);   // (dAlpha/dQbyp)
00665   // (dAlpha/dTraj) = (dAlpha/dQbyp) * (dQbyp/dTraj) + (dAlpha/dU1) * (dU1/dTraj) + (dAlpha/dU2) * (dU2/dTraj)
00666   AlgebraicMatrix JacTrajToAngle = JacQbypToAngleC * TrajToQbyp + JacOffsetToAngleC * TrajToOff1 + JacOffsetToAngleN * TrajToOff2;
00667   // (dCurv/dTraj) = (dCurv/dQbyp) * (dQbyp/dTraj) + (dCurv/dAlpha) * (dAlpha/dTraj) + (dCurv/dU) * (dU/dTraj)
00668   theInnerTrajectoryToCurvilinear = QbypToCurv * TrajToQbyp + AngleToCurv * JacTrajToAngle + OffsetToCurv * TrajToOff1;
00669   theInnerLocalToTrajectory = theInnerTrajectoryToCurvilinear.inverse(ierr) * allLocalToCurv[0];
00670   if (ierr) {
00671     edm::LogError("Alignment") << "@SUB=ReferenceTrajectory::addMaterialEffectsBrl"
00672                                << "Inversion 2 for fine broken lines failed: " << ierr;
00673     return false;
00674   }
00675 
00676   AlgebraicMatrix tempJacobian(allCurvatureChanges[0]);
00677   AlgebraicSymMatrix tempMSCov;
00678   AlgebraicSymMatrix tempMSCovProj;
00679   AlgebraicMatrix tempJacL, tempJacN;
00680   AlgebraicMatrix JacOffsetToMeas;
00681 
00682 // measurements from hits  
00683   for (unsigned int k = 0; k < allCurvlinJacobians.size(); ++k) {
00684 //  (dMeas/dU) = (dMeas/dLoc) * (dLoc/dCurv) * (dCurv/dU)
00685     JacOffsetToMeas = (allProjections[k] * allLocalToCurv[k].inverse(ierr) ) * OffsetToCurv;
00686     if (ierr) {
00687        edm::LogError("Alignment") << "@SUB=ReferenceTrajectory::addMaterialEffectsBrl"
00688                                    << "Inversion 3 for fine broken lines failed: " << ierr;
00689        return false;
00690     }
00691     theDerivatives[nMeasPerHit*k  ][offsetPar+2*k  ] =  JacOffsetToMeas[0][0];
00692     theDerivatives[nMeasPerHit*k  ][offsetPar+2*k+1] =  JacOffsetToMeas[0][1];
00693     theDerivatives[nMeasPerHit*k+1][offsetPar+2*k  ] =  JacOffsetToMeas[1][0];
00694     theDerivatives[nMeasPerHit*k+1][offsetPar+2*k+1] =  JacOffsetToMeas[1][1];
00695   }
00696 
00697 // measurement of MS kink  
00698   for (unsigned int k = 1; k < allCurvlinJacobians.size()-1; ++k) {
00699 // CHK 
00700     int iMsMeas = k-1; 
00701     int l    = k-1; // last hit
00702     int n    = k+1; // next hit
00703 
00704 // amount of multiple scattering in layer k  (angular error perp to direction)  
00705     tempMSCov = allDeltaParameterCovs[k].similarity(allLocalToCurv[k]);
00706     tempMSCovProj = tempMSCov.similarity(CurvToAngle);
00707     theMeasurementsCov[offsetMeas+nMeasPerHit*iMsMeas  ][offsetMeas+nMeasPerHit*iMsMeas  ] = tempMSCovProj[1][1];
00708     theMeasurementsCov[offsetMeas+nMeasPerHit*iMsMeas+1][offsetMeas+nMeasPerHit*iMsMeas+1] = tempMSCovProj[0][0];
00709 
00710 // transformation matices for offsets ( l <- k -> n )
00711     tempJacL = allCurvlinJacobians[k] * tempJacobian;
00712     JacCurvToOffsetL = CurvToOffset * tempJacL.inverse(ierr); // (dU'/dCurv') * (dCurv'/dCurv) @ last point
00713 
00714     if (ierr) {
00715        edm::LogError("Alignment") << "@SUB=ReferenceTrajectory::addMaterialEffectsBrl"
00716                                    << "Inversion 4 for fine broken lines failed: " << ierr;
00717        return false;
00718     }
00719     JacOffsetToOffsetL = JacCurvToOffsetL * OffsetToCurv; // J-: (dU'/dU)     = (dU'/dCurv') * (dCurv'/dCurv) * (dCurv/dU)
00720     JacAngleToOffsetL  = JacCurvToOffsetL * AngleToCurv;  // S-: (dU'/dAlpha) = (dU'/dCurv') * (dCurv'/dCurv) * (dCurv/dAlpha)
00721     JacQbypToOffsetL   = JacCurvToOffsetL * QbypToCurv;   // d-: (dU'/dQbyp)  = (dU'/dCurv') * (dCurv'/dCurv) * (dCurv/dQbyp)
00722     JacOffsetToAngleL  =-JacAngleToOffsetL.inverse(ierr); // W-
00723     if (ierr) {
00724        edm::LogError("Alignment") << "@SUB=ReferenceTrajectory::addMaterialEffectsBrl"
00725                                    << "Inversion 5 for fine broken lines failed: " << ierr;
00726        return false;
00727     }
00728     tempJacobian = tempJacobian * allCurvatureChanges[k];
00729     tempJacN = allCurvlinJacobians[n] * tempJacobian;
00730     JacCurvToOffsetN = CurvToOffset * tempJacN; // (dU'/dCurv') * (dCurv'/dCurv) @ next point
00731     JacOffsetToOffsetN = JacCurvToOffsetN * OffsetToCurv; // J+: (dU'/dU)     = (dU'/dCurv') * (dCurv'/dCurv) * (dCurv/dU)
00732     JacAngleToOffsetN  = JacCurvToOffsetN * AngleToCurv;  // S+: (dU'/dAlpha) = (dU'/dCurv') * (dCurv'/dCurv) * (dCurv/dAlpha)
00733     JacQbypToOffsetN   = JacCurvToOffsetN * QbypToCurv;   // d+: (dU'/dQbyp)  = (dU'/dCurv') * (dCurv'/dCurv) * (dCurv/dQbyp)
00734     JacOffsetToAngleN  = JacAngleToOffsetN.inverse(ierr); // W+
00735     if (ierr) {
00736        edm::LogError("Alignment") << "@SUB=ReferenceTrajectory::addMaterialEffectsBrl"
00737                                    << "Inversion 6 for fine broken lines failed: " << ierr;
00738        return false;
00739     }
00740     JacOffsetToAngleC = -(JacOffsetToAngleL * JacOffsetToOffsetL + JacOffsetToAngleN * JacOffsetToOffsetN);
00741     JacQbypToAngleC   = -(JacOffsetToAngleL * JacQbypToOffsetL   + JacOffsetToAngleN * JacQbypToOffsetN);
00742 
00743     // bending    
00744     theDerivatives[offsetMeas+nMeasPerHit*iMsMeas  ][              0] = JacQbypToAngleC[0][0];
00745     theDerivatives[offsetMeas+nMeasPerHit*iMsMeas+1][              0] = JacQbypToAngleC[1][0];
00746     // last layer
00747     theDerivatives[offsetMeas+nMeasPerHit*iMsMeas  ][offsetPar+2*l  ] = JacOffsetToAngleL[0][0];
00748     theDerivatives[offsetMeas+nMeasPerHit*iMsMeas  ][offsetPar+2*l+1] = JacOffsetToAngleL[0][1];
00749     theDerivatives[offsetMeas+nMeasPerHit*iMsMeas+1][offsetPar+2*l  ] = JacOffsetToAngleL[1][0];
00750     theDerivatives[offsetMeas+nMeasPerHit*iMsMeas+1][offsetPar+2*l+1] = JacOffsetToAngleL[1][1];
00751     // current layer
00752     theDerivatives[offsetMeas+nMeasPerHit*iMsMeas  ][offsetPar+2*k  ] = JacOffsetToAngleC[0][0];
00753     theDerivatives[offsetMeas+nMeasPerHit*iMsMeas  ][offsetPar+2*k+1] = JacOffsetToAngleC[0][1];
00754     theDerivatives[offsetMeas+nMeasPerHit*iMsMeas+1][offsetPar+2*k  ] = JacOffsetToAngleC[1][0];
00755     theDerivatives[offsetMeas+nMeasPerHit*iMsMeas+1][offsetPar+2*k+1] = JacOffsetToAngleC[1][1];
00756 
00757     // next layer
00758     theDerivatives[offsetMeas+nMeasPerHit*iMsMeas  ][offsetPar+2*n  ] = JacOffsetToAngleN[0][0];
00759     theDerivatives[offsetMeas+nMeasPerHit*iMsMeas  ][offsetPar+2*n+1] = JacOffsetToAngleN[0][1];
00760     theDerivatives[offsetMeas+nMeasPerHit*iMsMeas+1][offsetPar+2*n  ] = JacOffsetToAngleN[1][0];
00761     theDerivatives[offsetMeas+nMeasPerHit*iMsMeas+1][offsetPar+2*n+1] = JacOffsetToAngleN[1][1];
00762 
00763   }
00764 
00765   return true;
00766 }
00767 
00768 bool ReferenceTrajectory::addMaterialEffectsBrl(const std::vector<AlgebraicMatrix> &allProjections,
00769                                                 const std::vector<AlgebraicSymMatrix> &allDeltaParameterCovs,
00770                                                 const std::vector<AlgebraicMatrix> &allLocalToCurv,
00771                                                 const std::vector<double> &allSteps,
00772                                                 const GlobalTrajectoryParameters &gtp,
00773                                                 const double minStep)
00774 {
00775 //CHK: add material effects using broken lines
00776 //BrokenLinesCoarse: combine close by detectors, 
00777 //                   use approximate Jacobians from Steps (limit Qbyp -> 0),
00778 //                   bending only in RPhi (B=(0,0,Bz)), no energy loss correction
00779 
00780   int offsetPar = theNumberOfPars; 
00781   int offsetMeas = nMeasPerHit*allSteps.size();
00782   int ierr = 0;
00783 
00784   GlobalVector p = gtp.momentum();
00785   double cosLambda = sqrt((p.x()*p.x()+p.y()*p.y())/(p.x()*p.x()+p.y()*p.y()+p.z()*p.z()));
00786   double bFac = -gtp.magneticFieldInInverseGeV(gtp.position()).mag();
00787 
00788   // transformation from trajectory to curvilinear parameters at refTsos
00789   double delta (1.0/allSteps[1]);    
00790   theInnerTrajectoryToCurvilinear[0][0] = 1;
00791   theInnerTrajectoryToCurvilinear[1][2] = -delta;  
00792   theInnerTrajectoryToCurvilinear[1][4] =  delta;    
00793   theInnerTrajectoryToCurvilinear[2][0] = -0.5*bFac/delta;
00794   theInnerTrajectoryToCurvilinear[2][1] = -delta/cosLambda;
00795   theInnerTrajectoryToCurvilinear[2][3] =  delta/cosLambda;
00796   theInnerTrajectoryToCurvilinear[3][1] = 1;
00797   theInnerTrajectoryToCurvilinear[4][2] = 1;
00798   theInnerLocalToTrajectory = theInnerTrajectoryToCurvilinear.inverse(ierr) * allLocalToCurv[0];
00799   if (ierr) {
00800     edm::LogError("Alignment") << "@SUB=ReferenceTrajectory::addMaterialEffectsBrl"
00801                                << "Inversion 1 for coarse broken lines failed: " << ierr;
00802     return false;
00803   }
00804 
00805   AlgebraicMatrix CurvToAngle(2,5);  // dAlpha/dCurv
00806   CurvToAngle[1][1] = 1.;            // dalpha2/dlambda
00807   CurvToAngle[0][2] = cosLambda;     // dalpha1/dphi
00808   AlgebraicMatrix OffsetToCurv(5,2); // dCurv/dU
00809   OffsetToCurv[3][0] = 1.;           // dxt/du1
00810   OffsetToCurv[4][1] = 1.;           // dyt/du2
00811 
00812   AlgebraicSymMatrix tempMSCov;
00813   AlgebraicSymMatrix tempMSCovProj;
00814   AlgebraicMatrix JacOffsetToMeas;
00815 
00816   // combine closeby detectors into single plane
00817   std::vector<unsigned int> first(allSteps.size());
00818   std::vector<unsigned int> last (allSteps.size());
00819   std::vector<unsigned int> plane(allSteps.size());
00820   std::vector<double> sPlane(allSteps.size());
00821   unsigned int nPlane = 0;
00822   double sTot = 0;
00823 
00824   for (unsigned int k = 1; k < allSteps.size(); ++k) {
00825     sTot += allSteps[k];
00826     if (fabs(allSteps[k])>minStep) { nPlane += 1; first[nPlane] = k; }
00827     last[nPlane] = k;
00828     plane[k] = nPlane;
00829     sPlane[nPlane] += sTot;
00830   }
00831   if (nPlane < 2) return false; // pathological cases: need at least 2 planes
00832 
00833   theNumberOfMsPars = 2*(nPlane+1);
00834   theNumberOfMsMeas = 2*(nPlane-1);// unsigned underflow for nPlane == 0...
00835   for (unsigned int k = 0; k <= nPlane; ++k) { sPlane[k] /= (double) (last[k]-first[k]+1); }
00836 
00837   // measurements from hits
00838   sTot = 0; 
00839   for (unsigned int k = 0; k < allSteps.size(); ++k) {
00840     sTot += allSteps[k];
00841 //  (dMeas/dU) = (dMeas/dLoc) * (dLoc/dCurv) * (dCurv/dU)
00842     JacOffsetToMeas = (allProjections[k] * allLocalToCurv[k].inverse(ierr) ) * OffsetToCurv;
00843     if (ierr) {
00844       edm::LogError("Alignment") << "@SUB=ReferenceTrajectory::addMaterialEffectsBrl"
00845                                  << "Inversion 2 for coarse broken lines failed: " << ierr;
00846       return false;
00847     }
00848 
00849     unsigned int iPlane = plane[k];
00850     if (last[iPlane] == first[iPlane])
00851     { // single plane
00852       theDerivatives[nMeasPerHit*k  ][offsetPar+2*iPlane  ] =  JacOffsetToMeas[0][0];
00853       theDerivatives[nMeasPerHit*k  ][offsetPar+2*iPlane+1] =  JacOffsetToMeas[0][1];
00854       theDerivatives[nMeasPerHit*k+1][offsetPar+2*iPlane  ] =  JacOffsetToMeas[1][0];
00855       theDerivatives[nMeasPerHit*k+1][offsetPar+2*iPlane+1] =  JacOffsetToMeas[1][1];
00856     } else
00857     { // combined plane: (linear) interpolation
00858       unsigned int jPlane; // neighbor plane for interpolation
00859       if (fabs(sTot) < fabs(sPlane[iPlane])) { jPlane = (iPlane>0) ? iPlane - 1 : 1; }
00860       else  { jPlane = (iPlane<nPlane) ? iPlane + 1 : nPlane -1 ;}
00861       // interpolation weights
00862       double sDiff = sPlane[iPlane] - sPlane[jPlane];
00863       double iFrac = (sTot - sPlane[jPlane]) / sDiff;
00864       double jFrac = 1.0 - iFrac;
00865       theDerivatives[nMeasPerHit*k  ][offsetPar+2*iPlane  ] =  JacOffsetToMeas[0][0]*iFrac;
00866       theDerivatives[nMeasPerHit*k  ][offsetPar+2*iPlane+1] =  JacOffsetToMeas[0][1]*iFrac;
00867       theDerivatives[nMeasPerHit*k+1][offsetPar+2*iPlane  ] =  JacOffsetToMeas[1][0]*iFrac;
00868       theDerivatives[nMeasPerHit*k+1][offsetPar+2*iPlane+1] =  JacOffsetToMeas[1][1]*iFrac;
00869       theDerivatives[nMeasPerHit*k  ][offsetPar+2*jPlane  ] =  JacOffsetToMeas[0][0]*jFrac;
00870       theDerivatives[nMeasPerHit*k  ][offsetPar+2*jPlane+1] =  JacOffsetToMeas[0][1]*jFrac;
00871       theDerivatives[nMeasPerHit*k+1][offsetPar+2*jPlane  ] =  JacOffsetToMeas[1][0]*jFrac;
00872       theDerivatives[nMeasPerHit*k+1][offsetPar+2*jPlane+1] =  JacOffsetToMeas[1][1]*jFrac;
00873       // 2nd order neglected
00874       // theDerivatives[nMeasPerHit*k  ][                   0] = -0.5*bFac*sDiff*iFrac*sDiff*jFrac*cosLambda;
00875     }
00876   }
00877 
00878 // measurement of MS kink
00879   for (unsigned int i = 1; i < nPlane; ++i) {
00880 // CHK 
00881     int iMsMeas = i-1;
00882     int l    = i-1; // last hit
00883     int n    = i+1; // next hit
00884 
00885 // amount of multiple scattering in plane k
00886     for (unsigned int k = first[i]; k <= last[i]; ++k) {
00887       tempMSCov = allDeltaParameterCovs[k].similarity(allLocalToCurv[k]);
00888       tempMSCovProj = tempMSCov.similarity(CurvToAngle);
00889       theMeasurementsCov[offsetMeas+nMeasPerHit*iMsMeas  ][offsetMeas+nMeasPerHit*iMsMeas  ] += tempMSCovProj[0][0];
00890       theMeasurementsCov[offsetMeas+nMeasPerHit*iMsMeas+1][offsetMeas+nMeasPerHit*iMsMeas+1] += tempMSCovProj[1][1];
00891     }
00892 // broken line measurements for layer k, correlations between both planes neglected
00893     double stepK = sPlane[i] - sPlane[l];
00894     double stepN = sPlane[n] - sPlane[i];
00895     double deltaK (1.0/stepK);
00896     double deltaN (1.0/stepN);
00897     // bending (only in RPhi)
00898     theDerivatives[offsetMeas+nMeasPerHit*iMsMeas  ][              0] = -0.5*bFac*(stepK+stepN)*cosLambda;
00899     // last layer
00900     theDerivatives[offsetMeas+nMeasPerHit*iMsMeas  ][offsetPar+2*l  ] = deltaK;
00901     theDerivatives[offsetMeas+nMeasPerHit*iMsMeas+1][offsetPar+2*l+1] = deltaK;
00902     // current layer
00903     theDerivatives[offsetMeas+nMeasPerHit*iMsMeas  ][offsetPar+2*i  ] = -(deltaK + deltaN);
00904     theDerivatives[offsetMeas+nMeasPerHit*iMsMeas+1][offsetPar+2*i+1] = -(deltaK + deltaN);
00905     // next layer
00906     theDerivatives[offsetMeas+nMeasPerHit*iMsMeas  ][offsetPar+2*n  ] = deltaN;
00907     theDerivatives[offsetMeas+nMeasPerHit*iMsMeas+1][offsetPar+2*n+1] = deltaN;
00908   }
00909 
00910   return true;
00911 }