#include <RecoMuon/L3TrackFinder/interface/MuonCkfTrajectoryBuilder.h>
Public Member Functions | |
MuonCkfTrajectoryBuilder (const edm::ParameterSet &conf, const TrajectoryStateUpdator *updator, const Propagator *propagatorAlong, const Propagator *propagatorOpposite, const Propagator *propagatorProximity, const Chi2MeasurementEstimatorBase *estimator, const TransientTrackingRecHitBuilder *RecHitBuilder, const MeasurementTracker *measurementTracker, const TrajectoryFilter *filter) | |
virtual | ~MuonCkfTrajectoryBuilder () |
Protected Member Functions | |
void | collectMeasurement (const DetLayer *layer, const std::vector< const DetLayer * > &nl, const TrajectoryStateOnSurface ¤tState, std::vector< TM > &result, int &invalidHits, const Propagator *) const |
virtual void | findCompatibleMeasurements (const TempTrajectory &traj, std::vector< TrajectoryMeasurement > &result) const |
Protected Attributes | |
const Propagator * | theProximityPropagator |
double | theRescaleErrorIfFail |
bool | theUseSeedLayer |
Definition at line 6 of file MuonCkfTrajectoryBuilder.h.
MuonCkfTrajectoryBuilder::MuonCkfTrajectoryBuilder | ( | const edm::ParameterSet & | conf, | |
const TrajectoryStateUpdator * | updator, | |||
const Propagator * | propagatorAlong, | |||
const Propagator * | propagatorOpposite, | |||
const Propagator * | propagatorProximity, | |||
const Chi2MeasurementEstimatorBase * | estimator, | |||
const TransientTrackingRecHitBuilder * | RecHitBuilder, | |||
const MeasurementTracker * | measurementTracker, | |||
const TrajectoryFilter * | filter | |||
) |
Definition at line 16 of file MuonCkfTrajectoryBuilder.cc.
References edm::ParameterSet::getParameter(), theRescaleErrorIfFail, and theUseSeedLayer.
00024 : 00025 CkfTrajectoryBuilder(conf,updator,propagatorAlong,propagatorOpposite,estimator,RecHitBuilder,measurementTracker,filter), 00026 theProximityPropagator(propagatorProximity) 00027 { 00028 //and something specific to me ? 00029 theUseSeedLayer = conf.getParameter<bool>("useSeedLayer"); 00030 theRescaleErrorIfFail = conf.getParameter<double>("rescaleErrorIfFail"); 00031 }
MuonCkfTrajectoryBuilder::~MuonCkfTrajectoryBuilder | ( | ) | [virtual] |
void MuonCkfTrajectoryBuilder::collectMeasurement | ( | const DetLayer * | layer, | |
const std::vector< const DetLayer * > & | nl, | |||
const TrajectoryStateOnSurface & | currentState, | |||
std::vector< TM > & | result, | |||
int & | invalidHits, | |||
const Propagator * | prop | |||
) | const [protected] |
Definition at line 66 of file MuonCkfTrajectoryBuilder.cc.
References dumpMeasurements(), TransverseImpactPointExtrapolator::extrapolate(), TrajectoryStateOnSurface::globalMomentum(), TrajectoryStateOnSurface::globalPosition(), TrajectoryStateOnSurface::isValid(), LogDebug, LayerMeasurements::measurements(), middle, BaseCkfTrajectoryBuilder::theEstimator, BaseCkfTrajectoryBuilder::theLayerMeasurements, and tmp.
Referenced by findCompatibleMeasurements().
00066 { 00067 for (std::vector<const DetLayer*>::const_iterator il = nl.begin(); 00068 il != nl.end(); il++) { 00069 00070 TSOS stateToUse = currentState; 00071 00072 if (layer == (*il)){ 00073 LogDebug("CkfPattern")<<" self propagating in findCompatibleMeasurements.\n from: \n"<<stateToUse; 00074 //self navigation case 00075 // go to a middle point first 00076 TransverseImpactPointExtrapolator middle; 00077 GlobalPoint center(0,0,0); 00078 stateToUse = middle.extrapolate(stateToUse, center, *prop); 00079 00080 if (!stateToUse.isValid()) continue; 00081 LogDebug("CkfPattern")<<"to: "<<stateToUse; 00082 } 00083 00084 std::vector<TM> tmp = 00085 theLayerMeasurements->measurements((**il),currentState, *prop, *theEstimator); 00086 00087 LogDebug("CkfPattern")<<tmp.size()<<" measurements returned by LayerMeasurements"; 00088 00089 if ( !tmp.empty()) { 00090 // FIXME durty-durty-durty cleaning: never do that please ! 00091 /* for (vector<TM>::iterator it = tmp.begin(); it!=tmp.end(); ++it) 00092 {if (it->recHit()->det()==0) it=tmp.erase(it)--;}*/ 00093 00094 if ( result.empty()) result = tmp; 00095 else { 00096 // keep one dummy TM at the end, skip the others 00097 result.insert( result.end()-invalidHits, tmp.begin(), tmp.end()); 00098 } 00099 invalidHits++; 00100 } 00101 } 00102 00103 LogDebug("CkfPattern")<<"starting from:\n" 00104 <<"x: "<<currentState.globalPosition()<<"\n" 00105 <<"p: "<<currentState.globalMomentum()<<"\n" 00106 <<dumpMeasurements(result); 00107 }
void MuonCkfTrajectoryBuilder::findCompatibleMeasurements | ( | const TempTrajectory & | traj, | |
std::vector< TrajectoryMeasurement > & | result | |||
) | const [protected, virtual] |
Reimplemented from CkfTrajectoryBuilder.
Definition at line 112 of file MuonCkfTrajectoryBuilder.cc.
References alongMomentum, collectMeasurement(), PTrajectoryStateOnDet::detId(), GeometricSearchTracker::detLayer(), TempTrajectory::direction(), TempTrajectory::empty(), TrajectoryStateOnSurface::freeState(), g, MeasurementTracker::geometricSearchTracker(), MeasurementTracker::geomTracker(), i, id, TrackingGeometry::idToDet(), insideOut, edm::es::l(), TempTrajectory::lastLayer(), TempTrajectory::lastMeasurement(), LogDebug, Propagator::magneticField(), DetLayer::nextLayers(), outsideIn, TrajectoryMeasurement::recHit(), TrajectoryStateOnSurface::rescaleError(), TempTrajectory::seed(), python::multivaluedict::sort(), TrajectorySeed::startingState(), GeomDet::surface(), BaseCkfTrajectoryBuilder::theForwardPropagator, BaseCkfTrajectoryBuilder::theMeasurementTracker, theProximityPropagator, theRescaleErrorIfFail, theUseSeedLayer, TrajectoryStateTransform::transientState(), and TrajectoryMeasurement::updatedState().
00114 { 00115 int invalidHits = 0; 00116 00117 00118 std::vector<const DetLayer*> nl; 00119 00120 if (traj.empty()) 00121 { 00122 LogDebug("CkfPattern")<<"using JR patch for no measurement case"; 00123 //what if there are no measurement on the Trajectory 00124 00125 //set the currentState to be the one from the trajectory seed starting point 00126 PTrajectoryStateOnDet ptod =traj.seed().startingState(); 00127 DetId id(ptod.detId()); 00128 const GeomDet * g = theMeasurementTracker->geomTracker()->idToDet(id); 00129 const Surface * surface=&g->surface(); 00130 TrajectoryStateTransform tsTransform; 00131 TrajectoryStateOnSurface currentState(tsTransform.transientState(ptod,surface,theForwardPropagator->magneticField())); 00132 00133 //set the next layers to be that one the state is on 00134 const DetLayer * l=theMeasurementTracker->geometricSearchTracker()->detLayer(id); 00135 00136 if (theUseSeedLayer){ 00137 { 00138 //get the measurements on the layer first 00139 LogDebug("CkfPattern")<<"using the layer of the seed first."; 00140 nl.push_back(l); 00141 collectMeasurement(l,nl,currentState,result,invalidHits,theProximityPropagator); 00142 } 00143 00144 //if fails: try to rescale locally the state to find measurements 00145 if ((uint)invalidHits==result.size() && theRescaleErrorIfFail!=1.0 && result.size()!=0) 00146 { 00147 result.clear(); 00148 LogDebug("CkfPattern")<<"using a rescale by "<< theRescaleErrorIfFail <<" to find measurements."; 00149 TrajectoryStateOnSurface rescaledCurrentState = currentState; 00150 rescaledCurrentState.rescaleError(theRescaleErrorIfFail); 00151 invalidHits=0; 00152 collectMeasurement(l,nl,rescaledCurrentState,result,invalidHits,theProximityPropagator); 00153 } 00154 } 00155 00156 //if fails: go to next layers 00157 if (result.size()==0 || (uint)invalidHits==result.size()) 00158 { 00159 result.clear(); 00160 LogDebug("CkfPattern")<<"Need to go to next layer to get measurements"; 00161 //the following will "JUMP" the first layer measurements 00162 nl = l->nextLayers(*currentState.freeState(), traj.direction()); 00163 if (nl.size()==0){ 00164 LogDebug("CkfPattern")<<" there was no next layer with wellInside. Use the next with no check."; 00165 //means you did not get any compatible layer on the next 1/2 tracker layer. 00166 // use the next layers with no checking 00167 nl = l->nextLayers(((traj.direction()==alongMomentum)?insideOut:outsideIn)); 00168 } 00169 invalidHits=0; 00170 collectMeasurement(l,nl,currentState,result,invalidHits,theForwardPropagator); 00171 } 00172 00173 //if fails: this is on the next layers already, try rescaling locally the state 00174 if (result.size()!=0 && (uint)invalidHits==result.size() && theRescaleErrorIfFail!=1.0) 00175 { 00176 result.clear(); 00177 LogDebug("CkfPattern")<<"using a rescale by "<< theRescaleErrorIfFail <<" to find measurements on next layers."; 00178 TrajectoryStateOnSurface rescaledCurrentState = currentState; 00179 rescaledCurrentState.rescaleError(theRescaleErrorIfFail); 00180 invalidHits=0; 00181 collectMeasurement(l,nl,rescaledCurrentState, result,invalidHits,theForwardPropagator); 00182 } 00183 00184 } 00185 else //regular case 00186 { 00187 00188 TSOS currentState( traj.lastMeasurement().updatedState()); 00189 00190 nl = traj.lastLayer()->nextLayers( *currentState.freeState(), traj.direction()); 00191 if (nl.empty()){LogDebug("CkfPattern")<<" no next layers... going "<<traj.direction()<<"\n from: \n"<<currentState<<"\n from detId: "<<traj.lastMeasurement().recHit()->geographicalId().rawId(); return ;} 00192 00193 collectMeasurement(traj.lastLayer(),nl,currentState,result,invalidHits,theForwardPropagator); 00194 } 00195 00196 00197 // sort the final result, keep dummy measurements at the end 00198 if ( result.size() > 1) { 00199 sort( result.begin(), result.end()-invalidHits, TrajMeasLessEstim()); 00200 } 00201 00202 #ifdef DEBUG_INVALID 00203 bool afterInvalid = false; 00204 for (std::vector<TM>::const_iterator i=result.begin(); 00205 i!=result.end(); i++) { 00206 if ( ! i->recHit().isValid()) afterInvalid = true; 00207 if (afterInvalid && i->recHit().isValid()) { 00208 edm::LogError("CkfPattern") << "CkfTrajectoryBuilder error: valid hit after invalid!" ; 00209 } 00210 } 00211 #endif 00212 00213 //analyseMeasurements( result, traj); 00214 00215 }
const Propagator* MuonCkfTrajectoryBuilder::theProximityPropagator [mutable, protected] |
Definition at line 27 of file MuonCkfTrajectoryBuilder.h.
Referenced by findCompatibleMeasurements().
double MuonCkfTrajectoryBuilder::theRescaleErrorIfFail [protected] |
Definition at line 26 of file MuonCkfTrajectoryBuilder.h.
Referenced by findCompatibleMeasurements(), and MuonCkfTrajectoryBuilder().
bool MuonCkfTrajectoryBuilder::theUseSeedLayer [protected] |
Definition at line 25 of file MuonCkfTrajectoryBuilder.h.
Referenced by findCompatibleMeasurements(), and MuonCkfTrajectoryBuilder().