CMS 3D CMS Logo

/data/doxygen/doxygen-1.7.3/gen/CMSSW_4_2_8/src/TrackingTools/GsfTools/src/BasicMultiTrajectoryState.cc

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   theCombinedStateUp2Date( false)
00010 {
00011   theStates.reserve(tsvec.size());
00012   for (std::vector<TSOS>::const_iterator i=tsvec.begin(); i!=tsvec.end(); i++) {
00013     if (!i->isValid()) {
00014       throw cms::Exception("LogicError") << "MultiTrajectoryState constructed with invalid state";
00015     }
00016     if (i->hasError() != tsvec.front().hasError()) {
00017       throw cms::Exception("LogicError") << "MultiTrajectoryState mixes states with and without errors";
00018     }
00019     if ( &i->surface() != &tsvec.front().surface()) {
00020       throw cms::Exception("LogicError") << "MultiTrajectoryState mixes states with different surfaces";
00021     }
00022     if ( i->surfaceSide() != tsvec.front().surfaceSide()) {
00023       throw cms::Exception("LogicError") 
00024         << "MultiTrajectoryState mixes states defined before and after material";
00025     }
00026     if ( i->localParameters().pzSign()*tsvec.front().localParameters().pzSign()<0. ) {
00027       throw cms::Exception("LogicError") 
00028         << "MultiTrajectoryState mixes states with different signs of local p_z";
00029     }
00030     if ( i==tsvec.begin() ) {
00031       // only accept planes!!
00032       const BoundPlane* bp = dynamic_cast<const BoundPlane*>(&i->surface());
00033       if ( bp==0 )
00034         throw cms::Exception("LogicError") << "MultiTrajectoryState constructed on cylinder";
00035     }
00036     theStates.push_back( *i);
00037   }
00038 }
00039 
00040 void BasicMultiTrajectoryState::checkCombinedState() const
00041 {
00042   if (theCombinedStateUp2Date) return;
00043 
00044   theCombinedState = theCombiner.combine( theStates);
00045   theCombinedStateUp2Date = true;
00046   
00047 }
00048 
00049 double BasicMultiTrajectoryState::weight() const {
00050 
00051   if (theStates.empty()) {
00052     edm::LogError("BasicMultiTrajectoryState") 
00053       << "Asking for weight of empty MultiTrajectoryState, returning zero!";
00054     return 0.;
00055   }
00056 
00057   double sumw = 0.;
00058   for (std::vector<TSOS>::const_iterator it = theStates.begin(); it != theStates.end(); it++) {
00059     sumw += it->weight();
00060   }
00061   return sumw;
00062 }
00063 
00064 
00065 void BasicMultiTrajectoryState::rescaleError(double factor) {
00066 
00067   if (theStates.empty()) {
00068     edm::LogError("BasicMultiTrajectoryState") << "Trying to rescale errors of empty MultiTrajectoryState!";
00069     return;
00070   }
00071   
00072   for (std::vector<TSOS>::iterator it = theStates.begin(); it != theStates.end(); it++) {
00073     it->rescaleError(factor);
00074   }
00075   theCombinedStateUp2Date = false;
00076 }
00077 
00078 const MagneticField*
00079 BasicMultiTrajectoryState::magneticField () const
00080 {
00081   //
00082   // Magnetic field should be identical in all components:
00083   // avoid forcing the combination of states and take value from 1st component!
00084   //
00085   if (theStates.empty()) {
00086     edm::LogError("BasicMultiTrajectoryState") 
00087       << "Asking for magneticField of empty MultiTrajectoryState, returning null pointer!";
00088     return 0;
00089   }
00090   return theStates.front().magneticField();
00091 }
00092 
00093 SurfaceSide
00094 BasicMultiTrajectoryState::surfaceSide () const
00095 {
00096   //
00097   // SurfaceSide should be identical in all components:
00098   // avoid forcing the combination of states and take value from 1st component!
00099   //
00100   if (theStates.empty()) {
00101     edm::LogError("BasicMultiTrajectoryState") 
00102       << "Asking for magneticField of empty MultiTrajectoryState, returning atCenterOfSurface!";
00103     return atCenterOfSurface;
00104   }
00105   return theStates.front().surfaceSide();
00106 }
00107 
00108 void
00109 BasicMultiTrajectoryState::
00110 update( const LocalTrajectoryParameters& p,
00111         const Surface& aSurface,
00112         const MagneticField* field,
00113         const SurfaceSide side) 
00114 {
00115   throw cms::Exception("LogicError", 
00116                        "BasicMultiTrajectoryState::update(LocalTrajectoryParameters, Surface, ...) called even if canUpdateLocalParameters() is false");
00117 }
00118 
00119 void
00120 BasicMultiTrajectoryState::
00121 update( const LocalTrajectoryParameters& p,
00122         const LocalTrajectoryError& err,
00123         const Surface& aSurface,
00124         const MagneticField* field,
00125         const SurfaceSide side,
00126         double weight) 
00127 {
00128   throw cms::Exception("LogicError", 
00129                        "BasicMultiTrajectoryState::update(LocalTrajectoryParameters, LocalTrajectoryError, ...) called even if canUpdateLocalParameters() is false");
00130 }
00131