Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 #include "DataFormats/TrackCandidate/interface/TrackCandidate.h"
0002 #include "DataFormats/TrackingRecHit/interface/TrackingRecHitFwd.h"
0003 #include "DataFormats/TrackingRecHit/interface/InvalidTrackingRecHit.h"
0004 #include "RecoTracker/TransientTrackingRecHit/interface/TkTransientTrackingRecHitBuilder.h"
0005 #include "DataFormats/TrackReco/interface/Track.h"
0006 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0007 #include "Geometry/CommonDetUnit/interface/TrackingGeometry.h"
0008 #include "MagneticField/Engine/interface/MagneticField.h"
0009 #include "RecoTracker/TrackProducer/interface/DAFTrackProducerAlgorithm.h"
0010 #include "RecoTracker/SiTrackerMRHTools/interface/SiTrackerMultiRecHitUpdator.h"
0011 #include "RecoTracker/SiTrackerMRHTools/interface/MultiRecHitCollector.h"
0012 #include "TrackingTools/TrackFitters/interface/TrajectoryFitter.h"
0013 #include "TrackingTools/TrajectoryState/interface/TrajectoryStateOnSurface.h"
0014 #include "TrackingTools/TrajectoryState/interface/TrajectoryStateTransform.h"
0015 #include "TrackingTools/TransientTrackingRecHit/interface/TransientTrackingRecHitBuilder.h"
0016 #include "RecoTracker/TransientTrackingRecHit/interface/TkClonerImpl.h"
0017 #include "TrackingTools/PatternTools/interface/TSCBLBuilderNoMaterial.h"
0018 #include "TrackingTools/PatternTools/interface/TransverseImpactPointExtrapolator.h"
0019 #include "TrackingTools/TrackFitters/interface/TrajectoryStateWithArbitraryError.h"
0020 #include "TrackingTools/TransientTrackingRecHit/interface/TransientTrackingRecHit.h"
0021 #include "DataFormats/TrackerRecHit2D/interface/SiTrackerMultiRecHit.h"
0022 #include "DataFormats/TrackerRecHit2D/interface/OmniClusterRef.h"
0023 #include "DataFormats/TrackerRecHit2D/interface/TkCloner.h"
0024 #include "TrackingTools/PatternTools/interface/TrajAnnealing.h"
0025 #include "TrackingTools/TrackFitters/interface/TrajectoryStateCombiner.h"
0026 
0027 DAFTrackProducerAlgorithm::DAFTrackProducerAlgorithm(const edm::ParameterSet& conf)
0028     : conf_(conf), minHits_(conf.getParameter<int>("MinHits")) {}
0029 
0030 void DAFTrackProducerAlgorithm::runWithCandidate(const TrackingGeometry* theG,
0031                                                  const MagneticField* theMF,
0032                                                  const TrajTrackAssociationCollection& TTmap,
0033                                                  const MeasurementTrackerEvent* measTk,
0034                                                  const TrajectoryFitter* theFitter,
0035                                                  const TransientTrackingRecHitBuilder* builder,
0036                                                  const MultiRecHitCollector* measurementCollector,
0037                                                  const SiTrackerMultiRecHitUpdator* updator,
0038                                                  const reco::BeamSpot& bs,
0039                                                  AlgoProductCollection& algoResults,
0040                                                  TrajAnnealingCollection& trajann,
0041                                                  bool TrajAnnSaving_,
0042                                                  AlgoProductCollection& algoResultsBeforeDAF,
0043                                                  AlgoProductCollection& algoResultsAfterDAF) const {
0044   LogDebug("DAFTrackProducerAlgorithm") << "Size of map: " << TTmap.size() << "\n";
0045 
0046   int cont = 0;
0047   int nTracksChanged = 0;
0048 
0049   for (TrajTrackAssociationCollection::const_iterator itTTmap = TTmap.begin(); itTTmap != TTmap.end(); itTTmap++) {
0050     const edm::Ref<std::vector<Trajectory> > BeforeDAFTraj = itTTmap->key;
0051     std::vector<TrajectoryMeasurement> BeforeDAFTrajMeas = BeforeDAFTraj->measurements();
0052     const reco::TrackRef BeforeDAFTrack = itTTmap->val;
0053 
0054     float ndof = 0;
0055     Trajectory CurrentTraj;
0056 
0057     if (BeforeDAFTraj->isValid()) {
0058       LogDebug("DAFTrackProducerAlgorithm") << "The trajectory #" << cont + 1 << " is valid. \n";
0059 
0060       //getting the MultiRecHit collection and the trajectory with a first fit-smooth round
0061       std::pair<TransientTrackingRecHit::RecHitContainer, TrajectoryStateOnSurface> hits =
0062           collectHits(*BeforeDAFTraj, measurementCollector, &*measTk);
0063 
0064       //new initial fit
0065       CurrentTraj = fit(hits, theFitter, *BeforeDAFTraj);
0066 
0067       //starting the annealing program
0068       for (std::vector<double>::const_iterator ian = updator->getAnnealingProgram().begin();
0069            ian != updator->getAnnealingProgram().end();
0070            ian++) {
0071         if (CurrentTraj.isValid()) {
0072           LogDebug("DAFTrackProducerAlgorithm") << "Seed direction is " << CurrentTraj.seed().direction()
0073                                                 << ".Traj direction is " << CurrentTraj.direction() << std::endl;
0074 
0075           //updating MultiRecHits and fit-smooth again
0076           std::pair<TransientTrackingRecHit::RecHitContainer, TrajectoryStateOnSurface> curiterationhits =
0077               updateHits(CurrentTraj, updator, &*measTk, *ian);
0078           if (curiterationhits.first.size() < 3) {
0079             LogDebug("DAFTrackProducerAlgorithm")
0080                 << "Rejecting trajectory with " << curiterationhits.first.size() << " hits" << std::endl;
0081             CurrentTraj = Trajectory();
0082             break;
0083           }
0084 
0085           CurrentTraj = fit(curiterationhits, theFitter, CurrentTraj);
0086 
0087           //saving trajectory for each annealing cycle ...
0088           if (TrajAnnSaving_) {
0089             TrajAnnealing temp(CurrentTraj, *ian);
0090             trajann.push_back(temp);
0091           }
0092 
0093           LogDebug("DAFTrackProducerAlgorithm") << "done annealing value " << (*ian);
0094 
0095         } else
0096           break;
0097       }  //end annealing cycle
0098 
0099       int percOfHitsUnchangedAfterDAF =
0100           (1. * checkHits(*BeforeDAFTraj, CurrentTraj) / (1. * BeforeDAFTrajMeas.size())) * 100.;
0101       LogDebug("DAFTrackProducerAlgorithm")
0102           << "Ended annealing program with " << percOfHitsUnchangedAfterDAF << " unchanged." << std::endl;
0103 
0104       //computing the ndof keeping into account the weights
0105       ndof = calculateNdof(CurrentTraj);
0106 
0107       //checking if the trajectory has the minimum number of valid hits ( weight (>1e-6) )
0108       //in order to remove tracks with too many outliers.
0109       int goodHits = countingGoodHits(CurrentTraj);
0110 
0111       if (goodHits >= minHits_) {
0112         bool ok = buildTrack(CurrentTraj, algoResults, ndof, bs, &BeforeDAFTrack);
0113         // or filtered?
0114         if (ok)
0115           cont++;
0116 
0117         //saving tracks before and after DAF
0118         if ((100. - percOfHitsUnchangedAfterDAF) > 0.) {
0119           bool okBefore = buildTrack(*BeforeDAFTraj, algoResultsBeforeDAF, ndof, bs, &BeforeDAFTrack);
0120           bool okAfter = buildTrack(CurrentTraj, algoResultsAfterDAF, ndof, bs, &BeforeDAFTrack);
0121           if (okBefore && okAfter)
0122             nTracksChanged++;
0123         }
0124       } else {
0125         LogDebug("DAFTrackProducerAlgorithm") << "Rejecting trajectory with " << CurrentTraj.foundHits() << " hits";
0126       }
0127     }  //end run on track collection
0128     else {
0129       LogDebug("DAFTrackProducerAlgorithm") << "Rejecting empty trajectory" << std::endl;
0130     }
0131   }  //end run on track collection
0132 
0133   LogDebug("DAFTrackProducerAlgorithm") << "Number of Tracks found:   " << cont << "\n";
0134   LogDebug("DAFTrackProducerAlgorithm") << "Number of Tracks changed: " << nTracksChanged << "\n";
0135 }
0136 /*------------------------------------------------------------------------------------------------------*/
0137 std::pair<TransientTrackingRecHit::RecHitContainer, TrajectoryStateOnSurface> DAFTrackProducerAlgorithm::collectHits(
0138     const Trajectory vtraj,
0139     const MultiRecHitCollector* measurementCollector,
0140     const MeasurementTrackerEvent* measTk) const {
0141   LogDebug("DAFTrackProducerAlgorithm") << "Calling DAFTrackProducerAlgorithm::collectHits";
0142 
0143   //getting the traj measurements from the MeasurementCollector
0144   TransientTrackingRecHit::RecHitContainer hits;
0145   std::vector<TrajectoryMeasurement> collectedmeas = measurementCollector->recHits(vtraj, measTk);
0146 
0147   //if the MeasurementCollector is empty, make an "empty" pair
0148   //else taking the collected measured hits and building the pair
0149   if (collectedmeas.empty())
0150     return std::make_pair(TransientTrackingRecHit::RecHitContainer(), TrajectoryStateOnSurface());
0151 
0152   for (std::vector<TrajectoryMeasurement>::const_iterator iter = collectedmeas.begin(); iter != collectedmeas.end();
0153        iter++) {
0154     hits.push_back(iter->recHit());
0155   }
0156 
0157   //TrajectoryStateWithArbitraryError() == Creates a TrajectoryState with the same parameters
0158   //     as the input one, but with "infinite" errors, i.e. errors so big that they don't
0159   //     bias a fit starting with this state.
0160   //return std::make_pair(hits,TrajectoryStateWithArbitraryError()(collectedmeas.front().predictedState()));
0161 
0162   // I do not have to rescale the error because it is already rescaled in the fit code
0163   TrajectoryStateOnSurface initialStateFromTrack = collectedmeas.front().predictedState();
0164 
0165   LogDebug("DAFTrackProducerAlgorithm")
0166       << "Pair (hits, TSOS)  with TSOS predicted(collectedmeas.front().predictedState())";
0167   return std::make_pair(hits, initialStateFromTrack);
0168 }
0169 /*------------------------------------------------------------------------------------------------------*/
0170 std::pair<TransientTrackingRecHit::RecHitContainer, TrajectoryStateOnSurface> DAFTrackProducerAlgorithm::updateHits(
0171     const Trajectory vtraj,
0172     const SiTrackerMultiRecHitUpdator* updator,
0173     const MeasurementTrackerEvent* theMTE,
0174     double annealing) const {
0175   LogDebug("DAFTrackProducerAlgorithm") << "Calling DAFTrackProducerAlgorithm::updateHits";
0176   TransientTrackingRecHit::RecHitContainer hits;
0177   std::vector<TrajectoryMeasurement> vmeas = vtraj.measurements();
0178   std::vector<TrajectoryMeasurement>::reverse_iterator imeas;
0179   unsigned int hitcounter = 1;
0180 
0181   //I run inversely on the trajectory obtained and update the state
0182   for (imeas = vmeas.rbegin(); imeas != vmeas.rend(); imeas++, hitcounter++) {
0183     DetId id = imeas->recHit()->geographicalId();
0184     MeasurementDetWithData measDet = theMTE->idToDet(id);
0185 
0186     TrajectoryStateCombiner combiner;
0187     TrajectoryStateOnSurface combtsos;
0188     if (hitcounter == vmeas.size())
0189       combtsos = imeas->predictedState();  //fwd
0190     else if (hitcounter == 1)
0191       combtsos = imeas->backwardPredictedState();  //bwd
0192     else
0193       combtsos = combiner(imeas->backwardPredictedState(), imeas->predictedState());
0194 
0195     PrintHit(&*imeas->recHit(), combtsos);
0196     if (imeas->recHit()->isValid()) {
0197       TransientTrackingRecHit::RecHitPointer updated = updator->update(imeas->recHit(), combtsos, measDet, annealing);
0198       hits.push_back(updated);
0199     } else {
0200       hits.push_back(imeas->recHit());
0201     }
0202   }
0203 
0204   TrajectoryStateOnSurface updatedStateFromTrack = vmeas.back().predictedState();
0205 
0206   //return std::make_pair(hits,TrajectoryStateWithArbitraryError()(vmeas.back().updatedState()));
0207   LogDebug("DAFTrackProducerAlgorithm") << "Pair (hits, TSOS)  with TSOS predicted (vmeas.back().predictedState())";
0208 
0209   return std::make_pair(hits, updatedStateFromTrack);
0210 }
0211 /*------------------------------------------------------------------------------------------------------*/
0212 Trajectory DAFTrackProducerAlgorithm::fit(
0213     const std::pair<TransientTrackingRecHit::RecHitContainer, TrajectoryStateOnSurface>& hits,
0214     const TrajectoryFitter* theFitter,
0215     Trajectory vtraj) const {
0216   //creating a new trajectory starting from the direction of the seed of the input one and the hits
0217   Trajectory newVec = theFitter->fitOne(TrajectorySeed({}, {}, vtraj.seed().direction()), hits.first, hits.second);
0218 
0219   if (newVec.isValid())
0220     return newVec;
0221   else {
0222     LogDebug("DAFTrackProducerAlgorithm") << "Fit no valid.";
0223     return Trajectory();
0224   }
0225 }
0226 /*------------------------------------------------------------------------------------------------------*/
0227 bool DAFTrackProducerAlgorithm::buildTrack(const Trajectory vtraj,
0228                                            AlgoProductCollection& algoResults,
0229                                            float ndof,
0230                                            const reco::BeamSpot& bs,
0231                                            const reco::TrackRef* BeforeDAFTrack) const {
0232   LogDebug("DAFTrackProducerAlgorithm") << " BUILDER " << std::endl;
0233   ;
0234   TrajectoryStateOnSurface innertsos;
0235 
0236   if (vtraj.isValid()) {
0237     std::unique_ptr<Trajectory> theTraj(new Trajectory(vtraj));
0238 
0239     if (vtraj.direction() == alongMomentum) {
0240       //if (theTraj->direction() == oppositeToMomentum) {
0241       innertsos = vtraj.firstMeasurement().updatedState();
0242     } else {
0243       innertsos = vtraj.lastMeasurement().updatedState();
0244     }
0245 
0246     TSCBLBuilderNoMaterial tscblBuilder;
0247     TrajectoryStateClosestToBeamLine tscbl = tscblBuilder(*(innertsos.freeState()), bs);
0248 
0249     if (tscbl.isValid() == false)
0250       return false;
0251 
0252     GlobalPoint v = tscbl.trackStateAtPCA().position();
0253     math::XYZPoint pos(v.x(), v.y(), v.z());
0254     GlobalVector p = tscbl.trackStateAtPCA().momentum();
0255     math::XYZVector mom(p.x(), p.y(), p.z());
0256 
0257     //    LogDebug("TrackProducer") <<v<<p<<std::endl;
0258 
0259     auto algo = (*BeforeDAFTrack)->algo();
0260     std::unique_ptr<reco::Track> theTrack(new reco::Track(vtraj.chiSquared(),
0261                                                           ndof,  //in the DAF the ndof is not-integer
0262                                                           pos,
0263                                                           mom,
0264                                                           tscbl.trackStateAtPCA().charge(),
0265                                                           tscbl.trackStateAtPCA().curvilinearError(),
0266                                                           algo));
0267     theTrack->setQualityMask((*BeforeDAFTrack)->qualityMask());
0268 
0269     AlgoProduct aProduct{theTraj.release(), theTrack.release(), vtraj.direction(), 0};
0270     algoResults.push_back(aProduct);
0271 
0272     return true;
0273   } else {
0274     LogDebug("DAFTrackProducerAlgorithm") << " BUILDER NOT POSSIBLE: traj is not valid" << std::endl;
0275     ;
0276     return false;
0277   }
0278 }
0279 /*------------------------------------------------------------------------------------------------------*/
0280 int DAFTrackProducerAlgorithm::countingGoodHits(const Trajectory traj) const {
0281   int ngoodhits = 0;
0282   std::vector<TrajectoryMeasurement> vtm = traj.measurements();
0283 
0284   for (std::vector<TrajectoryMeasurement>::const_iterator tm = vtm.begin(); tm != vtm.end(); tm++) {
0285     //if the rechit is valid
0286     if (tm->recHit()->isValid()) {
0287       SiTrackerMultiRecHit const& mHit = dynamic_cast<SiTrackerMultiRecHit const&>(*tm->recHit());
0288       std::vector<const TrackingRecHit*> components = mHit.recHits();
0289 
0290       int iComp = 0;
0291 
0292       for (std::vector<const TrackingRecHit*>::const_iterator iter = components.begin(); iter != components.end();
0293            iter++, iComp++) {
0294         //if there is at least one component with weight higher than 1e-6 then the hit is not an outlier
0295         if (mHit.weight(iComp) > 1e-6) {
0296           ngoodhits++;
0297           iComp++;
0298           break;
0299         }
0300       }
0301     }
0302   }
0303 
0304   LogDebug("DAFTrackProducerAlgorithm") << "Original number of valid hits " << traj.foundHits()
0305                                         << " -> hit with good weight (>1e-6) are " << ngoodhits;
0306   return ngoodhits;
0307 }
0308 /*------------------------------------------------------------------------------------------------------*/
0309 
0310 void DAFTrackProducerAlgorithm::filter(const TrajectoryFitter* fitter,
0311                                        std::vector<Trajectory>& input,
0312                                        int minhits,
0313                                        std::vector<Trajectory>& output,
0314                                        const TransientTrackingRecHitBuilder* builder) const {
0315   if (input.empty())
0316     return;
0317 
0318   int ngoodhits = 0;
0319   std::vector<TrajectoryMeasurement> vtm = input[0].measurements();
0320   TransientTrackingRecHit::RecHitContainer hits;
0321 
0322   //count the number of non-outlier and non-invalid hits
0323   for (std::vector<TrajectoryMeasurement>::reverse_iterator tm = vtm.rbegin(); tm != vtm.rend(); tm++) {
0324     //if the rechit is valid
0325     if (tm->recHit()->isValid()) {
0326       SiTrackerMultiRecHit const& mHit = dynamic_cast<SiTrackerMultiRecHit const&>(*tm->recHit());
0327       std::vector<const TrackingRecHit*> components = mHit.recHits();
0328       int iComp = 0;
0329       bool isGood = false;
0330       for (std::vector<const TrackingRecHit*>::const_iterator iter = components.begin(); iter != components.end();
0331            iter++, iComp++) {
0332         //if there is at least one component with weight higher than 1e-6 then the hit is not an outlier
0333         if (mHit.weight(iComp) > 1e-6) {
0334           ngoodhits++;
0335           iComp++;
0336           isGood = true;
0337           break;
0338         }
0339       }
0340       if (isGood) {
0341         TkClonerImpl hc = static_cast<TkTransientTrackingRecHitBuilder const*>(builder)->cloner();
0342         auto tempHit = hc.makeShared(tm->recHit(), tm->updatedState());
0343         hits.push_back(tempHit);
0344       } else
0345         hits.push_back(std::make_shared<InvalidTrackingRecHit>(*tm->recHit()->det(), TrackingRecHit::missing));
0346     } else {
0347       TkClonerImpl hc = static_cast<TkTransientTrackingRecHitBuilder const*>(builder)->cloner();
0348       auto tempHit = hc.makeShared(tm->recHit(), tm->updatedState());
0349       hits.push_back(tempHit);
0350     }
0351   }
0352 
0353   LogDebug("DAFTrackProducerAlgorithm") << "Original number of valid hits " << input[0].foundHits()
0354                                         << "; after filtering " << ngoodhits;
0355   if (ngoodhits > input[0].foundHits())
0356     edm::LogError("DAFTrackProducerAlgorithm")
0357         << "Something wrong: the number of good hits from DAFTrackProducerAlgorithm::filter " << ngoodhits
0358         << " is higher than the original one " << input[0].foundHits();
0359 
0360   if (ngoodhits < minhits)
0361     return;
0362 
0363   TrajectoryStateOnSurface curstartingTSOS = input.front().lastMeasurement().updatedState();
0364   LogDebug("DAFTrackProducerAlgorithm") << "starting tsos for final refitting " << curstartingTSOS;
0365   //curstartingTSOS.rescaleError(100);
0366 
0367   output = fitter->fit(TrajectorySeed({}, {}, input.front().seed().direction()),
0368                        hits,
0369                        TrajectoryStateWithArbitraryError()(curstartingTSOS));
0370 
0371   LogDebug("DAFTrackProducerAlgorithm") << "After filtering " << output.size() << " trajectories";
0372 }
0373 /*------------------------------------------------------------------------------------------------------*/
0374 float DAFTrackProducerAlgorithm::calculateNdof(const Trajectory vtraj) const {
0375   if (!vtraj.isValid())
0376     return 0;
0377   float ndof = 0;
0378   const std::vector<TrajectoryMeasurement>& meas = vtraj.measurements();
0379   for (std::vector<TrajectoryMeasurement>::const_iterator iter = meas.begin(); iter != meas.end(); iter++) {
0380     if (iter->recHit()->isValid()) {
0381       SiTrackerMultiRecHit const& mHit = dynamic_cast<SiTrackerMultiRecHit const&>(*iter->recHit());
0382       std::vector<const TrackingRecHit*> components = mHit.recHits();
0383       int iComp = 0;
0384       for (std::vector<const TrackingRecHit*>::const_iterator iter2 = components.begin(); iter2 != components.end();
0385            iter2++, iComp++) {
0386         if ((*iter2)->isValid())
0387           ndof += ((*iter2)->dimension()) * mHit.weight(iComp);
0388       }
0389     }
0390   }
0391 
0392   return ndof - 5;
0393 }
0394 //------------------------------------------------------------------------------------------------
0395 int DAFTrackProducerAlgorithm::checkHits(Trajectory iInitTraj, const Trajectory iFinalTraj) const {
0396   std::vector<TrajectoryMeasurement> initmeasurements = iInitTraj.measurements();
0397   std::vector<TrajectoryMeasurement> finalmeasurements = iFinalTraj.measurements();
0398   std::vector<TrajectoryMeasurement>::iterator jmeas;
0399   int nSame = 0;
0400   int ihit = 0;
0401 
0402   if (initmeasurements.empty() || finalmeasurements.empty()) {
0403     LogDebug("DAFTrackProducerAlgorithm") << "Initial or Final Trajectory empty.";
0404     return 0;
0405   }
0406 
0407   if (initmeasurements.size() != finalmeasurements.size()) {
0408     LogDebug("DAFTrackProducerAlgorithm")
0409         << "Initial Trajectory size(" << initmeasurements.size() << " hits) "
0410         << "is different to final traj size (" << finalmeasurements.size() << ")! No checkHits possible! ";
0411     return 0;
0412   }
0413 
0414   for (jmeas = initmeasurements.begin(); jmeas != initmeasurements.end(); jmeas++) {
0415     const TrackingRecHit* initHit = jmeas->recHit()->hit();
0416     if (!initHit->isValid() && ihit == 0)
0417       continue;
0418 
0419     if (initHit->isValid()) {
0420       TrajectoryMeasurement imeas = finalmeasurements.at(ihit);
0421       const TrackingRecHit* finalHit = imeas.recHit()->hit();
0422       const TrackingRecHit* MaxWeightHit = nullptr;
0423       float maxweight = 0;
0424 
0425       const SiTrackerMultiRecHit* mrh = dynamic_cast<const SiTrackerMultiRecHit*>(finalHit);
0426       if (mrh) {
0427         std::vector<const TrackingRecHit*> components = mrh->recHits();
0428         std::vector<const TrackingRecHit*>::const_iterator icomp;
0429         int hitcounter = 0;
0430 
0431         for (icomp = components.begin(); icomp != components.end(); icomp++) {
0432           if ((*icomp)->isValid()) {
0433             double weight = mrh->weight(hitcounter);
0434             if (weight > maxweight) {
0435               MaxWeightHit = *icomp;
0436               maxweight = weight;
0437             }
0438           }
0439           hitcounter++;
0440         }
0441         if (!MaxWeightHit)
0442           continue;
0443 
0444         auto myref1 = reinterpret_cast<const BaseTrackerRecHit*>(initHit)->firstClusterRef();
0445         auto myref2 = reinterpret_cast<const BaseTrackerRecHit*>(MaxWeightHit)->firstClusterRef();
0446 
0447         if (myref1 == myref2) {
0448           nSame++;
0449         } else {
0450           LogDebug("DAFTrackProducerAlgorithm") << "diverso hit!" << std::endl;
0451           TrajectoryStateOnSurface dummState;
0452           LogTrace("DAFTrackProducerAlgorithm") << "  This hit was:\n ";
0453           PrintHit(initHit, dummState);
0454           LogTrace("DAFTrackProducerAlgorithm") << "  instead now is:\n ";
0455           PrintHit(MaxWeightHit, dummState);
0456         }
0457       }
0458     } else {
0459       TrajectoryMeasurement imeas = finalmeasurements.at(ihit);
0460       const TrackingRecHit* finalHit = imeas.recHit()->hit();
0461       if (!finalHit->isValid()) {
0462         nSame++;
0463       }
0464     }
0465 
0466     ihit++;
0467   }
0468 
0469   return nSame;
0470 }
0471 
0472 void DAFTrackProducerAlgorithm::PrintHit(const TrackingRecHit* const& hit, TrajectoryStateOnSurface& tsos) const {
0473   if (hit->isValid()) {
0474     LogTrace("DAFTrackProducerAlgorithm")
0475         << "  Valid Hit with DetId " << hit->geographicalId().rawId() << " and dim:" << hit->dimension()
0476         << " local position " << hit->localPosition() << " global position " << hit->globalPosition() << " and r "
0477         << hit->globalPosition().perp();
0478     if (tsos.isValid())
0479       LogTrace("DAFTrackProducerAlgorithm") << "  TSOS combtsos " << tsos.localPosition();
0480 
0481   } else {
0482     LogTrace("DAFTrackProducerAlgorithm") << "  Invalid Hit with DetId " << hit->geographicalId().rawId();
0483   }
0484 }