00001 #include "RecoTracker/TrackProducer/interface/MTFTrackProducerAlgorithm.h"
00002 #include "RecoTracker/SiTrackerMRHTools/interface/SiTrackerMultiRecHitUpdatorMTF.h"
00003 #include "RecoTracker/SiTrackerMRHTools/interface/MeasurementByLayerGrouper.h"
00004 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00005 #include "RecoTracker/SiTrackerMRHTools/interface/MultiTrackFilterHitCollector.h"
00006 #include "MagneticField/Engine/interface/MagneticField.h"
00007 #include "Geometry/CommonDetUnit/interface/TrackingGeometry.h"
00008
00009 #include "DataFormats/TrackCandidate/interface/TrackCandidate.h"
00010 #include "DataFormats/TrackingRecHit/interface/TrackingRecHitFwd.h"
00011 #include "DataFormats/TrackReco/interface/Track.h"
00012
00013 #include "TrackingTools/TrackFitters/interface/TrajectoryFitter.h"
00014 #include "TrackingTools/TrajectoryState/interface/TrajectoryStateOnSurface.h"
00015 #include "TrackingTools/TrajectoryState/interface/TrajectoryStateTransform.h"
00016 #include "TrackingTools/TransientTrackingRecHit/interface/TransientTrackingRecHitBuilder.h"
00017 #include "TrackingTools/TransientTrackingRecHit/interface/InvalidTransientRecHit.h"
00018 #include "TrackingTools/PatternTools/interface/TSCBLBuilderNoMaterial.h"
00019 #include "TrackingTools/PatternTools/interface/TransverseImpactPointExtrapolator.h"
00020 #include "TrackingTools/TrackFitters/interface/TrajectoryStateWithArbitraryError.h"
00021 #include "TrackingTools/TrackFitters/interface/TrajectoryStateCombiner.h"
00022 #include "Utilities/General/interface/CMSexception.h"
00023 #include "RecoTracker/TransientTrackingRecHit/interface/TSiTrackerMultiRecHit.h"
00024
00025 void MTFTrackProducerAlgorithm::runWithCandidate(const TrackingGeometry * theG,
00026 const MagneticField * theMF,
00027
00028 const std::vector<Trajectory>& theTrajectoryCollection,
00029
00030 const TrajectoryFitter * theFitter,
00031 const TransientTrackingRecHitBuilder* builder,
00032 const MultiTrackFilterHitCollector* measurementCollector,
00033 const SiTrackerMultiRecHitUpdatorMTF* updator,
00034 const reco::BeamSpot& bs,
00035 AlgoProductCollection& algoResults) const
00036 {
00037 edm::LogInfo("MTFTrackProducer") << "Number of Trajectories: " << theTrajectoryCollection.size() << "\n";
00038
00039 int cont = 0;
00040 float ndof = 0;
00041
00042 std::vector<Trajectory> mvtraj=theTrajectoryCollection;
00043
00044
00045
00046
00047 std::map<int, std::vector<TrajectoryMeasurement> > mvtm;
00048 int a=0;
00049
00050 for(std::vector<Trajectory>::iterator imvtraj=mvtraj.begin(); imvtraj!=mvtraj.end(); imvtraj++)
00051
00052 {
00053
00054 Trajectory *traj = &(*imvtraj);
00055 mvtm.insert( make_pair(a, traj->measurements()) );
00056 a++;
00057 }
00058
00059 LogDebug("MTFTrackProducerAlgorithm") << "after the cicle found " << mvtraj.size() << " trajectories" << std::endl;
00060 LogDebug("MTFTrackProducerAlgorithm") << "built a map of " << mvtm.size() << " elements (trajectories, yeah)" << std::endl;
00061
00062
00063 std::vector<std::pair<TransientTrackingRecHit::RecHitContainer, TrajectoryStateOnSurface> > vecmrhits;
00064 std::vector<Trajectory> transientmvtraj;
00065 int b=0;
00066 for(std::vector<Trajectory>::const_iterator im=mvtraj.begin(); im!=mvtraj.end(); im++)
00067
00068 {
00069
00070 LogDebug("MTFTrackProducerAlgorithm") << "about to collect hits for the trajectory number " << b << std::endl;
00071
00072
00073 std::pair<TransientTrackingRecHit::RecHitContainer, TrajectoryStateOnSurface> hits = collectHits(mvtm, measurementCollector, b);
00074 vecmrhits.push_back(hits);
00075
00076 LogDebug("MTFTrackProducerAlgorithm") << "hits collected";
00077
00078
00079 std::vector<Trajectory> vtraj(1,(*im));
00080
00081
00082 bool fitresult = fit(hits, theFitter, vtraj);
00083 LogDebug("MTFTrackProducerAlgorithm") << "fit done";
00084
00085
00086 if(fitresult == true)
00087 {
00088 LogDebug("MTFTrackProducerAlgorithm") << "fit was good for trajectory number" << b << "\n"
00089 << "the trajectory has size" << vtraj.size() << "\n"
00090 << "and its number of measurements is: " << vtraj.front().measurements().size() << "\n";
00091 Trajectory thetraj = vtraj.front();
00092
00093
00094 transientmvtraj.push_back(thetraj);
00095 }
00096
00097
00098
00099
00100
00101
00102
00103 b++;
00104 }
00105
00106 LogDebug("MTFTrackProducerAlgorithm") << "cicle finished";
00107
00108
00109 mvtraj.swap(transientmvtraj);
00110 LogDebug("MTFTrackProducerAlgorithm") << " with " << mvtraj.size() << "trajectories\n";
00111
00112
00113
00114
00115
00116 for(std::vector<double>::const_iterator ian = updator->getAnnealingProgram().begin(); ian != updator->getAnnealingProgram().end(); ian++)
00117
00118 {
00119
00120
00121 std::map<int, std::vector<TrajectoryMeasurement> > transientmvtm1;
00122 int b=0;
00123
00124 for(std::vector<Trajectory>::iterator imv=mvtraj.begin(); imv!=mvtraj.end(); imv++)
00125 {
00126 Trajectory *traj = &(*imv);
00127 transientmvtm1.insert( make_pair(b,traj->measurements()) );
00128 b++;
00129 }
00130
00131
00132 mvtm.swap(transientmvtm1);
00133
00134
00135 vecmrhits.clear();
00136
00137
00138 for (unsigned int d=0; d<mvtraj.size(); d++)
00139 {
00140
00141 std::pair<TransientTrackingRecHit::RecHitContainer, TrajectoryStateOnSurface>
00142 curiterationhits = updateHits(mvtm, measurementCollector, updator, *ian, builder, d);
00143
00144 vecmrhits.push_back(curiterationhits);
00145
00146 }
00147
00148
00149
00150 LogDebug("MTFTrackProducerAlgorithm") << "vector vecmrhits has size "
00151 << vecmrhits.size() << std::endl;
00152
00153
00154
00155
00156
00157
00158 std::vector<Trajectory> transientmvtrajone;
00159
00160 int n=0;
00161
00162 for(std::vector<Trajectory>::const_iterator j=mvtraj.begin(); j!=mvtraj.end(); j++)
00163 {
00164
00165
00166 std::vector<Trajectory> vtraj(1,(*j));
00167
00168 if(vtraj.size()){
00169
00170 LogDebug("MTFTrackProducerAlgorithm") << "Seed direction is " << vtraj.front().seed().direction()
00171 << "Traj direction is " << vtraj.front().direction();
00172
00173
00174
00175
00176
00177
00178
00179
00180 fit(vecmrhits[n], theFitter, vtraj);
00181
00182
00183
00184
00185
00186 }
00187
00188 else
00189
00190 {
00191 LogDebug("MTFTrackProducerAlgorithm") << "in map making skipping trajectory number: " << n <<"\n";
00192 continue;
00193 }
00194
00195
00196 Trajectory thetraj = vtraj.front();
00197
00198
00199
00200
00201
00202
00203
00204 transientmvtrajone.push_back(thetraj);
00205
00206
00207
00208
00209
00210
00211
00212
00213 n++;
00214 }
00215
00216
00217 mvtraj.swap(transientmvtrajone);
00218
00219
00220
00221 LogDebug("MTFTrackProducerAlgorithm") << "number of trajectories after annealing value "
00222 << (*ian)
00223 << " annealing step "
00224 << mvtraj.size() << std::endl;
00225 }
00226
00227
00228
00229
00230 LogDebug("MTFTrackProducerAlgorithm") << "Ended annealing program with " << mvtraj.size() << " trajectories" << std::endl;
00231
00232
00233
00234
00235
00236
00237 for(std::vector<Trajectory>::iterator it=mvtraj.begin(); it!=mvtraj.end();it++){
00238
00239 std::vector<Trajectory> vtraj(1,(*it));
00240
00241
00242
00243
00244
00245
00246
00247
00248
00249 ndof = calculateNdof(vtraj);
00250
00251 bool ok = buildTrack(vtraj, algoResults, ndof, bs);
00252 if(ok) cont++;
00253
00254 }
00255
00256 edm::LogInfo("TrackProducer") << "Number of Tracks found: " << cont << "\n";
00257
00258 }
00259
00260
00261
00262
00263
00264
00265 std::pair<TransientTrackingRecHit::RecHitContainer, TrajectoryStateOnSurface>
00266 MTFTrackProducerAlgorithm::collectHits(const std::map<int, std::vector<TrajectoryMeasurement> >& vtm,
00267 const MultiTrackFilterHitCollector* measurementCollector,
00268 int i) const{
00269
00270 TransientTrackingRecHit::RecHitContainer hits;
00271
00272
00273 std::vector<TrajectoryMeasurement> collectedmeas = measurementCollector->recHits(vtm,i,1.);
00274 LogDebug("MTFTrackProducerAlgorithm") << "hits collected by the MTF Measurement Collector " << std::endl
00275 << "trajectory number " << i << "has measurements" << collectedmeas.size();
00276
00277 if (collectedmeas.empty()) {
00278 LogDebug("MTFTrackProducerAlgorithm") << "In method collectHits, we had a problem and no measurements were collected "
00279 <<"for trajectory number " << i << std::endl;
00280 return std::make_pair(TransientTrackingRecHit::RecHitContainer(), TrajectoryStateOnSurface());
00281 }
00282
00283 for (std::vector<TrajectoryMeasurement>::const_iterator iter = collectedmeas.begin(); iter!=collectedmeas.end(); iter++){
00284 hits.push_back(iter->recHit());
00285 if(iter->recHit()->isValid())
00286 LogDebug("MTFTrackProducerAlgorithm") << "this MultiRecHit has size: " << iter->recHit()->recHits().size();
00287 }
00288
00289 return std::make_pair(hits,TrajectoryStateWithArbitraryError()(collectedmeas.front().predictedState()));
00290 }
00291
00292
00293 std::pair<TransientTrackingRecHit::RecHitContainer, TrajectoryStateOnSurface>
00294 MTFTrackProducerAlgorithm::updateHits(const std::map<int, std::vector<TrajectoryMeasurement> >& mapvtm,
00295 const MultiTrackFilterHitCollector* measurementCollector,
00296 const SiTrackerMultiRecHitUpdatorMTF* updator,
00297 double annealing,
00298 const TransientTrackingRecHitBuilder* builder,
00299 int i) const {
00300 using namespace std;
00301
00302 std::map< int, std::vector<TrajectoryMeasurement> >::const_iterator itmeas = mapvtm.find(i);
00303 LogDebug("SimpleMTFHitCollector") << "found the element "<< i
00304 << "in the map" << std::endl;
00305
00306 std::vector<TrajectoryMeasurement> meas = itmeas->second;
00307 TransientTrackingRecHit::RecHitContainer multirechits;
00308
00309 if (meas.empty()) {LogDebug("MTFTrackProducerAlgorithm::updateHits") << "Warning!!!The trajectory measurement vector is empty" << "\n";}
00310
00311 for(vector<TrajectoryMeasurement>::reverse_iterator imeas = meas.rbegin(); imeas != meas.rend(); imeas++)
00312
00313 {
00314 if( imeas->recHit()->isValid() )
00315
00316 {
00317
00318
00319 std::vector<const TrackingRecHit*> trechit = imeas->recHit()->recHits();
00320
00321
00322
00323
00324 LogDebug("MTFTrackProducerAlgorithm::updateHits") << "multirechit vector size: " << trechit.size() << "\n";
00325
00326 TransientTrackingRecHit::RecHitContainer hits;
00327
00328 for (vector<const TrackingRecHit*>::iterator itrechit=trechit.begin(); itrechit!=trechit.end(); itrechit++)
00329 {
00330
00331 hits.push_back( builder->build(*itrechit));
00332
00333
00334 }
00335
00336 MultiTrajectoryMeasurement multitm = measurementCollector->TSOSfinder(mapvtm, *imeas, i);
00337
00338 TrajectoryStateCombiner statecombiner;
00339 TrajectoryStateOnSurface combinedtsos = statecombiner.combine(imeas->predictedState(), imeas->backwardPredictedState());
00340 TrajectoryStateOnSurface predictedtsos = imeas->predictedState();
00341
00342
00343
00344 TransientTrackingRecHit::RecHitPointer mrh = updator->buildMultiRecHit(combinedtsos, hits, &multitm, annealing);
00345
00346
00347
00348
00349 multirechits.push_back(mrh);
00350
00351
00352 }
00353 else
00354 {
00355
00356 multirechits.push_back(imeas->recHit());
00357
00358 }
00359
00360
00361 }
00362
00363
00364 return std::make_pair(multirechits,TrajectoryStateWithArbitraryError()(itmeas->second.back().predictedState()));
00365
00366 }
00367
00368
00369
00370 bool MTFTrackProducerAlgorithm::fit(const std::pair<TransientTrackingRecHit::RecHitContainer, TrajectoryStateOnSurface>& hits,
00371 const TrajectoryFitter * theFitter,
00372 std::vector<Trajectory>& vtraj) const {
00373
00374 std::vector<Trajectory> newVec = theFitter->fit(TrajectorySeed(PTrajectoryStateOnDet(),
00375 BasicTrajectorySeed::recHitContainer(),
00376 vtraj.front().seed().direction()),
00377 hits.first,
00378 hits.second);
00379
00380
00381 if(newVec.size())
00382 {
00383 vtraj.reserve(newVec.size());
00384 vtraj.swap(newVec);
00385 LogDebug("MTFTrackProducerAlgorithm") << "swapped!" << std::endl;
00386 }
00387
00388
00389 else
00390 {
00391 LogDebug("MTFTrackProducerAlgorithm") <<" somewhwere, something went terribly wrong...in fitting or smoothing trajectory with measurements:"
00392 << vtraj.front().measurements().size()
00393 <<" was broken\n. We keep the old trajectory"
00394 << std::endl;
00395 return false;
00396 }
00397
00398 return true;
00399
00400 }
00401
00402
00403 bool MTFTrackProducerAlgorithm::buildTrack(const std::vector<Trajectory>& vtraj,
00404 AlgoProductCollection& algoResults,
00405 float ndof,
00406 const reco::BeamSpot& bs) const{
00407
00408 reco::Track * theTrack;
00409 Trajectory * theTraj;
00410
00411
00412
00413
00414
00415 LogDebug("MTFTrackProducerAlgorithm") <<" FITTER FOUND "<< vtraj.size() << " TRAJECTORIES" << std::endl;;
00416 TrajectoryStateOnSurface innertsos;
00417
00418 if (vtraj.size() != 0){
00419
00420 theTraj = new Trajectory( vtraj.front() );
00421
00422 if (theTraj->direction() == alongMomentum) {
00423
00424 innertsos = theTraj->firstMeasurement().updatedState();
00425
00426 } else {
00427 innertsos = theTraj->lastMeasurement().updatedState();
00428 }
00429
00430 TSCBLBuilderNoMaterial tscblBuilder;
00431 TrajectoryStateClosestToBeamLine tscbl = tscblBuilder(*(innertsos.freeState()),bs);
00432
00433 if (tscbl.isValid()==false) return false;
00434
00435 GlobalPoint v = tscbl.trackStateAtPCA().position();
00436 math::XYZPoint pos( v.x(), v.y(), v.z() );
00437 GlobalVector p = tscbl.trackStateAtPCA().momentum();
00438 math::XYZVector mom( p.x(), p.y(), p.z() );
00439
00440 LogDebug("TrackProducer") <<v<<p<<std::endl;
00441
00442 theTrack = new reco::Track(theTraj->chiSquared(),
00443 ndof,
00444 pos, mom, tscbl.trackStateAtPCA().charge(), tscbl.trackStateAtPCA().curvilinearError());
00445
00446
00447 LogDebug("TrackProducer") <<"track done\n";
00448
00449 AlgoProduct aProduct(theTraj,std::make_pair(theTrack,theTraj->direction()));
00450 LogDebug("TrackProducer") <<"track done1\n";
00451 algoResults.push_back(aProduct);
00452 LogDebug("TrackProducer") <<"track done2\n";
00453
00454 return true;
00455 }
00456 else return false;
00457 }
00458
00459 void MTFTrackProducerAlgorithm::filter(const TrajectoryFitter* fitter,
00460 std::vector<Trajectory>& input,
00461 int minhits,
00462 std::vector<Trajectory>& output) const {
00463 if (input.empty()) return;
00464
00465 int ngoodhits = 0;
00466
00467 std::vector<TrajectoryMeasurement> vtm = input[0].measurements();
00468
00469 TransientTrackingRecHit::RecHitContainer hits;
00470
00471
00472 for (std::vector<TrajectoryMeasurement>::reverse_iterator tm=vtm.rbegin(); tm!=vtm.rend();tm++){
00473
00474 if (tm->recHit()->isValid()) {
00475 TransientTrackingRecHit::ConstRecHitContainer components = tm->recHit()->transientHits();
00476 bool isGood = false;
00477 for (TransientTrackingRecHit::ConstRecHitContainer::iterator rechit = components.begin(); rechit != components.end(); rechit++){
00478
00479 if ((*rechit)->weight()>1e-6) {ngoodhits++; isGood = true; break;}
00480 }
00481 if (isGood) hits.push_back(tm->recHit()->clone(tm->updatedState()));
00482 else hits.push_back(InvalidTransientRecHit::build(tm->recHit()->det(), TrackingRecHit::missing));
00483 } else {
00484 hits.push_back(tm->recHit()->clone(tm->updatedState()));
00485 }
00486 }
00487
00488
00489 LogDebug("DAFTrackProducerAlgorithm") << "Original number of valid hits " << input[0].foundHits() << "; after filtering " << ngoodhits;
00490
00491 if (ngoodhits>input[0].foundHits()) edm::LogError("DAFTrackProducerAlgorithm") << "Something wrong: the number of good hits from DAFTrackProducerAlgorithm::filter " << ngoodhits << " is higher than the original one " << input[0].foundHits();
00492
00493 if (ngoodhits < minhits) return;
00494
00495 TrajectoryStateOnSurface curstartingTSOS = input.front().lastMeasurement().updatedState();
00496 LogDebug("DAFTrackProducerAlgorithm") << "starting tsos for final refitting " << curstartingTSOS ;
00497
00498
00499 output = fitter->fit(TrajectorySeed(PTrajectoryStateOnDet(),
00500 BasicTrajectorySeed::recHitContainer(),
00501 input.front().seed().direction()),
00502 hits,
00503 TrajectoryStateWithArbitraryError()(curstartingTSOS));
00504 LogDebug("DAFTrackProducerAlgorithm") << "After filtering " << output.size() << " trajectories";
00505
00506 }
00507
00508 float MTFTrackProducerAlgorithm::calculateNdof(const std::vector<Trajectory>& vtraj) const {
00509 if (vtraj.empty()) return 0;
00510 float ndof = 0;
00511 int nhits=0;
00512 const std::vector<TrajectoryMeasurement>& meas = vtraj.front().measurements();
00513 for (std::vector<TrajectoryMeasurement>::const_iterator iter = meas.begin(); iter != meas.end(); iter++){
00514 if (iter->recHit()->isValid()){
00515 TransientTrackingRecHit::ConstRecHitContainer components = iter->recHit()->transientHits();
00516 TransientTrackingRecHit::ConstRecHitContainer::const_iterator iter2;
00517 for (iter2 = components.begin(); iter2 != components.end(); iter2++){
00518 if ((*iter2)->isValid()){ndof+=((*iter2)->dimension())*(*iter2)->weight();
00519 LogDebug("DAFTrackProducerAlgorithm") << "hit dimension: "<<(*iter2)->dimension()
00520 <<" and weight: "<<(*iter2)->weight();
00521 nhits++;
00522 }
00523 }
00524 }
00525 }
00526
00527 LogDebug("DAFTrackProducerAlgorithm") <<"nhits: "<<nhits<< " ndof: "<<ndof-5;
00528 return ndof-5;
00529 }