Go to the documentation of this file.00001 #include "TrackingTools/GsfTools/interface/BasicMultiTrajectoryState.h"
00002
00003 #include "DataFormats/GeometrySurface/interface/BoundPlane.h"
00004 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00005
00006 using namespace SurfaceSideDefinition;
00007
00008 BasicMultiTrajectoryState::BasicMultiTrajectoryState( const std::vector<TSOS>& tsvec) :
00009 BasicTrajectoryState(tsvec.front().surface()) {
00010
00011
00012 const BoundPlane* bp = dynamic_cast<const BoundPlane*>(&tsvec.begin()->surface());
00013 if unlikely( bp==0 )
00014 throw cms::Exception("LogicError") << "MultiTrajectoryState constructed on cylinder";
00015
00016 theStates.reserve(tsvec.size());
00017 for (std::vector<TSOS>::const_iterator i=tsvec.begin(); i!=tsvec.end(); i++) {
00018 if unlikely(!i->isValid()) {
00019 throw cms::Exception("LogicError") << "MultiTrajectoryState constructed with invalid state";
00020 }
00021 if unlikely(i->hasError() != tsvec.front().hasError()) {
00022 throw cms::Exception("LogicError") << "MultiTrajectoryState mixes states with and without errors";
00023 }
00024 if unlikely( &i->surface() != &tsvec.front().surface()) {
00025 throw cms::Exception("LogicError") << "MultiTrajectoryState mixes states with different surfaces";
00026 }
00027 if unlikely( i->surfaceSide() != tsvec.front().surfaceSide()) {
00028 throw cms::Exception("LogicError")
00029 << "MultiTrajectoryState mixes states defined before and after material";
00030 }
00031 if unlikely( i->localParameters().pzSign()*tsvec.front().localParameters().pzSign()<0. ) {
00032 throw cms::Exception("LogicError")
00033 << "MultiTrajectoryState mixes states with different signs of local p_z";
00034 }
00035
00036 theStates.push_back( *i);
00037 }
00038 combine();
00039 }
00040
00041
00042
00043 void BasicMultiTrajectoryState::rescaleError(double factor) {
00044
00045 if unlikely(theStates.empty()) {
00046 edm::LogError("BasicMultiTrajectoryState") << "Trying to rescale errors of empty MultiTrajectoryState!";
00047 return;
00048 }
00049
00050 for (std::vector<TSOS>::iterator it = theStates.begin(); it != theStates.end(); it++) {
00051 it->rescaleError(factor);
00052 }
00053 combine();
00054 }
00055
00056 void
00057 BasicMultiTrajectoryState::combine() {
00058 const std::vector<TrajectoryStateOnSurface>& tsos = theStates;
00059
00060 if unlikely(tsos.empty()) {
00061 edm::LogError("MultiTrajectoryStateCombiner")
00062 << "Trying to collapse empty set of trajectory states!";
00063 return;
00064 }
00065
00066 double pzSign = tsos.front().localParameters().pzSign();
00067 for (std::vector<TrajectoryStateOnSurface>::const_iterator it = tsos.begin();
00068 it != tsos.end(); it++) {
00069 if unlikely(it->localParameters().pzSign() != pzSign) {
00070 edm::LogError("MultiTrajectoryStateCombiner")
00071 << "Trying to collapse trajectory states with different signs on p_z!";
00072 return;
00073 }
00074 }
00075
00076 if unlikely(tsos.size() == 1) {
00077 BasicTrajectoryState::update(tsos.front().localParameters(),
00078 tsos.front().localError(),
00079 tsos.front().surface(),
00080 tsos.front().magneticField(),
00081 tsos.front().surfaceSide(),
00082 tsos.front().weight()
00083 );
00084 return;
00085 }
00086
00087 double sumw = 0.;
00088
00089 AlgebraicVector5 mean;
00090 AlgebraicSymMatrix55 covarPart1, covarPart2;
00091 for (std::vector<TrajectoryStateOnSurface>::const_iterator it1 = tsos.begin();
00092 it1 != tsos.end(); it1++) {
00093 double weight = it1->weight();
00094 AlgebraicVector5 param = it1->localParameters().vector();
00095 sumw += weight;
00096 mean += weight * param;
00097 covarPart1 += weight * it1->localError().matrix();
00098 for (std::vector<TrajectoryStateOnSurface>::const_iterator it2 = it1 + 1;
00099 it2 != tsos.end(); it2++) {
00100 AlgebraicVector5 diff = param - it2->localParameters().vector();
00101 AlgebraicSymMatrix11 s = AlgebraicMatrixID();
00102 covarPart2 += (weight * it2->weight()) *
00103 ROOT::Math::Similarity(AlgebraicMatrix51(diff.Array(), 5), s);
00104
00105 }
00106 }
00107 double sumwI = 1.0/sumw;
00108 mean *= sumwI;
00109 covarPart1 *= sumwI; covarPart2 *= (sumwI*sumwI);
00110 AlgebraicSymMatrix55 covar = covarPart1 + covarPart2;
00111
00112 BasicTrajectoryState::update(LocalTrajectoryParameters(mean, pzSign),
00113 LocalTrajectoryError(covar),
00114 tsos.front().surface(),
00115 tsos.front().magneticField(),
00116 tsos.front().surfaceSide(),
00117 sumw);
00118 }
00119
00120 void
00121 BasicMultiTrajectoryState::
00122 update( const LocalTrajectoryParameters& p,
00123 const Surface& aSurface,
00124 const MagneticField* field,
00125 const SurfaceSide side)
00126 {
00127 throw cms::Exception("LogicError",
00128 "BasicMultiTrajectoryState::update(LocalTrajectoryParameters, Surface, ...) called even if canUpdateLocalParameters() is false");
00129 }
00130
00131 void
00132 BasicMultiTrajectoryState::
00133 update( const LocalTrajectoryParameters& p,
00134 const LocalTrajectoryError& err,
00135 const Surface& aSurface,
00136 const MagneticField* field,
00137 const SurfaceSide side,
00138 double weight)
00139 {
00140 throw cms::Exception("LogicError",
00141 "BasicMultiTrajectoryState::update(LocalTrajectoryParameters, LocalTrajectoryError, ...) called even if canUpdateLocalParameters() is false");
00142 }
00143