CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_2_9/src/RecoTracker/FinalTrackSelectors/src/SimpleTrackListMerger.cc

Go to the documentation of this file.
00001 //
00002 // Package:         RecoTracker/FinalTrackSelectors
00003 // Class:           SimpleTrackListMerger
00004 // 
00005 // Description:     TrackList Cleaner and Merger
00006 //
00007 // Original Author: Steve Wagner, stevew@pizero.colorado.edu
00008 // Created:         Sat Jan 14 22:00:00 UTC 2006
00009 //
00010 // $Author: stenson $
00011 // $Date: 2012/02/19 20:21:34 $
00012 // $Revision: 1.29 $
00013 //
00014 
00015 #include <memory>
00016 #include <string>
00017 #include <iostream>
00018 #include <cmath>
00019 #include <vector>
00020 
00021 #include "RecoTracker/FinalTrackSelectors/interface/SimpleTrackListMerger.h"
00022 
00023 #include "DataFormats/TrackerRecHit2D/interface/SiStripMatchedRecHit2DCollection.h"
00024 #include "DataFormats/TrackerRecHit2D/interface/SiStripRecHit2DCollection.h"
00025 #include "DataFormats/TrackerRecHit2D/interface/SiStripRecHit1DCollection.h"
00026 #include "DataFormats/TrajectorySeed/interface/TrajectorySeedCollection.h"
00027 #include "DataFormats/TrackCandidate/interface/TrackCandidateCollection.h"
00028 
00029 #include "FWCore/Framework/interface/ESHandle.h"
00030 
00031 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00032 
00033 #include "Geometry/CommonDetUnit/interface/GeomDetUnit.h"
00034 #include "Geometry/CommonDetUnit/interface/GeomDet.h"
00035 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
00036 #include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
00037 
00038 //#include "DataFormats/TrackReco/src/classes.h"
00039 
00040 #include "RecoTracker/TrackProducer/interface/ClusterRemovalRefSetter.h"
00041 
00042 
00043 namespace cms
00044 {
00045   // VI January 2012   to be migrated to omnicluster (or firstCluster)
00046   edm::ProductID clusterProduct( const TrackingRecHit *hit){
00047     edm::ProductID pID;
00048     //cast it into the proper class     and find productID
00049     DetId detid = hit->geographicalId(); 
00050     uint32_t subdet = detid.subdetId();
00051     if ((subdet == PixelSubdetector::PixelBarrel) || (subdet == PixelSubdetector::PixelEndcap)) {
00052       pID=reinterpret_cast<const SiPixelRecHit *>(hit)->cluster().id();
00053     } else {
00054       const std::type_info &type = typeid(*hit);
00055       if (type == typeid(SiStripRecHit2D)) {
00056         pID=reinterpret_cast<const SiStripRecHit2D *>(hit)->cluster().id();
00057       } else if (type == typeid(SiStripRecHit1D)) {
00058         pID=reinterpret_cast<const SiStripRecHit1D *>(hit)->cluster().id();     
00059       } else if (type == typeid(SiStripMatchedRecHit2D)) {
00060         const SiStripMatchedRecHit2D *mhit = reinterpret_cast<const SiStripMatchedRecHit2D *>(hit);
00061         pID=mhit->monoClusterRef().id();
00062       } else if (type == typeid(ProjectedSiStripRecHit2D)) {
00063         const ProjectedSiStripRecHit2D *phit = reinterpret_cast<const ProjectedSiStripRecHit2D *>(hit);
00064         pID=(&phit->originalHit())->cluster().id();
00065       } else throw cms::Exception("Unknown RecHit Type") << "RecHit of type " << type.name() << " not supported. (use c++filt to demangle the name)";
00066     }
00067         
00068     return pID;}
00069   
00070 
00071 
00072   SimpleTrackListMerger::SimpleTrackListMerger(edm::ParameterSet const& conf) : 
00073     conf_(conf)
00074   {
00075     copyExtras_ = conf_.getUntrackedParameter<bool>("copyExtras", true);
00076 
00077     produces<reco::TrackCollection>();
00078 
00079     makeReKeyedSeeds_ = conf_.getUntrackedParameter<bool>("makeReKeyedSeeds",false);
00080     if (makeReKeyedSeeds_){
00081       copyExtras_=true;
00082       produces<TrajectorySeedCollection>();
00083     }
00084 
00085     if (copyExtras_) {
00086         produces<reco::TrackExtraCollection>();
00087         produces<TrackingRecHitCollection>();
00088     }
00089     produces< std::vector<Trajectory> >();
00090     produces< TrajTrackAssociationCollection >();
00091     
00092 
00093   }
00094 
00095 
00096   // Virtual destructor needed.
00097   SimpleTrackListMerger::~SimpleTrackListMerger() { }  
00098 
00099   // Functions that gets called by framework every event
00100   void SimpleTrackListMerger::produce(edm::Event& e, const edm::EventSetup& es)
00101   {
00102     // retrieve producer name of input TrackCollection(s)
00103     std::string trackProducer1 = conf_.getParameter<std::string>("TrackProducer1");
00104     std::string trackProducer2 = conf_.getParameter<std::string>("TrackProducer2");
00105 
00106     double maxNormalizedChisq =  conf_.getParameter<double>("MaxNormalizedChisq");
00107     double minPT =  conf_.getParameter<double>("MinPT");
00108     unsigned int minFound = (unsigned int)conf_.getParameter<int>("MinFound");
00109     double epsilon =  conf_.getParameter<double>("Epsilon");
00110     bool use_sharesInput = true;
00111     if ( epsilon > 0.0 )use_sharesInput = false;
00112     double shareFrac =  conf_.getParameter<double>("ShareFrac");
00113     double foundHitBonus = conf_.getParameter<double>("FoundHitBonus");
00114     double lostHitPenalty = conf_.getParameter<double>("LostHitPenalty");
00115   
00116     bool promoteQuality = conf_.getParameter<bool>("promoteTrackQuality");
00117     bool allowFirstHitShare = conf_.getParameter<bool>("allowFirstHitShare");
00118 //
00119 
00120     // New track quality should be read from the file
00121     std::string qualityStr = conf_.getParameter<std::string>("newQuality");
00122     reco::TrackBase::TrackQuality qualityToSet;
00123     if (qualityStr != "") {
00124       qualityToSet = reco::TrackBase::qualityByName(conf_.getParameter<std::string>("newQuality"));
00125     }
00126     else 
00127       qualityToSet = reco::TrackBase::undefQuality;
00128 
00129     // extract tracker geometry
00130     //
00131     edm::ESHandle<TrackerGeometry> theG;
00132     es.get<TrackerDigiGeometryRecord>().get(theG);
00133 
00134 //    using namespace reco;
00135 
00136     // get Inputs 
00137     // if 1 input list doesn't exist, make an empty list, issue a warning, and continue
00138     // this allows SimpleTrackListMerger to be used as a cleaner only if handed just one list
00139     // if both input lists don't exist, will issue 2 warnings and generate an empty output collection
00140     //
00141     const reco::TrackCollection *TC1 = 0;
00142     static const reco::TrackCollection s_empty1, s_empty2;
00143     edm::Handle<reco::TrackCollection> trackCollection1;
00144     e.getByLabel(trackProducer1, trackCollection1);
00145     if (trackCollection1.isValid()) {
00146       TC1 = trackCollection1.product();
00147       //std::cout << "1st collection " << trackProducer1 << " has "<< TC1->size() << " tracks" << std::endl ;
00148     } else {
00149       TC1 = &s_empty1;
00150       edm::LogWarning("SimpleTrackListMerger") << "1st TrackCollection " << trackProducer1 << " not found; will only clean 2nd TrackCollection " << trackProducer2 ;
00151     }
00152     const reco::TrackCollection tC1 = *TC1;
00153 
00154     const reco::TrackCollection *TC2 = 0;
00155     edm::Handle<reco::TrackCollection> trackCollection2;
00156     e.getByLabel(trackProducer2, trackCollection2);
00157     if (trackCollection2.isValid()) {
00158       TC2 = trackCollection2.product();
00159       //std::cout << "2nd collection " << trackProducer2 << " has "<< TC2->size() << " tracks" << std::endl ;
00160     } else {
00161         TC2 = &s_empty2;
00162         edm::LogWarning("SimpleTrackListMerger") << "2nd TrackCollection " << trackProducer2 << " not found; will only clean 1st TrackCollection " << trackProducer1 ;
00163     }
00164     const reco::TrackCollection tC2 = *TC2;
00165 
00166     // Step B: create empty output collection
00167     outputTrks = std::auto_ptr<reco::TrackCollection>(new reco::TrackCollection);
00168     refTrks = e.getRefBeforePut<reco::TrackCollection>();      
00169 
00170     if (copyExtras_) {
00171         outputTrkExtras = std::auto_ptr<reco::TrackExtraCollection>(new reco::TrackExtraCollection);
00172         outputTrkExtras->reserve(TC1->size()+TC2->size());
00173         refTrkExtras    = e.getRefBeforePut<reco::TrackExtraCollection>();
00174         outputTrkHits   = std::auto_ptr<TrackingRecHitCollection>(new TrackingRecHitCollection);
00175         outputTrkHits->reserve((TC1->size()+TC2->size())*25);
00176         refTrkHits      = e.getRefBeforePut<TrackingRecHitCollection>();
00177         if (makeReKeyedSeeds_){
00178           outputSeeds = std::auto_ptr<TrajectorySeedCollection>(new TrajectorySeedCollection);
00179           outputSeeds->reserve(TC1->size()+TC2->size());
00180           refTrajSeeds = e.getRefBeforePut<TrajectorySeedCollection>();
00181         }
00182     }
00183 
00184     outputTrajs = std::auto_ptr< std::vector<Trajectory> >(new std::vector<Trajectory>()); 
00185     outputTrajs->reserve(TC1->size()+TC2->size());
00186     outputTTAss = std::auto_ptr< TrajTrackAssociationCollection >(new TrajTrackAssociationCollection());
00187     //outputTTAss->reserve(TC1->size()+TC2->size());//how do I reserve space for an association map?
00188 
00189   //
00190   //  no input tracks
00191   //
00192 
00193 //    if ( tC1.empty() ){
00194 //      LogDebug("RoadSearch") << "Found " << output.size() << " clouds.";
00195 //      e.put(output);
00196 //      return;  
00197 //    }
00198 
00199   //
00200   //  quality cuts first
00201   // 
00202     int i;
00203 
00204     std::vector<int> selected1; for (unsigned int i=0; i<tC1.size(); ++i){selected1.push_back(1);}
00205 
00206    if ( 0<tC1.size() ){
00207       i=-1;
00208       for (reco::TrackCollection::const_iterator track=tC1.begin(); track!=tC1.end(); track++){
00209         i++;
00210         if ((short unsigned)track->ndof() < 1){
00211           selected1[i]=0; 
00212           //std::cout << "L1Track "<< i << " rejected in SimpleTrackListMerger; ndof() < 1" << std::endl ;
00213           continue;
00214         }
00215         if (track->normalizedChi2() > maxNormalizedChisq){
00216           selected1[i]=0; 
00217           //std::cout << "L1Track "<< i << " rejected in SimpleTrackListMerger; normalizedChi2() > maxNormalizedChisq " << track->normalizedChi2() << " " << maxNormalizedChisq << std::endl ;
00218           continue;
00219         }
00220         if (track->found() < minFound){
00221           selected1[i]=0; 
00222           //std::cout << "L1Track "<< i << " rejected in SimpleTrackListMerger; found() < minFound " << track->found() << " " << minFound << std::endl ;
00223           continue;
00224         }
00225         if (track->pt() < minPT){
00226           selected1[i]=0; 
00227           //std::cout << "L1Track "<< i << " rejected in SimpleTrackListMerger; pt() < minPT " << track->pt() << " " << minPT << std::endl ;
00228           continue;
00229         }
00230       }//end loop over tracks
00231    }//end more than 0 track
00232 
00233 
00234     std::vector<int> selected2; for (unsigned int i=0; i<tC2.size(); ++i){selected2.push_back(1);}
00235 
00236    if ( 0<tC2.size() ){
00237       i=-1;
00238       for (reco::TrackCollection::const_iterator track=tC2.begin(); track!=tC2.end(); track++){
00239         i++;
00240         if ((short unsigned)track->ndof() < 1){
00241           selected2[i]=0; 
00242           //std::cout << "L2Track "<< i << " rejected in SimpleTrackListMerger; ndof() < 1" << std::endl ;
00243           continue;
00244         }
00245         if (track->normalizedChi2() > maxNormalizedChisq){
00246           selected2[i]=0; 
00247           //std::cout << "L2Track "<< i << " rejected in SimpleTrackListMerger; normalizedChi2() > maxNormalizedChisq " << track->normalizedChi2() << " " << maxNormalizedChisq << std::endl ;
00248           continue;
00249         }
00250         if (track->found() < minFound){
00251           selected2[i]=0; 
00252           //std::cout << "L2Track "<< i << " rejected in SimpleTrackListMerger; found() < minFound " << track->found() << " " << minFound << std::endl ;
00253           continue;
00254         }
00255         if (track->pt() < minPT){
00256           selected2[i]=0; 
00257           //std::cout << "L2Track "<< i << " rejected in SimpleTrackListMerger; pt() < minPT " << track->pt() << " " << minPT << std::endl ;
00258           continue;
00259         }
00260       }//end loop over tracks
00261    }//end more than 0 track
00262 
00263    std::map<reco::TrackCollection::const_iterator, std::vector<const TrackingRecHit*> > rh1;
00264    std::map<reco::TrackCollection::const_iterator, std::vector<const TrackingRecHit*> > rh2;
00265    for (reco::TrackCollection::const_iterator track=tC1.begin(); track!=tC1.end(); ++track){
00266      trackingRecHit_iterator itB = track->recHitsBegin();
00267      trackingRecHit_iterator itE = track->recHitsEnd();
00268      for (trackingRecHit_iterator it = itB;  it != itE; ++it) { 
00269        const TrackingRecHit* hit = &(**it);
00270        rh1[track].push_back(hit);
00271      }
00272    }
00273    for (reco::TrackCollection::const_iterator track=tC2.begin(); track!=tC2.end(); ++track){
00274      trackingRecHit_iterator jtB = track->recHitsBegin();
00275      trackingRecHit_iterator jtE = track->recHitsEnd();
00276      for (trackingRecHit_iterator jt = jtB;  jt != jtE; ++jt) { 
00277        const TrackingRecHit* hit = &(**jt);
00278        rh2[track].push_back(hit);
00279      }
00280    }
00281 
00282    if ( (0<tC1.size())&&(0<tC2.size()) ){
00283     i=-1;
00284     for (reco::TrackCollection::const_iterator track=tC1.begin(); track!=tC1.end(); ++track){
00285       i++; 
00286       if (!selected1[i])continue;
00287       std::vector<const TrackingRecHit*>& iHits = rh1[track]; 
00288       unsigned nh1 = iHits.size();
00289       int qualityMaskT1 = track->qualityMask();
00290       int j=-1;
00291       for (reco::TrackCollection::const_iterator track2=tC2.begin(); track2!=tC2.end(); ++track2){
00292         j++;
00293         if ((!selected2[j])||(!selected1[i]))continue;
00294         std::vector<const TrackingRecHit*>& jHits = rh2[track2]; 
00295         unsigned nh2 = jHits.size();
00296         int noverlap=0;
00297         int firstoverlap=0;
00298         for ( unsigned ih=0; ih<nh1; ++ih ) { 
00299           const TrackingRecHit* it = iHits[ih];
00300           if (it->isValid()){
00301             int jj=-1;
00302             for ( unsigned jh=0; jh<nh2; ++jh ) { 
00303               const TrackingRecHit* jt = jHits[jh];
00304               jj++;
00305               if (jt->isValid()){
00306                if (!use_sharesInput){
00307                 float delta = fabs ( it->localPosition().x()-jt->localPosition().x() ); 
00308                 if ((it->geographicalId()==jt->geographicalId())&&(delta<epsilon)) {
00309                   noverlap++;
00310                   if ( allowFirstHitShare && ( ih == 0 ) && ( jh == 0 ) ) firstoverlap=1;
00311                 }
00312                }else{
00313                 if ( it->sharesInput(jt,TrackingRecHit::some) ) {
00314                   noverlap++;
00315                   if ( allowFirstHitShare && ( ih == 0 ) && ( jh == 0 ) ) firstoverlap=1;
00316                 }
00317                }
00318               }
00319             }
00320           }
00321         }
00322         int newQualityMask = (qualityMaskT1 | track2->qualityMask()); // take OR of trackQuality 
00323         int nhit1 = track->numberOfValidHits();
00324         int nhit2 = track2->numberOfValidHits();
00325         if ( (noverlap-firstoverlap) > (std::min(nhit1,nhit2)-firstoverlap)*shareFrac ) {
00326           double score1 = foundHitBonus*nhit1 - lostHitPenalty*track->numberOfLostHits() - track->chi2();
00327           double score2 = foundHitBonus*nhit2 - lostHitPenalty*track2->numberOfLostHits() - track2->chi2();
00328           const double almostSame = 1.001;
00329           if ( score1 > almostSame * score2 ){
00330             selected2[j]=0; 
00331             selected1[i]=10+newQualityMask; // add 10 to avoid the case where mask = 1
00332           }else if ( score2 > almostSame * score1 ){
00333               selected1[i]=0; 
00334               selected2[j]=10+newQualityMask;  // add 10 to avoid the case where mask = 1
00335           }else{
00336             if (track->algo() <= track2->algo()) {
00337               selected2[j]=0;
00338               selected1[i]=10+newQualityMask; // add 10 to avoid the case where mask = 1
00339             }else{
00340               selected1[i]=0;
00341               selected2[j]=10+newQualityMask; // add 10 to avoid the case where mask = 1
00342             }
00343           }
00344         }//end got a duplicate
00345       }//end track2 loop
00346     }//end track loop
00347    }//end more than 1 track
00348 
00349   //
00350   //  output selected tracks - if any
00351   //
00352    trackRefs.resize(tC1.size()+tC2.size());
00353    std::vector<edm::RefToBase<TrajectorySeed> > seedsRefs(tC1.size()+tC2.size());
00354    size_t current = 0;
00355 
00356    if ( 0<tC1.size() ){
00357      i=0;
00358      for (reco::TrackCollection::const_iterator track=tC1.begin(); track!=tC1.end(); 
00359           ++track, ++current, ++i){
00360       if (!selected1[i]){
00361         trackRefs[current] = reco::TrackRef();
00362         continue;
00363       }
00364       const reco::Track & theTrack = * track;
00365       //fill the TrackCollection
00366       outputTrks->push_back( reco::Track( theTrack ) );
00367       if (selected1[i]>1 && promoteQuality){
00368         outputTrks->back().setQualityMask(selected1[i]-10);
00369         outputTrks->back().setQuality(qualityToSet);
00370       }
00371       if (copyExtras_) {
00372           //--------NEW----------
00373           edm::RefToBase<TrajectorySeed> origSeedRef = theTrack.seedRef();
00374           //creating a seed with rekeyed clusters if required
00375           if (makeReKeyedSeeds_){
00376             bool doRekeyOnThisSeed=false;
00377             
00378             edm::InputTag clusterRemovalInfos("");
00379             //grab on of the hits of the seed
00380             if (origSeedRef->nHits()!=0){
00381               TrajectorySeed::const_iterator firstHit=origSeedRef->recHits().first;
00382               const TrackingRecHit *hit = &*firstHit;
00383               if (firstHit->isValid()){
00384                 edm::ProductID  pID=clusterProduct(hit);
00385                 // the cluster collection either produced a removalInfo or mot
00386                 //get the clusterremoval info from the provenance: will rekey if this is found
00387                 edm::Handle<reco::ClusterRemovalInfo> CRIh;
00388                 edm::Provenance prov=e.getProvenance(pID);
00389                 clusterRemovalInfos=edm::InputTag(prov.moduleLabel(),
00390                                                   prov.productInstanceName(),
00391                                                   prov.processName());
00392                 doRekeyOnThisSeed=e.getByLabel(clusterRemovalInfos,CRIh);
00393               }//valid hit
00394             }//nhit!=0
00395             
00396             if (doRekeyOnThisSeed && !(clusterRemovalInfos==edm::InputTag("")))
00397               {
00398                 ClusterRemovalRefSetter refSetter(e,clusterRemovalInfos);
00399                 TrajectorySeed::recHitContainer  newRecHitContainer;
00400                 newRecHitContainer.reserve(origSeedRef->nHits());
00401                 TrajectorySeed::const_iterator iH=origSeedRef->recHits().first;
00402                 TrajectorySeed::const_iterator iH_end=origSeedRef->recHits().second;
00403                 for (;iH!=iH_end;++iH){
00404                   newRecHitContainer.push_back(*iH);
00405                   refSetter.reKey(&newRecHitContainer.back());
00406                 }
00407                 outputSeeds->push_back( TrajectorySeed( origSeedRef->startingState(),
00408                                                         newRecHitContainer,
00409                                                         origSeedRef->direction()));
00410               }
00411             //doRekeyOnThisSeed=true
00412             else{
00413               //just copy the one we had before
00414               outputSeeds->push_back( TrajectorySeed(*origSeedRef));
00415             }
00416             edm::Ref<TrajectorySeedCollection> pureRef(refTrajSeeds, outputSeeds->size()-1);
00417             origSeedRef=edm::RefToBase<TrajectorySeed>( pureRef);
00418           }//creating a new seed and rekeying it rechit clusters.
00419           //--------NEW----------
00420           // Fill TrackExtra collection
00421           outputTrkExtras->push_back( reco::TrackExtra( 
00422                         theTrack.outerPosition(), theTrack.outerMomentum(), theTrack.outerOk(),
00423                         theTrack.innerPosition(), theTrack.innerMomentum(), theTrack.innerOk(),
00424                         theTrack.outerStateCovariance(), theTrack.outerDetId(),
00425                         theTrack.innerStateCovariance(), theTrack.innerDetId(),
00426                         theTrack.seedDirection(), origSeedRef ) );
00427           seedsRefs[current]=origSeedRef;
00428           outputTrks->back().setExtra( reco::TrackExtraRef( refTrkExtras, outputTrkExtras->size() - 1) );
00429           reco::TrackExtra & tx = outputTrkExtras->back();
00430           tx.setResiduals(theTrack.residuals());
00431           // fill TrackingRecHits
00432           std::vector<const TrackingRecHit*>& iHits = rh1[track]; 
00433           unsigned nh1 = iHits.size();
00434           for ( unsigned ih=0; ih<nh1; ++ih ) { 
00435             const TrackingRecHit* hit = iHits[ih];
00436             //for( trackingRecHit_iterator hit = itB; hit != itE; ++hit ) {
00437             outputTrkHits->push_back( hit->clone() );
00438             tx.add( TrackingRecHitRef( refTrkHits, outputTrkHits->size() - 1) );
00439           }
00440       }
00441       trackRefs[current] = reco::TrackRef(refTrks, outputTrks->size() - 1);
00442     
00443 
00444     }//end faux loop over tracks
00445    }//end more than 0 track
00446 
00447    //Fill the trajectories, etc. for 1st collection
00448    edm::Handle< std::vector<Trajectory> > hTraj1;
00449    e.getByLabel(trackProducer1, hTraj1);
00450    edm::Handle< TrajTrackAssociationCollection > hTTAss1;
00451    e.getByLabel(trackProducer1, hTTAss1);
00452    refTrajs    = e.getRefBeforePut< std::vector<Trajectory> >();
00453 
00454    if (!hTraj1.failedToGet() && !hTTAss1.failedToGet()){
00455    for (size_t i = 0, n = hTraj1->size(); i < n; ++i) {
00456      edm::Ref< std::vector<Trajectory> > trajRef(hTraj1, i);
00457      TrajTrackAssociationCollection::const_iterator match = hTTAss1->find(trajRef);
00458      if (match != hTTAss1->end()) {
00459        const edm::Ref<reco::TrackCollection> &trkRef = match->val; 
00460        short oldKey = static_cast<short>(trkRef.key());
00461        if (trackRefs[oldKey].isNonnull()) {
00462          outputTrajs->push_back( *trajRef );
00463          //if making extras and the seeds at the same time, change the seed ref on the trajectory
00464          if (copyExtras_ && makeReKeyedSeeds_)
00465              outputTrajs->back().setSeedRef( seedsRefs[oldKey] );
00466          outputTTAss->insert ( edm::Ref< std::vector<Trajectory> >(refTrajs, outputTrajs->size() - 1), 
00467                                trackRefs[oldKey] );
00468        }
00469      }
00470    }
00471    }
00472 
00473    short offset = current; //save offset into trackRefs
00474 
00475    if ( 0<tC2.size() ){
00476     i=0;
00477     for (reco::TrackCollection::const_iterator track=tC2.begin(); track!=tC2.end();
00478          ++track, ++current, ++i){
00479       if (!selected2[i]){
00480         trackRefs[current] = reco::TrackRef();
00481         continue;
00482       }
00483       const reco::Track & theTrack = * track;
00484       //fill the TrackCollection
00485       outputTrks->push_back( reco::Track( theTrack ) );
00486       if (selected2[i]>1 && promoteQuality){
00487         outputTrks->back().setQualityMask(selected2[i]-10);
00488         outputTrks->back().setQuality(qualityToSet);
00489       }
00490       if (copyExtras_) {
00491           //--------NEW----------
00492           edm::RefToBase<TrajectorySeed> origSeedRef = theTrack.seedRef();
00493           //creating a seed with rekeyed clusters if required
00494           if (makeReKeyedSeeds_){
00495             bool doRekeyOnThisSeed=false;
00496             
00497             edm::InputTag clusterRemovalInfos("");
00498             //grab on of the hits of the seed
00499             if (origSeedRef->nHits()!=0){
00500               TrajectorySeed::const_iterator firstHit=origSeedRef->recHits().first;
00501               const TrackingRecHit *hit = &*firstHit;
00502               if (firstHit->isValid()){
00503                 edm::ProductID  pID=clusterProduct(hit);
00504                 // the cluster collection either produced a removalInfo or mot
00505                 //get the clusterremoval info from the provenance: will rekey if this is found
00506                 edm::Handle<reco::ClusterRemovalInfo> CRIh;
00507                 edm::Provenance prov=e.getProvenance(pID);
00508                 clusterRemovalInfos=edm::InputTag(prov.moduleLabel(),
00509                                                   prov.productInstanceName(),
00510                                                   prov.processName());
00511                 doRekeyOnThisSeed=e.getByLabel(clusterRemovalInfos,CRIh);
00512               }//valid hit
00513             }//nhit!=0
00514             
00515             if (doRekeyOnThisSeed && !(clusterRemovalInfos==edm::InputTag("")))
00516               {
00517                 ClusterRemovalRefSetter refSetter(e,clusterRemovalInfos);
00518                 TrajectorySeed::recHitContainer  newRecHitContainer;
00519                 newRecHitContainer.reserve(origSeedRef->nHits());
00520                 TrajectorySeed::const_iterator iH=origSeedRef->recHits().first;
00521                 TrajectorySeed::const_iterator iH_end=origSeedRef->recHits().second;
00522                 for (;iH!=iH_end;++iH){
00523                   newRecHitContainer.push_back(*iH);
00524                   refSetter.reKey(&newRecHitContainer.back());
00525                 }
00526                 outputSeeds->push_back( TrajectorySeed( origSeedRef->startingState(),
00527                                                         newRecHitContainer,
00528                                                         origSeedRef->direction()));
00529               }//doRekeyOnThisSeed=true
00530               else{
00531                 //just copy the one we had before
00532                 outputSeeds->push_back( TrajectorySeed(*origSeedRef));
00533               }
00534             edm::Ref<TrajectorySeedCollection> pureRef(refTrajSeeds, outputSeeds->size()-1);
00535             origSeedRef=edm::RefToBase<TrajectorySeed>( pureRef);
00536           }//creating a new seed and rekeying it rechit clusters.
00537           //--------NEW----------
00538           // Fill TrackExtra collection
00539           outputTrkExtras->push_back( reco::TrackExtra( 
00540                         theTrack.outerPosition(), theTrack.outerMomentum(), theTrack.outerOk(),
00541                         theTrack.innerPosition(), theTrack.innerMomentum(), theTrack.innerOk(),
00542                         theTrack.outerStateCovariance(), theTrack.outerDetId(),
00543                         theTrack.innerStateCovariance(), theTrack.innerDetId(),
00544                         theTrack.seedDirection(), origSeedRef ) );
00545           seedsRefs[current]=origSeedRef;
00546           outputTrks->back().setExtra( reco::TrackExtraRef( refTrkExtras, outputTrkExtras->size() - 1) );
00547           reco::TrackExtra & tx = outputTrkExtras->back();
00548           tx.setResiduals(theTrack.residuals());
00549           // fill TrackingRecHits
00550           std::vector<const TrackingRecHit*>& jHits = rh2[track]; 
00551           unsigned nh2 = jHits.size();
00552           for ( unsigned jh=0; jh<nh2; ++jh ) { 
00553             const TrackingRecHit* hit = jHits[jh];
00554             outputTrkHits->push_back( hit->clone() );
00555             tx.add( TrackingRecHitRef( refTrkHits, outputTrkHits->size() - 1) );
00556           }
00557       }
00558       trackRefs[current] = reco::TrackRef(refTrks, outputTrks->size() - 1);
00559 
00560     }//end faux loop over tracks
00561    }//end more than 0 track
00562 
00563    //Fill the trajectories, etc. for 2nd collection
00564    edm::Handle< std::vector<Trajectory> > hTraj2;
00565    e.getByLabel(trackProducer2, hTraj2);
00566    edm::Handle< TrajTrackAssociationCollection > hTTAss2;
00567    e.getByLabel(trackProducer2, hTTAss2);
00568 
00569    if (!hTraj2.failedToGet() && !hTTAss2.failedToGet()){
00570    for (size_t i = 0, n = hTraj2->size(); i < n; ++i) {
00571      edm::Ref< std::vector<Trajectory> > trajRef(hTraj2, i);
00572      TrajTrackAssociationCollection::const_iterator match = hTTAss2->find(trajRef);
00573      if (match != hTTAss2->end()) {
00574        const edm::Ref<reco::TrackCollection> &trkRef = match->val; 
00575        short oldKey = static_cast<short>(trkRef.key()) + offset;
00576        if (trackRefs[oldKey].isNonnull()) {
00577            outputTrajs->push_back( Trajectory(*trajRef) );
00578            //if making extras and the seeds at the same time, change the seed ref on the trajectory
00579            if (copyExtras_ && makeReKeyedSeeds_)
00580              outputTrajs->back().setSeedRef( seedsRefs[oldKey] );
00581            outputTTAss->insert ( edm::Ref< std::vector<Trajectory> >(refTrajs, outputTrajs->size() - 1), 
00582                    trackRefs[oldKey] ); 
00583        }
00584      }
00585    }}
00586 
00587     e.put(outputTrks);
00588     if (copyExtras_) {
00589         e.put(outputTrkExtras);
00590         e.put(outputTrkHits);
00591         if (makeReKeyedSeeds_)
00592           e.put(outputSeeds);
00593     }
00594     e.put(outputTrajs);
00595     e.put(outputTTAss);
00596     return;
00597 
00598   }//end produce
00599 }