Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:29:00

0001 #include "RecoTracker/TrackProducer/interface/KfTrackProducerBase.h"
0002 // system include files
0003 #include <memory>
0004 // user include files
0005 #include "DataFormats/TrackReco/interface/TrackResiduals.h"
0006 #include "FWCore/Framework/interface/Frameworkfwd.h"
0007 #include "FWCore/Framework/interface/MakerMacros.h"
0008 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0009 #include "TrackingTools/GeomPropagators/interface/Propagator.h"
0010 #include "TrackingTools/PatternTools/interface/Trajectory.h"
0011 
0012 #include "TrackingTools/PatternTools/interface/TrajTrackAssociation.h"
0013 
0014 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
0015 
0016 #include "TrackingTools/PatternTools/interface/ClusterRemovalRefSetter.h"
0017 #include "TrajectoryToResiduals.h"
0018 
0019 #include "RecoTracker/TransientTrackingRecHit/interface/Traj2TrackHits.h"
0020 
0021 void KfTrackProducerBase::putInEvt(edm::Event& evt,
0022                                    const Propagator* prop,
0023                                    const MeasurementTracker* measTk,
0024                                    std::unique_ptr<TrackingRecHitCollection>& selHits,
0025                                    std::unique_ptr<reco::TrackCollection>& selTracks,
0026                                    std::unique_ptr<reco::TrackExtraCollection>& selTrackExtras,
0027                                    std::unique_ptr<std::vector<Trajectory> >& selTrajectories,
0028                                    std::unique_ptr<std::vector<int> >& indecesInput,
0029                                    AlgoProductCollection& algoResults,
0030                                    TransientTrackingRecHitBuilder const* hitBuilder,
0031                                    const TrackerTopology* ttopo,
0032                                    int BeforeOrAfter) {
0033   TrackingRecHitRefProd rHits = evt.getRefBeforePut<TrackingRecHitCollection>();
0034   reco::TrackExtraRefProd rTrackExtras = evt.getRefBeforePut<reco::TrackExtraCollection>();
0035 
0036   edm::Ref<reco::TrackExtraCollection>::key_type idx = 0;
0037   edm::Ref<reco::TrackCollection>::key_type iTkRef = 0;
0038   edm::Ref<std::vector<Trajectory> >::key_type iTjRef = 0;
0039   std::map<unsigned int, unsigned int> tjTkMap;
0040 
0041   selTracks->reserve(algoResults.size());
0042   selTrackExtras->reserve(algoResults.size());
0043   if (trajectoryInEvent_)
0044     selTrajectories->reserve(algoResults.size());
0045 
0046   for (AlgoProductCollection::iterator i = algoResults.begin(); i != algoResults.end(); i++) {
0047     auto theTraj = (*i).trajectory;
0048     (*indecesInput).push_back((*i).indexInput);
0049     if (trajectoryInEvent_) {
0050       selTrajectories->push_back(*theTraj);
0051       iTjRef++;
0052     }
0053 
0054     auto theTrack = (*i).track;
0055 
0056     // Hits are going to be re-sorted along momentum few lines later.
0057     // Therefore the direction stored in the TrackExtra
0058     // has to be "alongMomentum" as well. Anyway, this direction can be differnt from the one of the orignal
0059     // seed! The name seedDirection() for the Track's method (and the corresponding data member) is
0060     // misleading and should be changed into something like "hitsDirection()". TO BE FIXED!
0061 
0062     PropagationDirection seedDir = alongMomentum;
0063 
0064     LogDebug("TrackProducer") << "In KfTrackProducerBase::putInEvt - seedDir=" << seedDir;
0065 
0066     selTracks->push_back(std::move(*theTrack));
0067     delete theTrack;
0068     iTkRef++;
0069 
0070     // Store indices in local map (starts at 0)
0071     if (trajectoryInEvent_)
0072       tjTkMap[iTjRef - 1] = iTkRef - 1;
0073 
0074     //sets the outermost and innermost TSOSs
0075     TrajectoryStateOnSurface outertsos;
0076     TrajectoryStateOnSurface innertsos;
0077     unsigned int innerId, outerId;
0078 
0079     // ---  NOTA BENE: the convention is to sort hits and measurements "along the momentum".
0080     // This is consistent with innermost and outermost labels only for tracks from LHC collision
0081     if (theTraj->direction() == alongMomentum) {
0082       outertsos = theTraj->lastMeasurement().updatedState();
0083       innertsos = theTraj->firstMeasurement().updatedState();
0084       outerId = theTraj->lastMeasurement().recHit()->geographicalId().rawId();
0085       innerId = theTraj->firstMeasurement().recHit()->geographicalId().rawId();
0086     } else {
0087       outertsos = theTraj->firstMeasurement().updatedState();
0088       innertsos = theTraj->lastMeasurement().updatedState();
0089       outerId = theTraj->firstMeasurement().recHit()->geographicalId().rawId();
0090       innerId = theTraj->lastMeasurement().recHit()->geographicalId().rawId();
0091     }
0092     // ---
0093     //build the TrackExtra
0094     GlobalPoint v = outertsos.globalParameters().position();
0095     GlobalVector p = outertsos.globalParameters().momentum();
0096     math::XYZVector outmom(p.x(), p.y(), p.z());
0097     math::XYZPoint outpos(v.x(), v.y(), v.z());
0098     v = innertsos.globalParameters().position();
0099     p = innertsos.globalParameters().momentum();
0100     math::XYZVector inmom(p.x(), p.y(), p.z());
0101     math::XYZPoint inpos(v.x(), v.y(), v.z());
0102 
0103     reco::TrackExtraRef teref = reco::TrackExtraRef(rTrackExtras, idx++);
0104     reco::Track& track = selTracks->back();
0105     track.setExtra(teref);
0106     //======= I want to set the second hitPattern here =============
0107     if (theSchool.isValid()) {
0108       edm::Handle<MeasurementTrackerEvent> mte;
0109       evt.getByToken(mteSrc_, mte);
0110       // NavigationSetter setter( *theSchool );
0111       setSecondHitPattern(theTraj, track, prop, &*mte, ttopo);
0112     }
0113     //==============================================================
0114 
0115     selTrackExtras->push_back(reco::TrackExtra(outpos,
0116                                                outmom,
0117                                                true,
0118                                                inpos,
0119                                                inmom,
0120                                                true,
0121                                                outertsos.curvilinearError(),
0122                                                outerId,
0123                                                innertsos.curvilinearError(),
0124                                                innerId,
0125                                                seedDir,
0126                                                theTraj->seedRef()));
0127 
0128     // FIXME will remove this obsolete config-param in a future PR
0129     assert(!useSplitting);
0130 
0131     reco::TrackExtra& tx = selTrackExtras->back();
0132     // ---  NOTA BENE: the convention is to sort hits and measurements "along the momentum".
0133     // This is consistent with innermost and outermost labels only for tracks from LHC collisions
0134     reco::TrackExtra::TrajParams trajParams;
0135     reco::TrackExtra::Chi2sFive chi2s;
0136     Traj2TrackHits t2t;
0137     auto ih = selHits->size();
0138     t2t(*theTraj, *selHits, trajParams, chi2s);
0139     auto ie = selHits->size();
0140     tx.setHits(rHits, ih, ie - ih);
0141     tx.setTrajParams(std::move(trajParams), std::move(chi2s));
0142     assert(tx.trajParams().size() == tx.recHitsSize());
0143     for (; ih < ie; ++ih) {
0144       auto const& hit = (*selHits)[ih];
0145       track.appendHitPattern(hit, *ttopo);
0146     }
0147 
0148     // ----
0149     tx.setResiduals(trajectoryToResiduals(*theTraj));
0150 
0151     delete theTraj;
0152   }
0153 
0154   // Now we can re-set refs to hits, as they have already been cloned
0155   if (rekeyClusterRefs_) {
0156     ClusterRemovalRefSetter refSetter(evt, clusterRemovalInfo_);
0157     for (TrackingRecHitCollection::iterator it = selHits->begin(), ed = selHits->end(); it != ed; ++it) {
0158       refSetter.reKey(&*it);
0159     }
0160   }
0161 
0162   LogTrace("TrackingRegressionTest") << "========== TrackProducer Info ===================";
0163   LogTrace("TrackingRegressionTest") << "number of finalTracks: " << selTracks->size();
0164   for (reco::TrackCollection::const_iterator it = selTracks->begin(); it != selTracks->end(); it++) {
0165     LogTrace("TrackingRegressionTest") << "track's n valid and invalid hit, chi2, pt, eta : " << it->found() << " , "
0166                                        << it->lost() << " , " << it->normalizedChi2() << " , " << it->pt() << " , "
0167                                        << it->eta();
0168   }
0169   LogTrace("TrackingRegressionTest") << "=================================================";
0170 
0171   selTracks->shrink_to_fit();
0172   selTrackExtras->shrink_to_fit();
0173   selHits->shrink_to_fit();
0174   indecesInput->shrink_to_fit();
0175   if (BeforeOrAfter == 1) {
0176     rTracks_ = evt.put(std::move(selTracks), "beforeDAF");
0177     evt.put(std::move(selTrackExtras), "beforeDAF");
0178     evt.put(std::move(indecesInput), "beforeDAF");
0179   } else if (BeforeOrAfter == 2) {
0180     rTracks_ = evt.put(std::move(selTracks), "afterDAF");
0181     evt.put(std::move(selTrackExtras), "afterDAF");
0182     evt.put(std::move(indecesInput), "afterDAF");
0183   } else {
0184     rTracks_ = evt.put(std::move(selTracks));
0185     evt.put(std::move(selTrackExtras));
0186     evt.put(std::move(selHits));
0187     evt.put(std::move(indecesInput));
0188   }
0189 
0190   if (trajectoryInEvent_ && BeforeOrAfter == 0) {
0191     selTrajectories->shrink_to_fit();
0192     edm::OrphanHandle<std::vector<Trajectory> > rTrajs = evt.put(std::move(selTrajectories));
0193 
0194     // Now Create traj<->tracks association map
0195     std::unique_ptr<TrajTrackAssociationCollection> trajTrackMap(new TrajTrackAssociationCollection(rTrajs, rTracks_));
0196     for (std::map<unsigned int, unsigned int>::iterator i = tjTkMap.begin(); i != tjTkMap.end(); i++) {
0197       edm::Ref<std::vector<Trajectory> > trajRef(rTrajs, (*i).first);
0198       edm::Ref<reco::TrackCollection> tkRef(rTracks_, (*i).second);
0199       trajTrackMap->insert(edm::Ref<std::vector<Trajectory> >(rTrajs, (*i).first),
0200                            edm::Ref<reco::TrackCollection>(rTracks_, (*i).second));
0201     }
0202     evt.put(std::move(trajTrackMap));
0203   }
0204 }