00001
00002
00003
00004
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 "Alignment/ReferenceTrajectories/interface/BeamSpotTransientTrackingRecHit.h"
00044 #include "Alignment/ReferenceTrajectories/interface/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
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
00113 std::auto_ptr<MaterialEffectsUpdator> aMaterialEffectsUpdator
00114 (this->createUpdator(materialEffects, mass));
00115 if (!aMaterialEffectsUpdator.get()) return false;
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
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
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 if (!tsctbl.isValid()) {
00155 edm::LogError("Alignment") << "@SUB=ReferenceTrajectory::construct"
00156 << "TrajectoryStateClostestToBeamLine invalid. Skip track.";
00157 return false;
00158 }
00159
00160 FreeTrajectoryState pcaFts = tsctbl.trackStateAtPCA();
00161 GlobalVector bd(beamSpot.dxdz(), beamSpot.dydz(), 1.0);
00162
00163
00164 AnalyticalPropagator propagator(magField);
00165 std::pair< TrajectoryStateOnSurface, double > tsosWithPath =
00166 propagator.propagateWithPath(pcaFts, refTsos.surface());
00167
00168 if (!tsosWithPath.first.isValid()) return false;
00169
00170 GlobalVector momDir(pcaFts.momentum());
00171 GlobalVector perpDir(bd.cross(momDir));
00172 BoundPlane::RotationType rotation(perpDir, bd);
00173
00174 BeamSpotGeomDet * bsGeom = new BeamSpotGeomDet(BoundPlane::build(bs, rotation, OpenBounds()));
00175
00176
00177 theRefTsos = TrajectoryStateOnSurface(pcaFts, bsGeom->surface());
00178
00179 TransientTrackingRecHit::ConstRecHitPointer bsRecHit =
00180 new BeamSpotTransientTrackingRecHit(beamSpot,
00181 bsGeom,
00182 theRefTsos.freeState()->momentum().phi());
00183 allRecHits.push_back(bsRecHit);
00184
00185 }
00186
00187
00188 TransientTrackingRecHit::ConstRecHitContainer::const_iterator itRecHit;
00189 for ( itRecHit = recHits.begin(); itRecHit != recHits.end(); ++itRecHit ) {
00190 const TransientTrackingRecHit::ConstRecHitPointer &hitPtr = *itRecHit;
00191 allRecHits.push_back(hitPtr);
00192 }
00193
00194 for ( itRecHit = allRecHits.begin(); itRecHit != allRecHits.end(); ++itRecHit ) {
00195
00196 const TransientTrackingRecHit::ConstRecHitPointer &hitPtr = *itRecHit;
00197 theRecHits.push_back(hitPtr);
00198
00199 if (0 == iRow) {
00200
00201
00202
00203 fullJacobian = AlgebraicMatrix(theParameters.num_row(), theParameters.num_row(), 1);
00204 allJacobians.push_back(fullJacobian);
00205 theTsosVec.push_back(theRefTsos);
00206 const JacobianLocalToCurvilinear startTrafo(hitPtr->det()->surface(), theRefTsos.localParameters(), *magField);
00207 const AlgebraicMatrix localToCurvilinear = asHepMatrix<5>(startTrafo.jacobian());
00208 if (materialEffects <= breakPoints) {
00209 theInnerTrajectoryToCurvilinear = asHepMatrix<5>(startTrafo.jacobian());
00210 theInnerLocalToTrajectory = AlgebraicMatrix(5, 5, 1);
00211 }
00212 allLocalToCurv.push_back(localToCurvilinear);
00213 allSteps.push_back(0.);
00214 allCurvlinJacobians.push_back(firstCurvlinJacobian);
00215
00216 } else {
00217
00218 AlgebraicMatrix nextJacobian;
00219 AlgebraicMatrix nextCurvlinJacobian;
00220 double nextStep = 0.;
00221 TrajectoryStateOnSurface nextTsos;
00222
00223 if (!this->propagate(previousHitPtr->det()->surface(), previousTsos,
00224 hitPtr->det()->surface(), nextTsos,
00225 nextJacobian, nextCurvlinJacobian, nextStep, propDir, magField)) {
00226 return false;
00227 }
00228
00229 allJacobians.push_back(nextJacobian);
00230 fullJacobian = nextJacobian * previousChangeInCurvature * fullJacobian;
00231 theTsosVec.push_back(nextTsos);
00232
00233 const JacobianLocalToCurvilinear startTrafo(hitPtr->det()->surface(), nextTsos.localParameters(), *magField);
00234 const AlgebraicMatrix localToCurvilinear = asHepMatrix<5>(startTrafo.jacobian());
00235 allLocalToCurv.push_back(localToCurvilinear);
00236 if (nextStep == 0.) {
00237 edm::LogError("Alignment") << "@SUB=ReferenceTrajectory::construct"
00238 << "step 0. from id " << previousHitPtr->geographicalId()
00239 << " to " << hitPtr->det()->geographicalId() << ".";
00240
00241 if (materialEffects == brokenLinesFine) {
00242 edm::LogError("Alignment") << "@SUB=ReferenceTrajectory::construct" << "Skip track.";
00243 return false;
00244 }
00245 }
00246 allSteps.push_back(nextStep);
00247 allCurvlinJacobians.push_back(nextCurvlinJacobian);
00248
00249 }
00250
00251
00252
00253 const TrajectoryStateOnSurface tmpTsos(theTsosVec.back().localParameters(), zeroErrors,
00254 theTsosVec.back().surface(), magField, surfaceSide);
00255 const TrajectoryStateOnSurface updatedTsos = aMaterialEffectsUpdator->updateState(tmpTsos, propDir);
00256
00257 if ( !updatedTsos.isValid() ) return false;
00258
00259 if ( theTsosVec.back().localParameters().charge() )
00260 {
00261 previousChangeInCurvature[0][0] = updatedTsos.localParameters().signedInverseMomentum()
00262 / theTsosVec.back().localParameters().signedInverseMomentum();
00263 }
00264
00265
00266 allDeltaParameterCovs.push_back( asHepMatrix<5>(updatedTsos.localError().matrix()) );
00267 allCurvatureChanges.push_back(previousChangeInCurvature);
00268
00269
00270 if ( useRecHit( hitPtr ) ) { allProjections.push_back( hitPtr->projectionMatrix() ); }
00271 else { allProjections.push_back( AlgebraicMatrix( 2, 5, 0) ); }
00272
00273
00274
00275 previousHitPtr = hitPtr;
00276 previousTsos = TrajectoryStateOnSurface(updatedTsos.globalParameters(),
00277 updatedTsos.surface(), surfaceSide);
00278
00279 if (materialEffects < brokenLinesCoarse) {
00280 this->fillDerivatives(allProjections.back(), fullJacobian, iRow);
00281 }
00282
00283 AlgebraicVector mixedLocalParams = asHepVector<5>(theTsosVec.back().localParameters().mixedFormatVector());
00284 this->fillTrajectoryPositions(allProjections.back(), mixedLocalParams, iRow);
00285 if ( useRecHit( hitPtr ) ) this->fillMeasurementAndError(hitPtr, iRow, updatedTsos);
00286
00287 iRow += nMeasPerHit;
00288 }
00289
00290 bool msOK = true;
00291 switch (materialEffects) {
00292 case none:
00293 break;
00294 case multipleScattering:
00295 case energyLoss:
00296 case combined:
00297 msOK = this->addMaterialEffectsCov(allJacobians, allProjections, allCurvatureChanges,
00298 allDeltaParameterCovs);
00299 break;
00300 case breakPoints:
00301 msOK = this->addMaterialEffectsBp(allJacobians, allProjections, allCurvatureChanges,
00302 allDeltaParameterCovs, allLocalToCurv);
00303 break;
00304 case brokenLinesCoarse:
00305 msOK = this->addMaterialEffectsBrl(allProjections, allDeltaParameterCovs, allLocalToCurv,
00306 allSteps, refTsos.globalParameters());
00307 break;
00308 case brokenLinesFine:
00309 msOK = this->addMaterialEffectsBrl(allCurvlinJacobians, allProjections, allCurvatureChanges,
00310 allDeltaParameterCovs, allLocalToCurv, refTsos.globalParameters());
00311 }
00312 if (!msOK) return false;
00313
00314 if (refTsos.hasError()) {
00315 AlgebraicSymMatrix parameterCov = asHepMatrix<5>(refTsos.localError().matrix());
00316 AlgebraicMatrix parDeriv;
00317 if (theNumberOfVirtualPars>0) {
00318 parDeriv = theDerivatives.sub( 1, nMeasPerHit*allJacobians.size(), 1, theParameters.num_row() );
00319 } else {
00320 parDeriv = theDerivatives;
00321 }
00322 theTrajectoryPositionCov = parameterCov.similarity(parDeriv);
00323 } else {
00324 theTrajectoryPositionCov = AlgebraicSymMatrix(theDerivatives.num_row(), 1);
00325 }
00326
00327 return true;
00328 }
00329
00330
00331
00332 MaterialEffectsUpdator*
00333 ReferenceTrajectory::createUpdator(MaterialEffects materialEffects, double mass) const
00334 {
00335 switch (materialEffects) {
00336
00337
00338 case none:
00339 case multipleScattering:
00340 return new MultipleScatteringUpdator(mass);
00341 case energyLoss:
00342 return new EnergyLossUpdator(mass);
00343 case combined:
00344 return new CombinedMaterialEffectsUpdator(mass);
00345 case breakPoints:
00346 return new CombinedMaterialEffectsUpdator(mass);
00347 case brokenLinesCoarse:
00348 case brokenLinesFine:
00349 return new CombinedMaterialEffectsUpdator(mass);
00350 }
00351
00352 return 0;
00353 }
00354
00355
00356
00357 bool ReferenceTrajectory::propagate(const BoundPlane &previousSurface, const TrajectoryStateOnSurface &previousTsos,
00358 const BoundPlane &newSurface, TrajectoryStateOnSurface &newTsos, AlgebraicMatrix &newJacobian,
00359 AlgebraicMatrix &newCurvlinJacobian, double &nextStep,
00360 const PropagationDirection propDir, const MagneticField *magField) const
00361 {
00362
00366
00367
00368
00369
00370 RKTestPropagator bPropagator(magField, propDir);
00371 Propagator &aPropagator = bPropagator;
00372 const std::pair<TrajectoryStateOnSurface, double> tsosWithPath =
00373 aPropagator.propagateWithPath(previousTsos, newSurface);
00374
00375
00376 if (!tsosWithPath.first.isValid()) return false;
00377
00378 nextStep = tsosWithPath.second;
00379
00380
00381 const AnalyticalCurvilinearJacobian aJacobian(previousTsos.globalParameters(),
00382 tsosWithPath.first.globalPosition(),
00383 tsosWithPath.first.globalMomentum(),
00384 tsosWithPath.second);
00385 const AlgebraicMatrix curvilinearJacobian = asHepMatrix<5,5>(aJacobian.jacobian());
00386
00387
00388 const JacobianLocalToCurvilinear startTrafo(previousSurface, previousTsos.localParameters(), *magField);
00389 const AlgebraicMatrix localToCurvilinear = asHepMatrix<5>(startTrafo.jacobian());
00390
00391
00392 const JacobianCurvilinearToLocal endTrafo(newSurface, tsosWithPath.first.localParameters(), *magField);
00393 const AlgebraicMatrix curvilinearToLocal = asHepMatrix<5>(endTrafo.jacobian());
00394
00395
00396
00397 newCurvlinJacobian = curvilinearJacobian;
00398 newJacobian = curvilinearToLocal * curvilinearJacobian * localToCurvilinear;
00399 newTsos = tsosWithPath.first;
00400
00401 return true;
00402 }
00403
00404
00405
00406 void ReferenceTrajectory::fillMeasurementAndError(const TransientTrackingRecHit::ConstRecHitPointer &hitPtr,
00407 unsigned int iRow,
00408 const TrajectoryStateOnSurface &updatedTsos)
00409 {
00410
00411
00412
00413
00414
00415
00416 TransientTrackingRecHit::ConstRecHitPointer newHitPtr(hitPtr->canImproveWithTrack() ?
00417 hitPtr->clone(updatedTsos) : hitPtr);
00418
00419 const LocalPoint localMeasurement = newHitPtr->localPosition();
00420 const LocalError localMeasurementCov = newHitPtr->localPositionError();
00421
00422 theMeasurements[iRow] = localMeasurement.x();
00423 theMeasurements[iRow+1] = localMeasurement.y();
00424 theMeasurementsCov[iRow][iRow] = localMeasurementCov.xx();
00425 theMeasurementsCov[iRow][iRow+1] = localMeasurementCov.xy();
00426 theMeasurementsCov[iRow+1][iRow+1] = localMeasurementCov.yy();
00427
00428
00429
00430
00431
00432
00433
00434 }
00435
00436
00437
00438 void ReferenceTrajectory::fillDerivatives(const AlgebraicMatrix &projection,
00439 const AlgebraicMatrix &fullJacobian,
00440 unsigned int iRow)
00441 {
00442
00443 const AlgebraicMatrix projectedJacobian(projection * fullJacobian);
00444 for (int i = 0; i < parameters().num_row(); ++i) {
00445 theDerivatives[iRow ][i] = projectedJacobian[0][i];
00446 theDerivatives[iRow+1][i] = projectedJacobian[1][i];
00447
00448
00449
00450
00451 }
00452 }
00453
00454
00455
00456 void ReferenceTrajectory::fillTrajectoryPositions(const AlgebraicMatrix &projection,
00457 const AlgebraicVector &mixedLocalParams,
00458 unsigned int iRow)
00459 {
00460
00461 const AlgebraicVector localPosition(projection * mixedLocalParams);
00462 theTrajectoryPositions[iRow] = localPosition[0];
00463 theTrajectoryPositions[iRow+1] = localPosition[1];
00464
00465
00466
00467
00468 }
00469
00470
00471
00472 bool ReferenceTrajectory::addMaterialEffectsCov(const std::vector<AlgebraicMatrix> &allJacobians,
00473 const std::vector<AlgebraicMatrix> &allProjections,
00474 const std::vector<AlgebraicSymMatrix> &allCurvatureChanges,
00475 const std::vector<AlgebraicSymMatrix> &allDeltaParameterCovs)
00476 {
00477
00478
00479
00480
00481 AlgebraicSymMatrix materialEffectsCov(nMeasPerHit * allJacobians.size(), 0);
00482
00483
00484 AlgebraicSymMatrix deltaMaterialEffectsCov;
00485
00486
00487 AlgebraicSymMatrix paramMaterialEffectsCov(allDeltaParameterCovs[0]);
00488
00489
00490
00491 AlgebraicMatrix tempParameterCov;
00492 AlgebraicMatrix tempMeasurementCov;
00493
00494 for (unsigned int k = 1; k < allJacobians.size(); ++k) {
00495
00496 paramMaterialEffectsCov = paramMaterialEffectsCov.similarity(allJacobians[k]);
00497
00498
00499 deltaMaterialEffectsCov = paramMaterialEffectsCov.similarity(allProjections[k]);
00500 materialEffectsCov[nMeasPerHit*k ][nMeasPerHit*k ] = deltaMaterialEffectsCov[0][0];
00501 materialEffectsCov[nMeasPerHit*k ][nMeasPerHit*k+1] = deltaMaterialEffectsCov[0][1];
00502 materialEffectsCov[nMeasPerHit*k+1][nMeasPerHit*k ] = deltaMaterialEffectsCov[1][0];
00503 materialEffectsCov[nMeasPerHit*k+1][nMeasPerHit*k+1] = deltaMaterialEffectsCov[1][1];
00504
00505
00506 paramMaterialEffectsCov += allDeltaParameterCovs[k];
00507
00508 tempParameterCov = paramMaterialEffectsCov;
00509
00510
00511 for (unsigned int l = k+1; l < allJacobians.size(); ++l) {
00512 tempParameterCov = allJacobians[l] * allCurvatureChanges[l] * tempParameterCov;
00513 tempMeasurementCov = allProjections[l] * tempParameterCov * allProjections[k].T();
00514
00515 materialEffectsCov[nMeasPerHit*l][nMeasPerHit*k] = tempMeasurementCov[0][0];
00516 materialEffectsCov[nMeasPerHit*k][nMeasPerHit*l] = tempMeasurementCov[0][0];
00517
00518 materialEffectsCov[nMeasPerHit*l][nMeasPerHit*k+1] = tempMeasurementCov[0][1];
00519 materialEffectsCov[nMeasPerHit*k+1][nMeasPerHit*l] = tempMeasurementCov[0][1];
00520
00521 materialEffectsCov[nMeasPerHit*l+1][nMeasPerHit*k] = tempMeasurementCov[1][0];
00522 materialEffectsCov[nMeasPerHit*k][nMeasPerHit*l+1] = tempMeasurementCov[1][0];
00523
00524 materialEffectsCov[nMeasPerHit*l+1][nMeasPerHit*k+1] = tempMeasurementCov[1][1];
00525 materialEffectsCov[nMeasPerHit*k+1][nMeasPerHit*l+1] = tempMeasurementCov[1][1];
00526 }
00527
00528
00529
00530 paramMaterialEffectsCov = paramMaterialEffectsCov.similarity(allCurvatureChanges[k]);
00531
00532 }
00533 theMeasurementsCov += materialEffectsCov;
00534
00535 return true;
00536 }
00537
00538
00539
00540 bool ReferenceTrajectory::addMaterialEffectsBp(const std::vector<AlgebraicMatrix> &allJacobians,
00541 const std::vector<AlgebraicMatrix> &allProjections,
00542 const std::vector<AlgebraicSymMatrix> &allCurvatureChanges,
00543 const std::vector<AlgebraicSymMatrix> &allDeltaParameterCovs,
00544 const std::vector<AlgebraicMatrix> &allLocalToCurv)
00545 {
00546
00547 int offsetPar = theNumberOfPars;
00548 int offsetMeas = nMeasPerHit * allJacobians.size();
00549 int ierr = 0;
00550
00551 AlgebraicMatrix tempJacobian;
00552 AlgebraicMatrix MSprojection(2,5);
00553 MSprojection[0][1] = 1;
00554 MSprojection[1][2] = 1;
00555 AlgebraicSymMatrix tempMSCov;
00556 AlgebraicSymMatrix tempMSCovProj;
00557 AlgebraicMatrix tempMSJacProj;
00558
00559 for (unsigned int k = 1; k < allJacobians.size(); ++k) {
00560
00561 int kbp = k-1;
00562 tempJacobian = allJacobians[k] * allCurvatureChanges[k];
00563 tempMSCov = allDeltaParameterCovs[k-1].similarity(allLocalToCurv[k-1]);
00564 tempMSCovProj = tempMSCov.similarity(MSprojection);
00565 theMeasurementsCov[offsetMeas+nMeasPerHit*kbp ][offsetMeas+nMeasPerHit*kbp ] = tempMSCovProj[0][0];
00566 theMeasurementsCov[offsetMeas+nMeasPerHit*kbp+1][offsetMeas+nMeasPerHit*kbp+1]= tempMSCovProj[1][1];
00567 theDerivatives[offsetMeas+nMeasPerHit*kbp ][offsetPar+2*kbp ] = 1.0;
00568 theDerivatives[offsetMeas+nMeasPerHit*kbp+1][offsetPar+2*kbp+1] = 1.0 ;
00569
00570 tempMSJacProj = (allProjections[k] * ( tempJacobian * allLocalToCurv[k-1].inverse(ierr) )) * MSprojection.T();
00571 if (ierr) {
00572 edm::LogError("Alignment") << "@SUB=ReferenceTrajectory::addMaterialEffectsBp"
00573 << "Inversion 1 for break points failed: " << ierr;
00574 return false;
00575 }
00576 theDerivatives[nMeasPerHit*k ][offsetPar+2*kbp ] = tempMSJacProj[0][0];
00577 theDerivatives[nMeasPerHit*k ][offsetPar+2*kbp+1] = tempMSJacProj[0][1];
00578 theDerivatives[nMeasPerHit*k+1][offsetPar+2*kbp ] = tempMSJacProj[1][0];
00579 theDerivatives[nMeasPerHit*k+1][offsetPar+2*kbp+1] = tempMSJacProj[1][1];
00580
00581 for (unsigned int l = k+1; l < allJacobians.size(); ++l) {
00582
00583 int kbp = k-1;
00584 tempJacobian = allJacobians[l] * allCurvatureChanges[l] * tempJacobian;
00585 tempMSJacProj = (allProjections[l] * ( tempJacobian * allLocalToCurv[k-1].inverse(ierr) )) * MSprojection.T();
00586 if (ierr) {
00587 edm::LogError("Alignment") << "@SUB=ReferenceTrajectory::addMaterialEffectsBp"
00588 << "Inversion 2 for break points failed: " << ierr;
00589 return false;
00590 }
00591 theDerivatives[nMeasPerHit*l ][offsetPar+2*kbp ] = tempMSJacProj[0][0];
00592 theDerivatives[nMeasPerHit*l ][offsetPar+2*kbp+1] = tempMSJacProj[0][1];
00593 theDerivatives[nMeasPerHit*l+1][offsetPar+2*kbp ] = tempMSJacProj[1][0];
00594 theDerivatives[nMeasPerHit*l+1][offsetPar+2*kbp+1] = tempMSJacProj[1][1];
00595
00596 }
00597
00598 }
00599
00600 return true;
00601 }
00602
00603
00604
00605 bool ReferenceTrajectory::addMaterialEffectsBrl(const std::vector<AlgebraicMatrix> &allCurvlinJacobians,
00606 const std::vector<AlgebraicMatrix> &allProjections,
00607 const std::vector<AlgebraicSymMatrix> &allCurvatureChanges,
00608 const std::vector<AlgebraicSymMatrix> &allDeltaParameterCovs,
00609 const std::vector<AlgebraicMatrix> &allLocalToCurv,
00610 const GlobalTrajectoryParameters >p)
00611 {
00612
00613
00614
00615
00616
00617
00618
00619
00620 int offsetPar = theNumberOfPars;
00621 int offsetMeas = nMeasPerHit*allCurvlinJacobians.size();
00622 int ierr = 0;
00623
00624 GlobalVector p = gtp.momentum();
00625 double cosLambda = sqrt((p.x()*p.x()+p.y()*p.y())/(p.x()*p.x()+p.y()*p.y()+p.z()*p.z()));
00626
00627
00628 AlgebraicMatrix QbypToCurv(5,1);
00629 QbypToCurv[0][0] = 1.;
00630 AlgebraicMatrix AngleToCurv(5,2);
00631 AngleToCurv[1][1] = 1.;
00632 AngleToCurv[2][0] = 1./cosLambda;
00633 AlgebraicMatrix CurvToAngle(2,5);
00634 CurvToAngle[1][1] = 1.;
00635 CurvToAngle[0][2] = cosLambda;
00636 AlgebraicMatrix OffsetToCurv(5,2);
00637 OffsetToCurv[3][0] = 1.;
00638 OffsetToCurv[4][1] = 1.;
00639 AlgebraicMatrix CurvToOffset(2,5);
00640 CurvToOffset[0][3] = 1.;
00641 CurvToOffset[1][4] = 1.;
00642
00643
00644 AlgebraicMatrix TrajToQbyp(1,5);
00645 TrajToQbyp[0][0] = 1.;
00646 AlgebraicMatrix TrajToOff1(2,5);
00647 TrajToOff1[0][1] = 1.;
00648 TrajToOff1[1][2] = 1.;
00649 AlgebraicMatrix TrajToOff2(2,5);
00650 TrajToOff2[0][3] = 1.;
00651 TrajToOff2[1][4] = 1.;
00652
00653 AlgebraicMatrix JacOffsetToAngleC, JacQbypToAngleC;
00654 AlgebraicMatrix JacCurvToOffsetL, JacOffsetToOffsetL, JacAngleToOffsetL, JacQbypToOffsetL, JacOffsetToAngleL;
00655 AlgebraicMatrix JacCurvToOffsetN, JacOffsetToOffsetN, JacAngleToOffsetN, JacQbypToOffsetN, JacOffsetToAngleN;
00656
00657
00658
00659 JacCurvToOffsetN = CurvToOffset * allCurvlinJacobians[1];
00660 JacOffsetToOffsetN = JacCurvToOffsetN * OffsetToCurv;
00661 JacAngleToOffsetN = JacCurvToOffsetN * AngleToCurv;
00662 JacQbypToOffsetN = JacCurvToOffsetN * QbypToCurv;
00663 JacOffsetToAngleN = JacAngleToOffsetN.inverse(ierr);
00664 if (ierr) {
00665 edm::LogError("Alignment") << "@SUB=ReferenceTrajectory::addMaterialEffectsBrl"
00666 << "Inversion 1 for fine broken lines failed: " << ierr;
00667 return false;
00668 }
00669 JacOffsetToAngleC = -(JacOffsetToAngleN * JacOffsetToOffsetN);
00670 JacQbypToAngleC = -(JacOffsetToAngleN * JacQbypToOffsetN);
00671
00672 AlgebraicMatrix JacTrajToAngle = JacQbypToAngleC * TrajToQbyp + JacOffsetToAngleC * TrajToOff1 + JacOffsetToAngleN * TrajToOff2;
00673
00674 theInnerTrajectoryToCurvilinear = QbypToCurv * TrajToQbyp + AngleToCurv * JacTrajToAngle + OffsetToCurv * TrajToOff1;
00675 theInnerLocalToTrajectory = theInnerTrajectoryToCurvilinear.inverse(ierr) * allLocalToCurv[0];
00676 if (ierr) {
00677 edm::LogError("Alignment") << "@SUB=ReferenceTrajectory::addMaterialEffectsBrl"
00678 << "Inversion 2 for fine broken lines failed: " << ierr;
00679 return false;
00680 }
00681
00682 AlgebraicMatrix tempJacobian(allCurvatureChanges[0]);
00683 AlgebraicSymMatrix tempMSCov;
00684 AlgebraicSymMatrix tempMSCovProj;
00685 AlgebraicMatrix tempJacL, tempJacN;
00686 AlgebraicMatrix JacOffsetToMeas;
00687
00688
00689 for (unsigned int k = 0; k < allCurvlinJacobians.size(); ++k) {
00690
00691 JacOffsetToMeas = (allProjections[k] * allLocalToCurv[k].inverse(ierr) ) * OffsetToCurv;
00692 if (ierr) {
00693 edm::LogError("Alignment") << "@SUB=ReferenceTrajectory::addMaterialEffectsBrl"
00694 << "Inversion 3 for fine broken lines failed: " << ierr;
00695 return false;
00696 }
00697 theDerivatives[nMeasPerHit*k ][offsetPar+2*k ] = JacOffsetToMeas[0][0];
00698 theDerivatives[nMeasPerHit*k ][offsetPar+2*k+1] = JacOffsetToMeas[0][1];
00699 theDerivatives[nMeasPerHit*k+1][offsetPar+2*k ] = JacOffsetToMeas[1][0];
00700 theDerivatives[nMeasPerHit*k+1][offsetPar+2*k+1] = JacOffsetToMeas[1][1];
00701 }
00702
00703
00704 for (unsigned int k = 1; k < allCurvlinJacobians.size()-1; ++k) {
00705
00706 int iMsMeas = k-1;
00707 int l = k-1;
00708 int n = k+1;
00709
00710
00711 tempMSCov = allDeltaParameterCovs[k].similarity(allLocalToCurv[k]);
00712 tempMSCovProj = tempMSCov.similarity(CurvToAngle);
00713 theMeasurementsCov[offsetMeas+nMeasPerHit*iMsMeas ][offsetMeas+nMeasPerHit*iMsMeas ] = tempMSCovProj[1][1];
00714 theMeasurementsCov[offsetMeas+nMeasPerHit*iMsMeas+1][offsetMeas+nMeasPerHit*iMsMeas+1] = tempMSCovProj[0][0];
00715
00716
00717 tempJacL = allCurvlinJacobians[k] * tempJacobian;
00718 JacCurvToOffsetL = CurvToOffset * tempJacL.inverse(ierr);
00719
00720 if (ierr) {
00721 edm::LogError("Alignment") << "@SUB=ReferenceTrajectory::addMaterialEffectsBrl"
00722 << "Inversion 4 for fine broken lines failed: " << ierr;
00723 return false;
00724 }
00725 JacOffsetToOffsetL = JacCurvToOffsetL * OffsetToCurv;
00726 JacAngleToOffsetL = JacCurvToOffsetL * AngleToCurv;
00727 JacQbypToOffsetL = JacCurvToOffsetL * QbypToCurv;
00728 JacOffsetToAngleL =-JacAngleToOffsetL.inverse(ierr);
00729 if (ierr) {
00730 edm::LogError("Alignment") << "@SUB=ReferenceTrajectory::addMaterialEffectsBrl"
00731 << "Inversion 5 for fine broken lines failed: " << ierr;
00732 return false;
00733 }
00734 tempJacobian = tempJacobian * allCurvatureChanges[k];
00735 tempJacN = allCurvlinJacobians[n] * tempJacobian;
00736 JacCurvToOffsetN = CurvToOffset * tempJacN;
00737 JacOffsetToOffsetN = JacCurvToOffsetN * OffsetToCurv;
00738 JacAngleToOffsetN = JacCurvToOffsetN * AngleToCurv;
00739 JacQbypToOffsetN = JacCurvToOffsetN * QbypToCurv;
00740 JacOffsetToAngleN = JacAngleToOffsetN.inverse(ierr);
00741 if (ierr) {
00742 edm::LogError("Alignment") << "@SUB=ReferenceTrajectory::addMaterialEffectsBrl"
00743 << "Inversion 6 for fine broken lines failed: " << ierr;
00744 return false;
00745 }
00746 JacOffsetToAngleC = -(JacOffsetToAngleL * JacOffsetToOffsetL + JacOffsetToAngleN * JacOffsetToOffsetN);
00747 JacQbypToAngleC = -(JacOffsetToAngleL * JacQbypToOffsetL + JacOffsetToAngleN * JacQbypToOffsetN);
00748
00749
00750 theDerivatives[offsetMeas+nMeasPerHit*iMsMeas ][ 0] = JacQbypToAngleC[0][0];
00751 theDerivatives[offsetMeas+nMeasPerHit*iMsMeas+1][ 0] = JacQbypToAngleC[1][0];
00752
00753 theDerivatives[offsetMeas+nMeasPerHit*iMsMeas ][offsetPar+2*l ] = JacOffsetToAngleL[0][0];
00754 theDerivatives[offsetMeas+nMeasPerHit*iMsMeas ][offsetPar+2*l+1] = JacOffsetToAngleL[0][1];
00755 theDerivatives[offsetMeas+nMeasPerHit*iMsMeas+1][offsetPar+2*l ] = JacOffsetToAngleL[1][0];
00756 theDerivatives[offsetMeas+nMeasPerHit*iMsMeas+1][offsetPar+2*l+1] = JacOffsetToAngleL[1][1];
00757
00758 theDerivatives[offsetMeas+nMeasPerHit*iMsMeas ][offsetPar+2*k ] = JacOffsetToAngleC[0][0];
00759 theDerivatives[offsetMeas+nMeasPerHit*iMsMeas ][offsetPar+2*k+1] = JacOffsetToAngleC[0][1];
00760 theDerivatives[offsetMeas+nMeasPerHit*iMsMeas+1][offsetPar+2*k ] = JacOffsetToAngleC[1][0];
00761 theDerivatives[offsetMeas+nMeasPerHit*iMsMeas+1][offsetPar+2*k+1] = JacOffsetToAngleC[1][1];
00762
00763
00764 theDerivatives[offsetMeas+nMeasPerHit*iMsMeas ][offsetPar+2*n ] = JacOffsetToAngleN[0][0];
00765 theDerivatives[offsetMeas+nMeasPerHit*iMsMeas ][offsetPar+2*n+1] = JacOffsetToAngleN[0][1];
00766 theDerivatives[offsetMeas+nMeasPerHit*iMsMeas+1][offsetPar+2*n ] = JacOffsetToAngleN[1][0];
00767 theDerivatives[offsetMeas+nMeasPerHit*iMsMeas+1][offsetPar+2*n+1] = JacOffsetToAngleN[1][1];
00768
00769 }
00770
00771 return true;
00772 }
00773
00774 bool ReferenceTrajectory::addMaterialEffectsBrl(const std::vector<AlgebraicMatrix> &allProjections,
00775 const std::vector<AlgebraicSymMatrix> &allDeltaParameterCovs,
00776 const std::vector<AlgebraicMatrix> &allLocalToCurv,
00777 const std::vector<double> &allSteps,
00778 const GlobalTrajectoryParameters >p,
00779 const double minStep)
00780 {
00781
00782
00783
00784
00785
00786 int offsetPar = theNumberOfPars;
00787 int offsetMeas = nMeasPerHit*allSteps.size();
00788 int ierr = 0;
00789
00790 GlobalVector p = gtp.momentum();
00791 double cosLambda = sqrt((p.x()*p.x()+p.y()*p.y())/(p.x()*p.x()+p.y()*p.y()+p.z()*p.z()));
00792 double bFac = -gtp.magneticFieldInInverseGeV(gtp.position()).mag();
00793
00794
00795 double delta (1.0/allSteps[1]);
00796 theInnerTrajectoryToCurvilinear[0][0] = 1;
00797 theInnerTrajectoryToCurvilinear[1][2] = -delta;
00798 theInnerTrajectoryToCurvilinear[1][4] = delta;
00799 theInnerTrajectoryToCurvilinear[2][0] = -0.5*bFac/delta;
00800 theInnerTrajectoryToCurvilinear[2][1] = -delta/cosLambda;
00801 theInnerTrajectoryToCurvilinear[2][3] = delta/cosLambda;
00802 theInnerTrajectoryToCurvilinear[3][1] = 1;
00803 theInnerTrajectoryToCurvilinear[4][2] = 1;
00804 theInnerLocalToTrajectory = theInnerTrajectoryToCurvilinear.inverse(ierr) * allLocalToCurv[0];
00805 if (ierr) {
00806 edm::LogError("Alignment") << "@SUB=ReferenceTrajectory::addMaterialEffectsBrl"
00807 << "Inversion 1 for coarse broken lines failed: " << ierr;
00808 return false;
00809 }
00810
00811 AlgebraicMatrix CurvToAngle(2,5);
00812 CurvToAngle[1][1] = 1.;
00813 CurvToAngle[0][2] = cosLambda;
00814 AlgebraicMatrix OffsetToCurv(5,2);
00815 OffsetToCurv[3][0] = 1.;
00816 OffsetToCurv[4][1] = 1.;
00817
00818 AlgebraicSymMatrix tempMSCov;
00819 AlgebraicSymMatrix tempMSCovProj;
00820 AlgebraicMatrix JacOffsetToMeas;
00821
00822
00823 std::vector<unsigned int> first(allSteps.size());
00824 std::vector<unsigned int> last (allSteps.size());
00825 std::vector<unsigned int> plane(allSteps.size());
00826 std::vector<double> sPlane(allSteps.size());
00827 unsigned int nPlane = 0;
00828 double sTot = 0;
00829
00830 for (unsigned int k = 1; k < allSteps.size(); ++k) {
00831 sTot += allSteps[k];
00832 if (fabs(allSteps[k])>minStep) { nPlane += 1; first[nPlane] = k; }
00833 last[nPlane] = k;
00834 plane[k] = nPlane;
00835 sPlane[nPlane] += sTot;
00836 }
00837 if (nPlane < 2) return false;
00838
00839 theNumberOfVirtualPars = 2*(nPlane+1);
00840 theNumberOfVirtualMeas = 2*(nPlane-1);
00841 for (unsigned int k = 0; k <= nPlane; ++k) { sPlane[k] /= (double) (last[k]-first[k]+1); }
00842
00843
00844 sTot = 0;
00845 for (unsigned int k = 0; k < allSteps.size(); ++k) {
00846 sTot += allSteps[k];
00847
00848 JacOffsetToMeas = (allProjections[k] * allLocalToCurv[k].inverse(ierr) ) * OffsetToCurv;
00849 if (ierr) {
00850 edm::LogError("Alignment") << "@SUB=ReferenceTrajectory::addMaterialEffectsBrl"
00851 << "Inversion 2 for coarse broken lines failed: " << ierr;
00852 return false;
00853 }
00854
00855 unsigned int iPlane = plane[k];
00856 if (last[iPlane] == first[iPlane])
00857 {
00858 theDerivatives[nMeasPerHit*k ][offsetPar+2*iPlane ] = JacOffsetToMeas[0][0];
00859 theDerivatives[nMeasPerHit*k ][offsetPar+2*iPlane+1] = JacOffsetToMeas[0][1];
00860 theDerivatives[nMeasPerHit*k+1][offsetPar+2*iPlane ] = JacOffsetToMeas[1][0];
00861 theDerivatives[nMeasPerHit*k+1][offsetPar+2*iPlane+1] = JacOffsetToMeas[1][1];
00862 } else
00863 {
00864 unsigned int jPlane;
00865 if (fabs(sTot) < fabs(sPlane[iPlane])) { jPlane = (iPlane>0) ? iPlane - 1 : 1; }
00866 else { jPlane = (iPlane<nPlane) ? iPlane + 1 : nPlane -1 ;}
00867
00868 double sDiff = sPlane[iPlane] - sPlane[jPlane];
00869 double iFrac = (sTot - sPlane[jPlane]) / sDiff;
00870 double jFrac = 1.0 - iFrac;
00871 theDerivatives[nMeasPerHit*k ][offsetPar+2*iPlane ] = JacOffsetToMeas[0][0]*iFrac;
00872 theDerivatives[nMeasPerHit*k ][offsetPar+2*iPlane+1] = JacOffsetToMeas[0][1]*iFrac;
00873 theDerivatives[nMeasPerHit*k+1][offsetPar+2*iPlane ] = JacOffsetToMeas[1][0]*iFrac;
00874 theDerivatives[nMeasPerHit*k+1][offsetPar+2*iPlane+1] = JacOffsetToMeas[1][1]*iFrac;
00875 theDerivatives[nMeasPerHit*k ][offsetPar+2*jPlane ] = JacOffsetToMeas[0][0]*jFrac;
00876 theDerivatives[nMeasPerHit*k ][offsetPar+2*jPlane+1] = JacOffsetToMeas[0][1]*jFrac;
00877 theDerivatives[nMeasPerHit*k+1][offsetPar+2*jPlane ] = JacOffsetToMeas[1][0]*jFrac;
00878 theDerivatives[nMeasPerHit*k+1][offsetPar+2*jPlane+1] = JacOffsetToMeas[1][1]*jFrac;
00879
00880
00881 }
00882 }
00883
00884
00885 for (unsigned int i = 1; i < nPlane; ++i) {
00886
00887 int iMsMeas = i-1;
00888 int l = i-1;
00889 int n = i+1;
00890
00891
00892 for (unsigned int k = first[i]; k <= last[i]; ++k) {
00893 tempMSCov = allDeltaParameterCovs[k].similarity(allLocalToCurv[k]);
00894 tempMSCovProj = tempMSCov.similarity(CurvToAngle);
00895 theMeasurementsCov[offsetMeas+nMeasPerHit*iMsMeas ][offsetMeas+nMeasPerHit*iMsMeas ] += tempMSCovProj[0][0];
00896 theMeasurementsCov[offsetMeas+nMeasPerHit*iMsMeas+1][offsetMeas+nMeasPerHit*iMsMeas+1] += tempMSCovProj[1][1];
00897 }
00898
00899 double stepK = sPlane[i] - sPlane[l];
00900 double stepN = sPlane[n] - sPlane[i];
00901 double deltaK (1.0/stepK);
00902 double deltaN (1.0/stepN);
00903
00904 theDerivatives[offsetMeas+nMeasPerHit*iMsMeas ][ 0] = -0.5*bFac*(stepK+stepN)*cosLambda;
00905
00906 theDerivatives[offsetMeas+nMeasPerHit*iMsMeas ][offsetPar+2*l ] = deltaK;
00907 theDerivatives[offsetMeas+nMeasPerHit*iMsMeas+1][offsetPar+2*l+1] = deltaK;
00908
00909 theDerivatives[offsetMeas+nMeasPerHit*iMsMeas ][offsetPar+2*i ] = -(deltaK + deltaN);
00910 theDerivatives[offsetMeas+nMeasPerHit*iMsMeas+1][offsetPar+2*i+1] = -(deltaK + deltaN);
00911
00912 theDerivatives[offsetMeas+nMeasPerHit*iMsMeas ][offsetPar+2*n ] = deltaN;
00913 theDerivatives[offsetMeas+nMeasPerHit*iMsMeas+1][offsetPar+2*n+1] = deltaN;
00914 }
00915
00916 return true;
00917 }