CMS 3D CMS Logo

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