CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_2_SLHC2/src/RecoLocalMuon/DTSegment/src/DTSegmentUpdator.cc

Go to the documentation of this file.
00001 
00012 /* This Class Header */
00013 #include "RecoLocalMuon/DTSegment/src/DTSegmentUpdator.h"
00014 
00015 /* Collaborating Class Header */
00016 
00017 //mene
00018 #include "DataFormats/DTRecHit/interface/DTChamberRecSegment2D.h"
00019 
00020 #include "DataFormats/DTRecHit/interface/DTRecSegment2D.h"
00021 #include "DataFormats/DTRecHit/interface/DTSLRecSegment2D.h"
00022 #include "DataFormats/DTRecHit/interface/DTChamberRecSegment2D.h"
00023 #include "DataFormats/DTRecHit/interface/DTRecSegment4D.h"
00024 #include "DataFormats/DTRecHit/interface/DTRecHit1D.h"
00025 #include "DataFormats/GeometryCommonDetAlgo/interface/ErrorFrameTransformer.h"
00026 
00027 #include "RecoLocalMuon/DTSegment/src/DTSegmentCand.h"
00028 #include "RecoLocalMuon/DTRecHit/interface/DTRecHitAlgoFactory.h"
00029 #include "RecoLocalMuon/DTSegment/src/DTLinearFit.h"
00030 
00031 #include "Geometry/Records/interface/MuonGeometryRecord.h"
00032 #include "Geometry/DTGeometry/interface/DTGeometry.h"
00033 #include "Geometry/DTGeometry/interface/DTLayer.h"
00034 
00035 #include "FWCore/Utilities/interface/Exception.h"
00036 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00037 #include "FWCore/Framework/interface/EventSetup.h"
00038 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00039 
00040 /* C++ Headers */
00041 #include <string>
00042 
00043 using namespace std;
00044 using namespace edm;
00045 
00046 /* ====================================================================== */
00047 
00049 DTSegmentUpdator::DTSegmentUpdator(const ParameterSet& config) :
00050   theFitter(new DTLinearFit()) ,
00051   vdrift_4parfit(config.getParameter<bool>("performT0_vdriftSegCorrection")),
00052   T0_hit_resolution(config.getParameter<double>("hit_afterT0_resolution")),
00053   perform_delta_rejecting(config.getParameter<bool>("perform_delta_rejecting")),
00054   debug(config.getUntrackedParameter<bool>("debug",false)) 
00055 {  
00056   string theAlgoName = config.getParameter<string>("recAlgo");
00057   theAlgo = DTRecHitAlgoFactory::get()->create(theAlgoName, 
00058                                                config.getParameter<ParameterSet>("recAlgoConfig"));
00059 
00060   if(debug)
00061     cout << "[DTSegmentUpdator] Constructor called" << endl;
00062 }
00063 
00065 DTSegmentUpdator::~DTSegmentUpdator() {
00066   delete theFitter;
00067 }
00068 
00069 /* Operations */ 
00070 
00071 void DTSegmentUpdator::setES(const EventSetup& setup){
00072   setup.get<MuonGeometryRecord>().get(theGeom);
00073   theAlgo->setES(setup);
00074 }
00075 
00076 void DTSegmentUpdator::update(DTRecSegment4D* seg, const bool calcT0) const {
00077 
00078   if(debug) cout << "[DTSegmentUpdator] Starting to update the segment" << endl;
00079 
00080   const bool hasPhi = seg->hasPhi();
00081   const bool hasZed = seg->hasZed();
00082 
00083   //reject the bad hits (due to delta rays)
00084   if(perform_delta_rejecting && hasPhi) rejectBadHits(seg->phiSegment());
00085 
00086   int step = (hasPhi && hasZed) ? 3 : 2;
00087   if(calcT0) step = 4;
00088 
00089   GlobalPoint  pos = theGeom->idToDet(seg->geographicalId())->toGlobal(seg->localPosition());
00090   GlobalVector dir = theGeom->idToDet(seg->geographicalId())->toGlobal(seg->localDirection());
00091 
00092   if(calcT0) calculateT0corr(seg);
00093 
00094   if(hasPhi) updateHits(seg->phiSegment(),pos,dir,step);
00095   if(hasZed) updateHits(seg->zSegment()  ,pos,dir,step);
00096 
00097   fit(seg);
00098 }
00099 
00100 void DTSegmentUpdator::update(DTRecSegment2D* seg) const {
00101   GlobalPoint pos = (theGeom->idToDet(seg->geographicalId()))->toGlobal(seg->localPosition());
00102   GlobalVector dir = (theGeom->idToDet(seg->geographicalId()))->toGlobal(seg->localDirection());
00103 
00104   updateHits(seg,pos,dir);
00105   fit(seg);
00106 }
00107 
00108 void DTSegmentUpdator::fit(DTRecSegment4D* seg)  const {
00109   // after the update must refit the segments
00110   if(seg->hasPhi()) fit(seg->phiSegment());
00111   if(seg->hasZed()) fit(seg->zSegment());
00112 
00113   const DTChamber* theChamber = theGeom->chamber(seg->chamberId());
00114 
00115   if(seg->hasPhi() && seg->hasZed() ) {
00116 
00117     DTChamberRecSegment2D *segPhi=seg->phiSegment();
00118     DTSLRecSegment2D *segZed=seg->zSegment();
00119 
00120     // NB Phi seg is already in chamber ref
00121     LocalPoint posPhiInCh = segPhi->localPosition();
00122     LocalVector dirPhiInCh= segPhi->localDirection();
00123 
00124     // Zed seg is in SL one
00125     const DTSuperLayer* zSL = theChamber->superLayer(segZed->superLayerId());
00126     LocalPoint zPos(segZed->localPosition().x(), 
00127                     (zSL->toLocal(theChamber->toGlobal(segPhi->localPosition()))).y(),
00128                     0.);
00129 
00130     LocalPoint posZInCh = theChamber->toLocal(zSL->toGlobal(zPos));
00131 
00132     LocalVector dirZInCh = theChamber->toLocal(zSL->toGlobal(segZed->localDirection()));
00133 
00134     LocalPoint posZAt0 = posZInCh + dirZInCh*(-posZInCh.z())/cos(dirZInCh.theta());
00135 
00136     // given the actual definition of chamber refFrame, (with z poiniting to IP),
00137     // the zed component of direction is negative.
00138     LocalVector dir=LocalVector(dirPhiInCh.x()/fabs(dirPhiInCh.z()),dirZInCh.y()/fabs(dirZInCh.z()),-1.);
00139 
00140     seg->setPosition(LocalPoint(posPhiInCh.x(),posZAt0.y(),0.));
00141     seg->setDirection(dir.unit());
00142 
00143     AlgebraicSymMatrix mat(4);
00144 
00145     // set cov matrix
00146     mat[0][0] = segPhi->parametersError()[0][0]; //sigma (dx/dz)
00147     mat[0][2] = segPhi->parametersError()[0][1]; //cov(dx/dz,x)
00148     mat[2][2] = segPhi->parametersError()[1][1]; //sigma (x)
00149 
00150     seg->setCovMatrix(mat);
00151     seg->setCovMatrixForZed(posZInCh);
00152 
00153   }
00154 
00155   else if (seg->hasPhi()) {
00156     DTChamberRecSegment2D *segPhi=seg->phiSegment();
00157 
00158     seg->setPosition(segPhi->localPosition());
00159     seg->setDirection(segPhi->localDirection());
00160 
00161     AlgebraicSymMatrix mat(4);
00162     // set cov matrix
00163     mat[0][0] = segPhi->parametersError()[0][0]; //sigma (dx/dz)
00164     mat[0][2] = segPhi->parametersError()[0][1]; //cov(dx/dz,x)
00165     mat[2][2] = segPhi->parametersError()[1][1]; //sigma (x)
00166 
00167     seg->setCovMatrix(mat);
00168   }
00169 
00170   else if (seg->hasZed()) {
00171     DTSLRecSegment2D *segZed = seg->zSegment();
00172 
00173     // Zed seg is in SL one
00174     GlobalPoint glbPosZ = ( theGeom->superLayer(segZed->superLayerId()) )->toGlobal(segZed->localPosition());
00175     LocalPoint posZInCh = ( theGeom->chamber(segZed->superLayerId().chamberId()) )->toLocal(glbPosZ);
00176 
00177     GlobalVector glbDirZ = (theGeom->superLayer(segZed->superLayerId()) )->toGlobal(segZed->localDirection());
00178     LocalVector dirZInCh = (theGeom->chamber(segZed->superLayerId().chamberId()) )->toLocal(glbDirZ);
00179 
00180     LocalPoint posZAt0 = posZInCh+
00181       dirZInCh*(-posZInCh.z())/cos(dirZInCh.theta());
00182 
00183     seg->setPosition(posZAt0);
00184     seg->setDirection(dirZInCh);
00185 
00186     AlgebraicSymMatrix mat(4);
00187     // set cov matrix
00188     seg->setCovMatrix(mat);
00189     seg->setCovMatrixForZed(posZInCh);
00190   }
00191 }
00192 
00193 bool DTSegmentUpdator::fit(DTSegmentCand* seg) const {
00194   if (!seg->good()) return false;
00195 
00196   vector<float> x;
00197   vector<float> y;
00198   vector<float> sigy;
00199 
00200   DTSegmentCand::AssPointCont hits=seg->hits();
00201   for (DTSegmentCand::AssPointCont::const_iterator iter=hits.begin();
00202        iter!=hits.end(); ++iter) {
00203     LocalPoint pos = (*iter).first->localPosition((*iter).second);
00204     x.push_back(pos.z()); 
00205     y.push_back(pos.x());
00206     sigy.push_back(sqrt((*iter).first->localPositionError().xx()));
00207   }
00208 
00209   LocalPoint pos;
00210   LocalVector dir;
00211   AlgebraicSymMatrix covMat(2);
00212 
00213   double chi2 = 0.;
00214   fit(x,y,sigy,pos,dir,covMat,chi2);
00215 
00216   seg->setPosition(pos);
00217   seg->setDirection(dir);
00218 
00219   //cout << "pos " << segPosition<< endl;
00220   //cout << "dir " << segDirection<< endl;
00221 
00222   seg->setCovMatrix(covMat);
00223   // cout << "Mat " << covMat << endl;
00224 
00225   seg->setChi2(chi2);
00226   return true;
00227 }
00228 
00229 void DTSegmentUpdator::fit(DTRecSegment2D* seg) const {
00230   // WARNING: since this method is called both with a 2D and a 2DPhi as argument
00231   // seg->geographicalId() can be a superLayerId or a chamberId 
00232 
00233   vector<float> x;
00234   vector<float> y;
00235   vector<float> sigy;
00236 
00237   vector<DTRecHit1D> hits=seg->specificRecHits();
00238   for (vector<DTRecHit1D>::const_iterator hit=hits.begin();
00239        hit!=hits.end(); ++hit) {
00240 
00241     // I have to get the hits position (the hit is in the layer rf) in SL frame...
00242     GlobalPoint glbPos = ( theGeom->layer( hit->wireId().layerId() ) )->toGlobal(hit->localPosition());
00243     LocalPoint pos = ( theGeom->idToDet(seg->geographicalId()) )->toLocal(glbPos);
00244 
00245     x.push_back(pos.z()); 
00246     y.push_back(pos.x());
00247 
00248     // Get local error in SL frame
00249     //RB: is it right in this way? 
00250     ErrorFrameTransformer tran;
00251     GlobalError glbErr =
00252       tran.transform( hit->localPositionError(),(theGeom->layer( hit->wireId().layerId() ))->surface());
00253     LocalError slErr =
00254       tran.transform( glbErr, (theGeom->idToDet(seg->geographicalId()))->surface());
00255 
00256     sigy.push_back(sqrt(slErr.xx()));
00257   }
00258 
00259   LocalPoint pos;
00260   LocalVector dir;
00261   AlgebraicSymMatrix covMat(2);
00262   double chi2 = 0.;
00263 
00264   fit(x,y,sigy,pos,dir,covMat,chi2);
00265 
00266   seg->setPosition(pos);
00267   seg->setDirection(dir);
00268 
00269   //cout << "pos " << segPosition << endl;
00270   //cout << "dir " << segDirection << endl;
00271 
00272   seg->setCovMatrix(covMat);
00273   // cout << "Mat " << mat << endl;
00274 
00275   seg->setChi2(chi2);
00276 }
00277 
00278 void DTSegmentUpdator::fit(const vector<float>& x,
00279                            const vector<float>& y, 
00280                            const vector<float>& sigy,
00281                            LocalPoint& pos,
00282                            LocalVector& dir,
00283                            AlgebraicSymMatrix& covMatrix,
00284                            double& chi2)  const {
00285 
00286   float slope     = 0.;
00287   float intercept = 0.;
00288   float covss     = 0.;
00289   float covii     = 0.;
00290   float covsi     = 0.;
00291 
00292   // do the fit
00293   theFitter->fit(x,y,x.size(),sigy,slope,intercept,covss,covii,covsi);
00294   // cout << "slope " << slope << endl;
00295   // cout << "intercept " << intercept << endl;
00296 
00297   // intercept is the x() in chamber frame when the segment cross the chamber
00298   // plane (at z()=0), the y() is not measured, so let's put the center of the
00299   // chamber.
00300   pos = LocalPoint(intercept,0.,0.);
00301 
00302   //  slope is dx()/dz(), while dy()/dz() is by definition 0, finally I want the
00303   //  segment to point outward, so opposite to local z
00304   dir = LocalVector(-slope,0.,-1.).unit();
00305 
00306   covMatrix = AlgebraicSymMatrix(2);
00307   covMatrix[0][0] = covss; // this is var(dy/dz)
00308   covMatrix[1][1] = covii; // this is var(y)
00309   covMatrix[1][0] = covsi; // this is cov(dy/dz,y)
00310 
00311   /* Calculate chi2. */
00312   chi2 = 0.;
00313   for(unsigned int i=0; i<x.size() ; ++i) {
00314     double resid= y[i] - (intercept + slope*x[i]);
00315     chi2 += (resid/sigy[i])*(resid/sigy[i]);
00316   }
00317 }
00318 
00319 // The GlobalPoint and the GlobalVector can be either the glb position and the direction
00320 // of the 2D-segment itself or the glb position and direction of the 4D segment
00321 void DTSegmentUpdator::updateHits(DTRecSegment2D* seg, GlobalPoint &gpos,
00322                                   GlobalVector &gdir, const int step) const{
00323 
00324   // it is not necessary to have DTRecHit1D* to modify the obj in the container
00325   // but I have to be carefully, since I cannot make a copy before the iteration!
00326 
00327   vector<DTRecHit1D> toBeUpdatedRecHits = seg->specificRecHits();
00328   vector<DTRecHit1D> updatedRecHits;
00329 
00330   for (vector<DTRecHit1D>::iterator hit= toBeUpdatedRecHits.begin(); 
00331        hit!=toBeUpdatedRecHits.end(); ++hit) {
00332 
00333     const DTLayer* layer = theGeom->layer( hit->wireId().layerId() );
00334 
00335     LocalPoint segPos=layer->toLocal(gpos);
00336     LocalVector segDir=layer->toLocal(gdir);
00337 
00338     // define impact angle needed by the step 2
00339     const float angle = atan(segDir.x()/-segDir.z());
00340 
00341     // define the local position (extr.) of the segment. Needed by the third step 
00342     LocalPoint segPosAtLayer=segPos+segDir*(-segPos.z())/cos(segDir.theta());
00343 
00344     DTRecHit1D newHit1D = (*hit);
00345 
00346     bool ok = true;
00347 
00348     if (step == 2) {
00349       ok = theAlgo->compute(layer,*hit,angle,newHit1D);
00350 
00351     } else if (step == 3) {
00352 
00353       LocalPoint hitPos(hit->localPosition().x(),+segPosAtLayer.y(),0.);
00354 
00355       GlobalPoint glbpos= theGeom->layer( hit->wireId().layerId() )->toGlobal(hitPos);
00356 
00357       newHit1D.setPosition(hitPos);
00358 
00359       ok = theAlgo->compute(layer,*hit,angle,glbpos,newHit1D);
00360 
00361     } else if (step == 4) {
00362 
00363       //const double vminf = seg->vDrift();   //  vdrift correction are recorded in the segment    
00364       double vminf =0.;
00365       if(vdrift_4parfit) vminf = seg->vDrift();   // use vdrift recorded in the segment only if vdrift_4parfit=True
00366 
00367       double cminf = 0.;
00368       if(seg->ist0Valid()) cminf = - seg->t0()*0.00543;
00369 
00370       //cout << "In updateHits: t0 = " << seg->t0() << endl;
00371       //cout << "In updateHits: vminf = " << vminf << endl;
00372       //cout << "In updateHits: cminf = " << cminf << endl;
00373 
00374       const float xwire = layer->specificTopology().wirePosition(hit->wireId().wire());
00375       const float distance = fabs(hit->localPosition().x() - xwire);
00376 
00377       const int ilc = ( hit->lrSide() == DTEnums::Left ) ? 1 : -1;
00378 
00379       const double dy_corr = (vminf*ilc*distance-cminf*ilc ); 
00380 
00381       LocalPoint point(hit->localPosition().x() + dy_corr, +segPosAtLayer.y(), 0.);
00382 
00383       LocalError error(T0_hit_resolution*T0_hit_resolution,0.,0.);
00384 
00385       newHit1D.setPositionAndError(point, error);
00386 
00387       //FIXME: check that the hit is still inside the cell
00388       ok = true;
00389 
00390     } else throw cms::Exception("DTSegmentUpdator")<<" updateHits called with wrong step " << endl;
00391 
00392     if (ok) updatedRecHits.push_back(newHit1D);
00393     else {
00394       LogError("DTSegmentUpdator")<<"DTSegmentUpdator::updateHits failed update" << endl;
00395       throw cms::Exception("DTSegmentUpdator")<<"updateHits failed update"<<endl;
00396     }
00397 
00398   }
00399   seg->update(updatedRecHits);
00400 }
00401 
00402 void DTSegmentUpdator::rejectBadHits(DTChamberRecSegment2D* phiSeg) const {
00403 
00404   vector<float> x;
00405   vector<float> y;
00406   
00407   if(debug) cout << " Inside the segment updator, now loop on hits:   ( x == z_loc , y == x_loc) " << endl;
00408  
00409   vector<DTRecHit1D> hits = phiSeg->specificRecHits();
00410   for (vector<DTRecHit1D>::const_iterator hit=hits.begin();
00411        hit!=hits.end(); ++hit) {
00412 
00413     // I have to get the hits position (the hit is in the layer rf) in SL frame...
00414     GlobalPoint glbPos = ( theGeom->layer( hit->wireId().layerId() ) )->toGlobal(hit->localPosition());
00415     LocalPoint pos = ( theGeom->idToDet(phiSeg->geographicalId()) )->toLocal(glbPos);
00416 
00417     x.push_back(pos.z()); 
00418     y.push_back(pos.x());
00419   }
00420 
00421   if(debug){
00422     cout << " end of segment! " << endl;
00423     cout << " size = Number of Hits: " << x.size() << "  " << y.size() << endl;
00424   }
00425   
00426   // Perform the 2 par fit:
00427   float par[2]={0.,0.}; // q , m
00428 
00429   //variables to perform the fit:
00430   float Sx = 0.;
00431   float Sy = 0.;
00432   float Sx2 = 0.;
00433   float Sy2 = 0.;
00434   float Sxy = 0.;
00435 
00436   size_t N =  x.size();
00437         
00438   for(size_t i = 0; i < N;++i){
00439     Sx += x.at(i);
00440     Sy += y.at(i);
00441     Sx2 += x.at(i)*x.at(i);
00442     Sy2 += y.at(i)*y.at(i);
00443     Sxy += x.at(i)*y.at(i);
00444                 
00445   }
00446         
00447   const float delta = N*Sx2 - Sx*Sx;
00448   par[0] = ( Sx2*Sy - Sx*Sxy )/delta;
00449   par[1] = ( N*Sxy - Sx*Sy )/delta;
00450 
00451   if(debug) cout << "fit 2 parameters done ----> par0: "<< par[0] << "  par1: "<< par[1] << endl;
00452 
00453   // Calc residuals:
00454   float residuals[N];
00455         
00456   for(size_t i = 0; i < N;++i)
00457     residuals[i] = 0;
00458         
00459   for(size_t i = 0; i < N;++i)          
00460     residuals[i] = y.at(i) - par[1]*x.at(i) - par[0];
00461         
00462   if(debug) cout << " Residuals computed! "<<  endl;
00463                 
00464                 
00465   // Perform bad hit rejecting -- update hits
00466   vector<DTRecHit1D> updatedRecHits;
00467         
00468   float mean_residual = 0.; //mean of the absolute values of residuals
00469         
00470   for (size_t i = 0; i < N; ++i)
00471     mean_residual += fabs(residuals[i]);
00472         
00473   mean_residual = mean_residual/(N - 2);        
00474         
00475   if(debug) cout << " mean_residual: "<< mean_residual << endl;
00476         
00477   int i = 0;
00478         
00479   for (vector<DTRecHit1D>::const_iterator hit=hits.begin();
00480        hit!=hits.end(); ++hit) {
00481                 
00482     DTRecHit1D newHit1D = (*hit);
00483 
00484     if(fabs(residuals[i])/mean_residual < 1.5){
00485                                         
00486       updatedRecHits.push_back(newHit1D);
00487       if(debug) cout << " accepted "<< i+1 << "th hit" <<"  Irej: " << fabs(residuals[i])/mean_residual << endl;
00488       ++i;
00489     }
00490     else {
00491       if(debug) cout << " rejected "<< i+1 << "th hit" <<"  Irej: " << fabs(residuals[i])/mean_residual << endl;
00492       ++i;
00493       continue;
00494     }
00495   }
00496         
00497   phiSeg->update(updatedRecHits);       
00498   
00499   //final check!
00500   if(debug){ 
00501   
00502     vector<float> x_upd;
00503     vector<float> y_upd;
00504  
00505     cout << " Check the update action: " << endl;
00506  
00507     vector<DTRecHit1D> hits_upd = phiSeg->specificRecHits();
00508     for (vector<DTRecHit1D>::const_iterator hit=hits_upd.begin();
00509          hit!=hits_upd.end(); ++hit) {
00510 
00511       // I have to get the hits position (the hit is in the layer rf) in SL frame...
00512       GlobalPoint glbPos = ( theGeom->layer( hit->wireId().layerId() ) )->toGlobal(hit->localPosition());
00513       LocalPoint pos = ( theGeom->idToDet(phiSeg->geographicalId()) )->toLocal(glbPos);
00514 
00515       x_upd.push_back(pos.z()); 
00516       y_upd.push_back(pos.x());
00517 
00518       cout << " x_upd: "<< pos.z() << "  y_upd: "<< pos.x() << endl;
00519 
00520 
00521     }
00522   
00523     cout << " end of segment! " << endl;
00524     cout << " size = Number of Hits: " << x_upd.size() << "  " << y_upd.size() << endl;
00525     
00526   }// end debug
00527   
00528   return;
00529 } //end DTSegmentUpdator::rejectBadHits
00530 
00531 void DTSegmentUpdator::calculateT0corr(DTRecSegment4D* seg) const {
00532   if(seg->hasPhi()) calculateT0corr(seg->phiSegment());
00533   if(seg->hasZed()) calculateT0corr(seg->zSegment());
00534 }
00535 
00536 void DTSegmentUpdator::calculateT0corr(DTRecSegment2D* seg) const {
00537   // WARNING: since this method is called both with a 2D and a 2DPhi as argument
00538   // seg->geographicalId() can be a superLayerId or a chamberId 
00539 
00540   vector<double> d_drift;
00541   vector<float> x;
00542   vector<float> y;
00543   vector<int> lc;
00544 
00545   vector<DTRecHit1D> hits=seg->specificRecHits();
00546 
00547   DTWireId wireId;
00548   int nptfit = 0;
00549 
00550   for (vector<DTRecHit1D>::const_iterator hit=hits.begin();
00551        hit!=hits.end(); ++hit) {
00552 
00553     // I have to get the hits position (the hit is in the layer rf) in SL frame...
00554     GlobalPoint glbPos = ( theGeom->layer( hit->wireId().layerId() ) )->toGlobal(hit->localPosition());
00555     LocalPoint pos = ( theGeom->idToDet(seg->geographicalId()) )->toLocal(glbPos);
00556 
00557     const DTLayer* layer = theGeom->layer( hit->wireId().layerId() );
00558     float xwire = layer->specificTopology().wirePosition(hit->wireId().wire());
00559     float distance = fabs(hit->localPosition().x() - xwire);
00560 
00561     int ilc = ( hit->lrSide() == DTEnums::Left ) ? 1 : -1;
00562 
00563     nptfit++;
00564     x.push_back(pos.z()); 
00565     y.push_back(pos.x());
00566     lc.push_back(ilc);
00567     d_drift.push_back(distance);
00568 
00569     // cout << " d_drift "<<distance  <<" npt= " <<npt<<endl;
00570   }
00571 
00572   double chi2fit = 0.;
00573   float cminf    = 0.;
00574   double vminf   = 0.;
00575 
00576   if ( nptfit > 2 ) {
00577     //NB chi2fit is normalized
00578     Fit4Var(x,y,lc,d_drift,nptfit,cminf,vminf,chi2fit);
00579 
00580     double t0cor = -999.;
00581     if(cminf > -998.) t0cor = - cminf/0.00543 ; // in ns
00582 
00583     //cout << "In calculateT0corr: t0 = " << t0cor << endl;
00584     //cout << "In calculateT0corr: vminf = " << vminf << endl;
00585     //cout << "In calculateT0corr: cminf = " << cminf << endl;
00586     //cout << "In calculateT0corr: chi2 = " << chi2fit << endl;
00587 
00588     seg->setT0(t0cor);          // time  and
00589     seg->setVdrift(vminf);   //  vdrift correction are recorded in the segment    
00590   }
00591 }
00592 
00593 
00594 void DTSegmentUpdator::Fit4Var(const vector<float>& xfit,
00595                                const vector<float>& yfit,
00596                                const vector<int>& lfit,
00597                                const vector<double>& tfit,
00598                                const int nptfit,
00599                                float& cminf,
00600                                double& vminf,
00601                                double& chi2fit) const { 
00602 
00603   const double sigma = 0.0295;// errors can be inserted .just load them/that is the usual TB resolution value for DT chambers 
00604   double aminf = 0.;
00605   double bminf = 0.;
00606   int nppar = 0;
00607   double sx = 0.;
00608   double  sx2 = 0.;
00609   double sy = 0.;
00610   double sxy = 0.;
00611   double sl = 0.;
00612   double sl2 = 0.;
00613   double sly = 0.;
00614   double slx = 0.;
00615   double st = 0.;
00616   double st2 = 0.;
00617   double slt = 0.;
00618   double sltx = 0.;
00619   double slty = 0.;
00620   double chi2fitN2 = -1. ;
00621   double chi2fit3 = -1.;
00622   double chi2fitN3 = -1. ;
00623   double chi2fitN4 = -1.;
00624   float bminf3 = bminf;
00625   float aminf3 = aminf;
00626   float cminf3 = cminf;
00627   int nppar2 = 0;
00628   int nppar3 = 0;
00629   int nppar4 = 0;
00630 
00631   cminf = -999.;
00632   vminf = 0.;
00633 
00634   for (int j=0; j<nptfit; j++){
00635     sx  = sx + xfit[j];       
00636     sy  = sy + yfit[j];
00637     sx2 = sx2 + xfit[j]*xfit[j];
00638     sxy = sxy + xfit[j]*yfit[j];
00639     sl  = sl + lfit[j];       
00640     sl2 = sl2 + lfit[j]*lfit[j];
00641     sly = sly + lfit[j]*yfit[j];
00642     slx = slx + lfit[j]*xfit[j];
00643     st = st + tfit[j];
00644     st2 = st2 + tfit[j] * tfit[j];
00645     slt = slt + lfit[j] * tfit[j];
00646     sltx = sltx + lfit[j] * tfit[j]*xfit[j];
00647     slty = slty + lfit[j] * tfit[j]*yfit[j];
00648 
00649   } //end loop
00650 
00651   const double delta = nptfit*sx2 - sx*sx;
00652 
00653   double a = 0.;
00654   double b = 0.;               
00655 
00656   if (delta!=0){   //
00657     a = (sx2*sy - sx*sxy)/delta;
00658     b = (nptfit*sxy - sx*sy)/delta;
00659 
00660     //  cout << " NPAR=2 : slope = "<<b<< "    intercept = "<<a <<endl;
00661     for (int j=0; j<nptfit; j++){
00662       const double ypred = a + b*xfit[j];
00663       const double dy = (yfit[j] - ypred)/sigma;
00664       chi2fit = chi2fit + dy*dy;
00665     } //end loop chi2
00666   }
00667 
00668   bminf = b;
00669   aminf = a;
00670 
00671   nppar = 2; 
00672   nppar2 = nppar; 
00673 
00674   chi2fitN2 = chi2fit/(nptfit-2);
00675 
00676   // cout << "dt0 = 0chi2fit = " << chi2fit << "  slope = "<<b<<endl;
00677 
00678   if (nptfit >= 3) {
00679 
00680     const double d1 = sy;
00681     const double d2 = sxy;
00682     const double d3 = sly;
00683     const double c1 = sl;
00684     const double c2 = slx;
00685     const double c3 = sl2;
00686     const double b1 = sx;
00687     const double b2 = sx2;
00688     const double b3 = slx;
00689     const double a1 = nptfit;
00690     const double a2 = sx;
00691     const double a3 = sl;
00692 
00693     //these parameters are not used in the 4-variables fit
00694     const double b4 = b2*a1-b1*a2;
00695     const double c4 = c2*a1-c1*a2;
00696     const double d4 = d2*a1-d1*a2;
00697     const double b5 = a1*b3-a3*b1;
00698     const double c5 = a1*c3-a3*c1;
00699     const double d5 = a1*d3-d1*a3;
00700     const double a6 = slt;
00701     const double b6 = sltx;
00702     const double c6 = st;
00703     const double v6 = st2;      
00704     const double d6 = slty;
00705 
00706     if (((c5*b4-c4*b5)*b4*a1)!=0) {
00707       nppar = 3;
00708       chi2fit = 0.;
00709       cminf = (d5*b4-d4*b5)/(c5*b4-c4*b5);
00710       bminf = d4/b4 -cminf *c4/b4;
00711       aminf = (d1/a1 -cminf*c1/a1 -bminf*b1/a1);
00712 
00713       for (int j=0; j<nptfit; j++){
00714         const double ypred = aminf + bminf*xfit[j];
00715         const double dy = (yfit[j]-cminf*lfit[j] - ypred)/sigma;
00716         chi2fit = chi2fit + dy*dy;
00717 
00718       } //end loop chi2
00719       chi2fit3 = chi2fit;
00720       if (nptfit>3)
00721         chi2fitN3 = chi2fit /(nptfit-3);
00722 
00723     }
00724     else {
00725       cminf = -999.;
00726       bminf = b;
00727       aminf = a;
00728       chi2fit3 = chi2fit;
00729       chi2fitN3 = chi2fit /(nptfit-2);
00730     }
00731 
00732     bminf3 = bminf;
00733     aminf3 = aminf;
00734     cminf3 = cminf;
00735     nppar3 = nppar;
00736 
00737     if (debug) {
00738       cout << "dt0= 0 : slope 2 = " << b << " pos in  = " << a << " chi2fitN2 = " << chi2fitN2
00739            << " nppar = " << nppar2 << " nptfit = " << nptfit << endl;
00740       cout << "dt0 = 0 : slope 3 = " << bminf << " pos out = " << aminf << " chi2fitN3 = "
00741            << chi2fitN3 << " nppar = " << nppar3 << " T0_ev ns = " << cminf/0.00543 << endl;
00742     } 
00743 
00744     //***********************************
00745     //     cout << " vdrift_4parfit "<< vdrift_4parfit<<endl;
00746     if( nptfit>=5) { 
00747       const double det = (a1*a1*(b2*v6 - b6*b6) - a1*(a2*a2*v6 - 2*a2*a6*b6 + a6*a6*b2 + b2*c6*c6 + b3*(b3*v6 - 2*b6*c6))
00748                           + a2*a2*c6*c6 + 2*a2*(a3*(b3*v6 - b6*c6) - a6*b3*c6) + a3*a3*(b6*b6 - b2*v6)
00749                           + a6*(2*a3*(b2*c6 - b3*b6) + a6*b3*b3)); 
00750 
00751       // the dv/vdrift correction may be computed  under vdrift_4parfit request;
00752       if (det != 0) { 
00753         nppar = 4;
00754         chi2fit = 0.;
00755         // computation of   a, b, c e v
00756         aminf = (a1*(a2*(b6*d6 - v6*d2) + a6*(b6*d2 - b2*d6) + d1*(b2*v6 - b6*b6)) - a2*(b3*(c6*d6 - v6*d3)
00757                  + c6*(b6*d3 - c6*d2)) + a3*(b2*(c6*d6 - v6*d3) + b3*(v6*d2 - b6*d6) + b6*(b6*d3 - c6*d2))
00758                  + a6*(b2*c6*d3 + b3*(b3*d6 - b6*d3 - c6*d2)) - d1*(b2*c6*c6 + b3*(b3*v6 - 2*b6*c6)))/det;
00759         bminf = - (a1*a1*(b6*d6 - v6*d2) - a1*(a2*(a6*d6 - v6*d1) - a6*a6*d2 + a6*b6*d1 + b3*(c6*d6 - v6*d3)
00760                  + c6*(b6*d3 - c6*d2)) + a2*(a3*(c6*d6 - v6*d3) + c6*(a6*d3 - c6*d1)) + a3*a3*(v6*d2 - b6*d6)
00761                  + a3*(a6*(b3*d6 + b6*d3 - 2*c6*d2) - d1*(b3*v6 - b6*c6)) - a6*b3*(a6*d3 - c6*d1))/det;
00762         cminf = -(a1*(b2*(c6*d6 - v6*d3) + b3*(v6*d2 - b6*d6) + b6*(b6*d3 - c6*d2)) + a2*a2*(v6*d3 - c6*d6)
00763                  + a2*(a3*(b6*d6 - v6*d2) + a6*(b3*d6 - 2*b6*d3 + c6*d2) - d1*(b3*v6 - b6*c6))
00764                  + a3*(d1*(b2*v6 - b6*b6) - a6*(b2*d6 - b6*d2)) + a6*(a6*(b2*d3 - b3*d2) - d1*(b2*c6 - b3*b6)))/det;
00765         vminf = - (a1*a1*(b2*d6 - b6*d2) - a1*(a2*a2*d6 - a2*(a6*d2 + b6*d1) + a6*b2*d1 + b2*c6*d3
00766                  + b3*(b3*d6 - b6*d3 - c6*d2)) + a2*a2*c6*d3 + a2*(a3*(2*b3*d6 - b6*d3 - c6*d2) - b3*(a6*d3 + c6*d1))
00767                  + a3*a3*(b6*d2 - b2*d6) + a3*(a6*(b2*d3 - b3*d2) + d1*(b2*c6 - b3*b6)) + a6*b3*b3*d1)/det;
00768 
00769         //  chi 2
00770         for (int j=0; j<nptfit; j++) {
00771           const double ypred = aminf + bminf*xfit[j];
00772           const double dy = (yfit[j]+vminf*lfit[j]*tfit[j]-cminf*lfit[j] -ypred)/sigma; 
00773           chi2fit = chi2fit + dy*dy;
00774 
00775         } //end loop chi2
00776         if (nptfit<=nppar){ 
00777           chi2fitN4=-1;
00778           //            cout << "nptfit " << nptfit << " nppar " << nppar << endl;
00779         }
00780         else{
00781           chi2fitN4= chi2fit / (nptfit-nppar); 
00782         }
00783       }
00784       else {
00785         vminf = 0.;
00786 
00787         if (nptfit <= nppar) chi2fitN4=-1;
00788         else chi2fitN4  = chi2fit / (nptfit-nppar); 
00789       }
00790 
00791       if (fabs(vminf) >= 0.29) {
00792         // for safety and for code construction..dont accept correction on dv/vdrift greater then 0.09
00793         vminf = 0.;
00794         cminf = cminf3;
00795         aminf = aminf3;
00796         bminf = bminf3;
00797         nppar = 3;
00798         chi2fit = chi2fit3;
00799       }
00800 
00801     }  //end if vdrift
00802 
00803      if(!vdrift_4parfit){         //if not required explicitly leave the t0 and track step as at step 3
00804                                   // just update vdrift value vmin for storing in the segments for monitoring
00805        cminf = cminf3;
00806        aminf = aminf3;
00807        bminf = bminf3;
00808        nppar = 3;
00809        chi2fit = chi2fit3;
00810      }
00811 
00812     nppar4 = nppar;
00813 
00814   }  //end nptfit >=3
00815 
00816   if (debug) {
00817     cout << "   dt0= 0 : slope 4  = " << bminf << " pos out = " << aminf <<" chi2fitN4 = " << chi2fitN4
00818          << "  nppar= " << nppar4 << " T0_ev ns= " << cminf/0.00543 <<" delta v = " << vminf <<endl;
00819     cout << nptfit << " nptfit " << " end  chi2fit = " << chi2fit/ (nptfit-nppar ) << " T0_ev ns= " << cminf/0.00543 << " delta v = " << vminf <<endl;
00820   }
00821 
00822   if ( fabs(vminf) >= 0.09 && debug ) {  //checks only vdrift less then 10 % accepted
00823     cout << "vminf gt 0.09 det=  " << endl;
00824     cout << "dt0= 0 : slope 4 = "<< bminf << " pos out = " << aminf << " chi2fitN4 = " << chi2fitN4
00825          << " T0_ev ns = " << cminf/0.00543 << " delta v = "<< vminf << endl;
00826     cout << "dt0 = 0 : slope 2 = "<< b << " pos in = " << a <<" chi2fitN2 = " << chi2fitN2
00827          << " nppar = " << nppar-1 << " nptfit = " << nptfit <<endl;
00828     cout << "dt0 = 0 : slope 3 = " << bminf << " pos out = " << aminf << " chi2fitN3 = "
00829          << chi2fitN3 << " T0_ev ns = " << cminf/0.00543 << endl;
00830     cout << nptfit   <<" nptfit "<< "   end  chi2fit = " << chi2fit << "T0_ev ns= " << cminf/0.00543 << "delta v = "<< vminf <<endl;        
00831   }
00832 
00833   if (nptfit != nppar) chi2fit = chi2fit / (nptfit-nppar);
00834 }