CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_4/src/RecoVertex/AdaptiveVertexFit/src/AdaptiveVertexFitter.cc

Go to the documentation of this file.
00001 #include "RecoVertex/AdaptiveVertexFit/interface/AdaptiveVertexFitter.h"
00002 #include "DataFormats/GeometryCommonDetAlgo/interface/GlobalError.h"
00003 #include "RecoVertex/VertexTools/interface/AnnealingSchedule.h"
00004 #include "RecoVertex/VertexTools/interface/GeometricAnnealing.h"
00005 #include "RecoVertex/VertexTools/interface/VertexTrackFactory.h"
00006 #include "RecoVertex/VertexPrimitives/interface/VertexException.h"
00007 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00008 
00009 #include <algorithm>
00010 
00011 using namespace edm;
00012 
00013 // #define STORE_WEIGHTS
00014 #ifdef STORE_WEIGHTS
00015 #include <dataharvester/Writer.h>
00016 #endif
00017 
00018 using namespace std;
00019 
00020 namespace {
00021   struct CompareTwoTracks {
00022     int operator() ( const reco::TransientTrack & a, const reco::TransientTrack & b ) {
00023             return ( a.impactPointState().globalMomentum().perp() >
00024                      b.impactPointState().globalMomentum().perp() ) ;
00025     };
00026   };
00027   // typedef ReferenceCountingPointer<VertexTrack<5> > RefCountedVertexTrack;
00028   typedef AdaptiveVertexFitter::RefCountedVertexTrack RefCountedVertexTrack;
00029 
00030   GlobalError fitError()
00031   {
00032     static GlobalError err;
00033     static bool once=true;
00034     if (once)
00035     {
00036       // that's how we model the lin pt error for the initial seed!
00037       static const float initialError = 10000;
00038       AlgebraicSymMatrix33 ret;
00039       ret(0,0)=initialError;
00040       ret(1,1)=initialError;
00041       ret(2,2)=initialError;
00042       err=GlobalError ( ret );
00043       once=false;
00044     }
00045     return err;
00046   }
00047 
00048   GlobalError linPointError()
00049   {
00050     static GlobalError err;
00051     static bool once=true;
00052     if (once)
00053     {
00054       // that's how we model the error of the linearization point.
00055       // for track weighting!
00056       AlgebraicSymMatrix33 ret;
00057       ret(0,0)=.3; ret(1,1)=.3; ret(2,2)=3.;
00058       // ret(0,0)=1e-7; ret(1,1)=1e-7; ret(2,2)=1e-7;
00059       err=GlobalError ( ret );
00060       once=false;
00061     }
00062     return err;
00063   }
00064 
00065   struct DistanceToRefPoint {
00066     DistanceToRefPoint ( const GlobalPoint & ref ) : theRef ( ref ) {}
00067 
00068     bool operator() ( const RefCountedVertexTrack & v1, const RefCountedVertexTrack & v2 )
00069     {
00070       return ( distance ( v1 ) < distance ( v2 ) );
00071     }
00072 
00073     float distance ( const RefCountedVertexTrack & v1 )
00074     {
00075       return ( v1->linearizedTrack()->track().initialFreeState().position() - theRef ).mag2();
00076     }
00077 
00078     private:
00079       GlobalPoint theRef;
00080   };
00081 
00082   #ifdef STORE_WEIGHTS
00083   map < RefCountedLinearizedTrackState, int > ids;
00084   int iter=0;
00085 
00086   int getId ( const RefCountedLinearizedTrackState & r )
00087   {
00088     static int ctr=1;
00089     if ( ids.count(r) == 0 )
00090     {
00091       ids[r]=ctr++;
00092     }
00093     return ids[r];
00094   }
00095   #endif
00096 }
00097 
00098 AdaptiveVertexFitter::AdaptiveVertexFitter(
00099       const AnnealingSchedule & ann,
00100       const LinearizationPointFinder & linP,
00101       const VertexUpdator<5> & updator,
00102       const VertexTrackCompatibilityEstimator<5> & crit,
00103       const VertexSmoother<5> & smoother,
00104       const AbstractLTSFactory<5> & ltsf ) :
00105     theNr(0),
00106     theLinP(linP.clone()), theUpdator( updator.clone()),
00107     theSmoother ( smoother.clone() ), theAssProbComputer( ann.clone() ),
00108     theComp ( crit.clone() ), theLinTrkFactory ( ltsf.clone() ),
00109     gsfIntermediarySmoothing_(false)
00110 {
00111   setParameters();
00112 }
00113 
00114 void AdaptiveVertexFitter::setWeightThreshold ( float w )
00115 {
00116   theWeightThreshold=w;
00117 }
00118 
00119 AdaptiveVertexFitter::AdaptiveVertexFitter
00120                         (const AdaptiveVertexFitter & o ) :
00121     theMaxShift ( o.theMaxShift ), theMaxLPShift ( o.theMaxLPShift ),
00122     theMaxStep ( o.theMaxStep ), theWeightThreshold ( o.theWeightThreshold ),
00123     theNr ( o.theNr ),
00124     theLinP ( o.theLinP->clone() ), theUpdator ( o.theUpdator->clone() ),
00125     theSmoother ( o.theSmoother->clone() ),
00126     theAssProbComputer ( o.theAssProbComputer->clone() ),
00127     theComp ( o.theComp->clone() ),
00128     theLinTrkFactory ( o.theLinTrkFactory->clone() ),
00129     gsfIntermediarySmoothing_(o.gsfIntermediarySmoothing_)
00130 {}
00131 
00132 AdaptiveVertexFitter::~AdaptiveVertexFitter()
00133 {
00134   delete theLinP;
00135   delete theUpdator;
00136   delete theSmoother;
00137   delete theAssProbComputer;
00138   delete theComp;
00139   delete theLinTrkFactory;
00140 }
00141 
00142 void AdaptiveVertexFitter::setParameters( double maxshift, double maxlpshift, 
00143                                           unsigned maxstep, double weightthreshold )
00144 {
00145   theMaxShift = maxshift;
00146   theMaxLPShift = maxlpshift;
00147   theMaxStep = maxstep;
00148   theWeightThreshold=weightthreshold;
00149 }
00150 
00151 
00152 void AdaptiveVertexFitter::setParameters ( const edm::ParameterSet & s )
00153 {
00154     setParameters ( s.getParameter<double>("maxshift"),
00155                     s.getParameter<double>("maxlpshift"),
00156                     s.getParameter<int>("maxstep"),
00157                     s.getParameter<double>("weightthreshold") );
00158 }
00159 
00160 CachingVertex<5>
00161 AdaptiveVertexFitter::vertex(const vector<reco::TransientTrack> & unstracks) const
00162 {
00163   if ( unstracks.size() < 2 )
00164   {
00165     LogError("RecoVertex|AdaptiveVertexFitter")
00166       << "Supplied fewer than two tracks. Vertex is invalid.";
00167     return CachingVertex<5>(); // return invalid vertex
00168   };
00169   vector < reco::TransientTrack > tracks = unstracks;
00170   sort ( tracks.begin(), tracks.end(), CompareTwoTracks() );
00171   // Linearization Point
00172   GlobalPoint linP = theLinP->getLinearizationPoint(tracks);
00173   // Initial vertex seed, with a very large error matrix
00174   VertexState lseed (linP, linPointError() );
00175   vector<RefCountedVertexTrack> vtContainer = linearizeTracks(tracks, lseed);
00176 
00177   VertexState seed (linP, fitError() );
00178   return fit(vtContainer, seed, false);
00179 }
00180 
00181 CachingVertex<5>
00182 AdaptiveVertexFitter::vertex(const vector<RefCountedVertexTrack> & tracks) const
00183 {
00184   if ( tracks.size() < 2 )
00185   {
00186     LogError("RecoVertex|AdaptiveVertexFitter")
00187       << "Supplied fewer than two tracks. Vertex is invalid.";
00188     return CachingVertex<5>(); // return invalid vertex
00189   };
00190   // Initial vertex seed, with a very small weight matrix
00191   GlobalPoint linP = tracks[0]->linearizedTrack()->linearizationPoint();
00192   VertexState seed (linP, fitError() );
00193   return fit(tracks, seed, false);
00194 }
00195 
00196 CachingVertex<5>
00197 AdaptiveVertexFitter::vertex(const vector<RefCountedVertexTrack> & tracks, const reco::BeamSpot & spot ) const
00198 {
00199   if ( tracks.size() < 1 )
00200   {
00201     LogError("RecoVertex|AdaptiveVertexFitter")
00202       << "Supplied no tracks. Vertex is invalid.";
00203     return CachingVertex<5>(); // return invalid vertex
00204   };
00205   VertexState beamSpotState(spot);
00206   return fit(tracks, beamSpotState, true );
00207 }
00208 
00209 
00210 
00214 CachingVertex<5>
00215 AdaptiveVertexFitter::vertex(const vector<reco::TransientTrack> & tracks,
00216                   const GlobalPoint& linPoint) const
00217 {
00218   if ( tracks.size() < 2 )
00219   {
00220     LogError("RecoVertex|AdaptiveVertexFitter")
00221       << "Supplied fewer than two tracks. Vertex is invalid.";
00222     return CachingVertex<5>(); // return invalid vertex
00223   };
00224   // Initial vertex seed, with a very large error matrix
00225   VertexState seed (linPoint, linPointError() );
00226   vector<RefCountedVertexTrack> vtContainer = linearizeTracks(tracks, seed);
00227   VertexState fitseed (linPoint, fitError() );
00228   return fit(vtContainer, fitseed, false);
00229 }
00230 
00231 
00236 CachingVertex<5> 
00237 AdaptiveVertexFitter::vertex(const vector<reco::TransientTrack> & unstracks,
00238                                const reco::BeamSpot& beamSpot) const
00239 {
00240   if ( unstracks.size() < 1 )
00241   {
00242     LogError("RecoVertex|AdaptiveVertexFitter")
00243       << "Supplied no tracks. Vertex is invalid.";
00244     return CachingVertex<5>(); // return invalid vertex
00245   };
00246 
00247   VertexState beamSpotState(beamSpot);
00248   vector<RefCountedVertexTrack> vtContainer;
00249 
00250   vector < reco::TransientTrack > tracks = unstracks;
00251   sort ( tracks.begin(), tracks.end(), CompareTwoTracks() );
00252 
00253   if (tracks.size() > 1) {
00254     // Linearization Point search if there are more than 1 track
00255     GlobalPoint linP = theLinP->getLinearizationPoint(tracks);
00256     VertexState lpState(linP, linPointError() );
00257     vtContainer = linearizeTracks(tracks, lpState);
00258   } else {
00259     // otherwise take the beamspot position.
00260     vtContainer = linearizeTracks(tracks, beamSpotState);
00261   }
00262 
00263   return fit(vtContainer, beamSpotState, true);
00264 }
00265 
00266 
00272 CachingVertex<5> AdaptiveVertexFitter::vertex(const vector<reco::TransientTrack> & tracks,
00273                   const GlobalPoint& priorPos,
00274                   const GlobalError& priorError) const
00275 
00276 {
00277   if ( tracks.size() < 1 )
00278   {
00279     LogError("RecoVertex|AdaptiveVertexFitter")
00280       << "Supplied no tracks. Vertex is invalid.";
00281     return CachingVertex<5>(); // return invalid vertex
00282   };
00283   VertexState seed (priorPos, priorError);
00284   vector<RefCountedVertexTrack> vtContainer = linearizeTracks(tracks, seed);
00285   return fit( vtContainer, seed, true );
00286 }
00287 
00288 
00293 CachingVertex<5> AdaptiveVertexFitter::vertex(
00294                 const vector<RefCountedVertexTrack> & tracks,
00295                   const GlobalPoint& priorPos,
00296                   const GlobalError& priorError) const
00297 {
00298   if ( tracks.size() < 1 )
00299   {
00300     LogError("RecoVertex|AdaptiveVertexFitter")
00301       << "Supplied no tracks. Vertex is invalid.";
00302     return CachingVertex<5>(); // return invalid vertex
00303   };
00304   VertexState seed (priorPos, priorError);
00305   return fit(tracks, seed, true);
00306 }
00307 
00308 
00316 vector<AdaptiveVertexFitter::RefCountedVertexTrack>
00317 AdaptiveVertexFitter::linearizeTracks(const vector<reco::TransientTrack> & tracks,
00318                                       const VertexState & seed ) const
00319 {
00320   const GlobalPoint & linP ( seed.position() );
00321   vector<RefCountedLinearizedTrackState> lTracks;
00322   for(vector<reco::TransientTrack>::const_iterator i = tracks.begin();
00323       i != tracks.end(); ++i )
00324   {
00325     try {
00326       RefCountedLinearizedTrackState lTrData
00327         = theLinTrkFactory->linearizedTrackState(linP, *i);
00328       lTracks.push_back(lTrData);
00329     } catch ( exception & e ) {
00330       LogInfo("RecoVertex/AdaptiveVertexFitter") 
00331         << "Exception " << e.what() << " in ::linearizeTracks."
00332         << "Your future vertex has just lost a track.";
00333     };
00334   }
00335   return weightTracks(lTracks, seed );
00336 }
00337 
00343 vector<AdaptiveVertexFitter::RefCountedVertexTrack>
00344 AdaptiveVertexFitter::reLinearizeTracks(
00345                                 const vector<RefCountedVertexTrack> & tracks,
00346                                 const CachingVertex<5> & vertex ) const
00347 {
00348   VertexState seed = vertex.vertexState();
00349   GlobalPoint linP = seed.position();
00350   vector<RefCountedLinearizedTrackState> lTracks;
00351   for(vector<RefCountedVertexTrack>::const_iterator i = tracks.begin();
00352     i != tracks.end(); i++)
00353   {
00354     try {
00355       RefCountedLinearizedTrackState lTrData
00356         = theLinTrkFactory->linearizedTrackState( linP, (**i).linearizedTrack()->track() );
00357       /*
00358       RefCountedLinearizedTrackState lTrData =
00359               (**i).linearizedTrack()->stateWithNewLinearizationPoint(linP);
00360               */
00361       lTracks.push_back(lTrData);
00362     } catch ( exception & e ) {
00363       LogInfo("RecoVertex/AdaptiveVertexFitter") 
00364         << "Exception " << e.what() << " in ::relinearizeTracks. "
00365         << "Will not relinearize this track.";
00366       lTracks.push_back ( (**i).linearizedTrack() );
00367     };
00368   };
00369   return reWeightTracks(lTracks, vertex );
00370 }
00371 
00372 AdaptiveVertexFitter * AdaptiveVertexFitter::clone() const
00373 {
00374   return new AdaptiveVertexFitter( * this );
00375 }
00376 
00377 double AdaptiveVertexFitter::getWeight ( float chi2 ) const
00378 {
00379   double weight = theAssProbComputer->weight(chi2);
00380 
00381   if ( weight > 1.0 )
00382   {
00383     LogInfo("RecoVertex/AdaptiveVertexFitter") << "Weight " << weight << " > 1.0!";
00384     weight=1.0;
00385   };
00386 
00387   if ( weight < 1e-20 )
00388   {
00389     // LogInfo("RecoVertex/AdaptiveVertexFitter") << "Weight " << weight << " < 0.0!";
00390     weight=1e-20;
00391   };
00392   return weight;
00393 }
00394 
00395 vector<AdaptiveVertexFitter::RefCountedVertexTrack>
00396 AdaptiveVertexFitter::reWeightTracks(
00397                     const vector<RefCountedLinearizedTrackState> & lTracks,
00398                     const CachingVertex<5> & vertex ) const
00399 {
00400   VertexState seed = vertex.vertexState();
00401   // cout << "[AdaptiveVertexFitter] now reweight around " << seed.position() << endl;
00402   theNr++;
00403   // GlobalPoint pos = seed.position();
00404 
00405   vector<RefCountedVertexTrack> finalTracks;
00406   VertexTrackFactory<5> vTrackFactory;
00407   #ifdef STORE_WEIGHTS
00408   iter++;
00409   #endif
00410   for(vector<RefCountedLinearizedTrackState>::const_iterator i
00411         = lTracks.begin(); i != lTracks.end(); i++)
00412   {
00413     double weight=0.;
00414     // cout << "[AdaptiveVertexFitter] estimate " << endl;
00415     pair < bool, double > chi2Res ( false, 0. );
00416     try {
00417       chi2Res =  theComp->estimate ( vertex, *i );
00418     } catch ( ... ) {};
00419     // cout << "[AdaptiveVertexFitter] /estimate " << endl;
00420     if (!chi2Res.first) {
00421       // cout << "[AdaptiveVertexFitter] aie... vertex candidate is at  " << vertex.position() << endl;
00422       LogInfo("AdaptiveVertexFitter" ) << "When reweighting, chi2<0. Will add this track with w=0.";
00423       // edm::LogInfo("AdaptiveVertexFitter" ) << "pt=" << (**i).track().pt();
00424     }else {
00425       weight = getWeight ( chi2Res.second );
00426     }
00427     
00428     RefCountedVertexTrack vTrData
00429        = vTrackFactory.vertexTrack(*i, seed, weight );
00430 
00431     #ifdef STORE_WEIGHTS
00432     map < string, dataharvester::MultiType > m;
00433     m["chi2"]=chi2;
00434     m["w"]=theAssProbComputer->weight(chi2);
00435     m["T"]=theAssProbComputer->currentTemp();
00436     m["n"]=iter;
00437     m["pos"]="reweight";
00438     m["id"]=getId ( *i );
00439     dataharvester::Writer::file("w.txt").save ( m );
00440     #endif
00441 
00442     finalTracks.push_back(vTrData);
00443   }
00444   sort ( finalTracks.begin(), finalTracks.end(), 
00445          DistanceToRefPoint ( vertex.position() ) );
00446   // cout << "[AdaptiveVertexFitter] /now reweight" << endl;
00447   return finalTracks;
00448 }
00449 
00450 vector<AdaptiveVertexFitter::RefCountedVertexTrack>
00451 AdaptiveVertexFitter::weightTracks(
00452                     const vector<RefCountedLinearizedTrackState> & lTracks,
00453                     const VertexState & seed ) const
00454 {
00455   theNr++;
00456   CachingVertex<5> seedvtx ( seed, vector<RefCountedVertexTrack> (), 0. );
00459   theAssProbComputer->resetAnnealing();
00460 
00461   vector<RefCountedVertexTrack> finalTracks;
00462   VertexTrackFactory<5> vTrackFactory;
00463   #ifdef STORE_WEIGHTS
00464   iter++;
00465   #endif
00466   for(vector<RefCountedLinearizedTrackState>::const_iterator i
00467         = lTracks.begin(); i != lTracks.end(); i++)
00468   {
00469 
00470     double weight = 0.;
00471     pair<bool, double> chi2Res = theComp->estimate ( seedvtx, *i );
00472     if (!chi2Res.first) {
00473       // cout << "[AdaptiveVertexFitter] Aiee! " << endl;
00474       LogInfo ("AdaptiveVertexFitter" ) << "When weighting a track, chi2 calculation failed;"
00475                                            << " will add with w=0.";
00476     } else {
00477       weight = getWeight ( chi2Res.second );
00478     }
00479     RefCountedVertexTrack vTrData
00480        = vTrackFactory.vertexTrack(*i, seed, weight );
00481     #ifdef STORE_WEIGHTS
00482     map < string, dataharvester::MultiType > m;
00483     m["chi2"]=chi2;
00484     m["w"]=theAssProbComputer->weight(chi2);
00485     m["T"]=theAssProbComputer->currentTemp();
00486     m["n"]=iter;
00487     m["id"]=getId ( *i );
00488     m["pos"]="weight";
00489     dataharvester::Writer::file("w.txt").save ( m );
00490     #endif
00491     finalTracks.push_back(vTrData);
00492   }
00493   return finalTracks;
00494 }
00495 
00501 vector<AdaptiveVertexFitter::RefCountedVertexTrack>
00502 AdaptiveVertexFitter::reWeightTracks(
00503                             const vector<RefCountedVertexTrack> & tracks,
00504                             const CachingVertex<5> & seed) const
00505 {
00506   vector<RefCountedLinearizedTrackState> lTracks;
00507   for(vector<RefCountedVertexTrack>::const_iterator i = tracks.begin();
00508     i != tracks.end(); i++)
00509   {
00510     lTracks.push_back((**i).linearizedTrack());
00511   }
00512 
00513   return reWeightTracks(lTracks, seed);
00514 }
00515 
00516 
00517 /*
00518  * The method where the vertex fit is actually done!
00519  */
00520 
00521 CachingVertex<5>
00522 AdaptiveVertexFitter::fit( const vector<RefCountedVertexTrack> & tracks,
00523                           const VertexState & priorSeed,
00524                           bool withPrior) const
00525 {
00526   // cout << "[AdaptiveVertexFit] fit with " << tracks.size() << endl;
00527   theAssProbComputer->resetAnnealing();
00528 
00529   vector<RefCountedVertexTrack> initialTracks;
00530   GlobalPoint priorVertexPosition = priorSeed.position();
00531   GlobalError priorVertexError = priorSeed.error();
00532 
00533   CachingVertex<5> returnVertex( priorVertexPosition,priorVertexError,
00534                               initialTracks,0);
00535   if (withPrior)
00536   {
00537     returnVertex = CachingVertex<5>(priorVertexPosition,priorVertexError,
00538                     priorVertexPosition,priorVertexError,initialTracks,0);
00539   }
00540 
00541   // vector<RefCountedVertexTrack> globalVTracks = tracks;
00542   // sort the tracks, according to distance to seed!
00543   vector<RefCountedVertexTrack> globalVTracks ( tracks.size() );
00544 
00545   partial_sort_copy ( tracks.begin(), tracks.end(), 
00546       globalVTracks.begin(), globalVTracks.end(), DistanceToRefPoint ( priorSeed.position() ) );
00547 
00548   // main loop through all the VTracks
00549   int lpStep = 0; int step = 0;
00550 
00551   CachingVertex<5> initialVertex = returnVertex;
00552 
00553   GlobalPoint newPosition = priorVertexPosition;
00554   GlobalPoint previousPosition = newPosition;
00555 
00556   int ns_trks=0; // number of significant tracks.
00557   // If we have only two significant tracks, we return an invalid vertex
00558 
00559   // cout << "[AdaptiveVertexFit] start " << tracks.size() << endl;
00560   /*
00561   for ( vector< RefCountedVertexTrack >::const_iterator 
00562         i=globalVTracks.begin(); i!=globalVTracks.end() ; ++i )
00563   {
00564     cout << "  " << (**i).linearizedTrack()->track().initialFreeState().momentum() << endl;
00565   }*/
00566   do {
00567     ns_trks=0;
00568     CachingVertex<5> fVertex = initialVertex;
00569     // cout << "[AdaptiveVertexFit] step " << step << " at " << fVertex.position() << endl;
00570     if ((previousPosition - newPosition).transverse() > theMaxLPShift)
00571     {
00572       // relinearize and reweight.
00573       // (reLinearizeTracks also reweights tracks)
00574       // cout << "[AdaptiveVertexFit] relinearize at " << returnVertex.position() << endl;
00575       if (gsfIntermediarySmoothing_) returnVertex = theSmoother->smooth(returnVertex);
00576       globalVTracks = reLinearizeTracks( globalVTracks, returnVertex );
00577       lpStep++;
00578     } else if (step) {
00579       // reweight, if it is not the first step
00580       // cout << "[AdaptiveVertexFit] reweight at " << returnVertex.position() << endl;
00581       if (gsfIntermediarySmoothing_) returnVertex = theSmoother->smooth(returnVertex);
00582       globalVTracks = reWeightTracks( globalVTracks, returnVertex );
00583     }
00584     // cout << "[AdaptiveVertexFit] relinarized, reweighted" << endl;
00585     // update sequentially the vertex estimate
00586     CachingVertex<5> nVertex;
00587     for(vector<RefCountedVertexTrack>::const_iterator i
00588           = globalVTracks.begin(); i != globalVTracks.end(); i++)
00589     {
00590       if ((**i).weight() > 0.) nVertex = theUpdator->add( fVertex, *i );
00591         else  nVertex = fVertex;
00592       if (nVertex.isValid()) {
00593         if ( (**i).weight() >= theWeightThreshold )
00594         {
00595           ns_trks++;
00596         };
00597 
00598         if ( fabs ( nVertex.position().z() ) > 10000. ||
00599              nVertex.position().perp()>120.)
00600         {
00601           // were more than 100 m off!!
00602           LogInfo ("AdaptiveVertexFitter" ) << "Vertex candidate just took off to " << nVertex.position()
00603                                             << "! Will discard this update!";
00604 //          //<< "track pt was " << (**i).linearizedTrack()->track().pt()
00605 //                                           << "track momentum was " << (**i).linearizedTrack()->track().initialFreeState().momentum()
00606 //                                           << "track position was " << (**i).linearizedTrack()->track().initialFreeState().position()
00607 //                                           << "track chi2 was " << (**i).linearizedTrack()->track().chi2()
00608 //                                           << "track ndof was " << (**i).linearizedTrack()->track().ndof()
00609 //                                           << "track w was " << (**i).weight()
00610 //                                           << "track schi2 was " << (**i).smoothedChi2();
00611         } else {
00612                 fVertex = nVertex;
00613         }
00614       } else {
00615         LogInfo("RecoVertex/AdaptiveVertexFitter") 
00616           << "The updator returned an invalid vertex when adding track "
00617           << i-globalVTracks.begin() 
00618                 << ".\n Your vertex might just have lost one good track.";
00619       };
00620     }
00621     previousPosition = newPosition;
00622     newPosition = fVertex.position();
00623     returnVertex = fVertex;
00624     theAssProbComputer->anneal();
00625     step++;
00626     if ( step >= theMaxStep ) break;
00627 
00628   } while (
00629       // repeat as long as
00630       // - vertex moved too much or
00631       // - we're not yet annealed
00632          ( ((previousPosition - newPosition).mag() > theMaxShift) ||
00633            (!(theAssProbComputer->isAnnealed()) ) ) ) ;
00634 
00635   if ( theWeightThreshold > 0. &&  ns_trks < 2 && !withPrior ) 
00636   {
00637     LogDebug("AdaptiveVertexFitter") 
00638       << "fewer than two significant tracks (w>" << theWeightThreshold << ")."
00639       << " Fitted vertex is invalid.";
00640     return CachingVertex<5>(); // return invalid vertex
00641   }
00642 
00643   #ifdef STORE_WEIGHTS
00644   map < string, dataharvester::MultiType > m;
00645   m["chi2"]=chi2;
00646   m["w"]=theAssProbComputer->weight(chi2);
00647   m["T"]=theAssProbComputer->currentTemp();
00648   m["n"]=iter;
00649   m["id"]=getId ( *i );
00650   m["pos"]="final";
00651   dataharvester::Writer::file("w.txt").save ( m );
00652   #endif
00653   // cout << "[AdaptiveVertexFit] /fit" << endl;
00654   return theSmoother->smooth( returnVertex );
00655 }