GFKalman.cxx
Go to the documentation of this file.
1 /* Copyright 2008-2010, Technische Universitaet Muenchen,
2  Authors: Christian Hoeppner & Sebastian Neubert
3 
4  This file is part of GENFIT.
5 
6  GENFIT is free software: you can redistribute it and/or modify
7  it under the terms of the GNU Lesser General Public License as published
8  by the Free Software Foundation, either version 3 of the License, or
9  (at your option) any later version.
10 
11  GENFIT is distributed in the hope that it will be useful,
12  but WITHOUT ANY WARRANTY; without even the implied warranty of
13  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14  GNU Lesser General Public License for more details.
15 
16  You should have received a copy of the GNU Lesser General Public License
17  along with GENFIT. If not, see <http://www.gnu.org/licenses/>.
18 */
20 
21 #include <iostream>
22 
23 #include "TMath.h"
24 #include "TRandom.h"
25 #include "TDatabasePDG.h"
26 
27 #include "larreco/Genfit/GFTrack.h"
31 
32 #include "cetlib_except/exception.h"
34 #include "art_root_io/TFileService.h"
35 
36 #define COVEXC "cov_is_zero"
37 
38 
39 genf::GFKalman::GFKalman():fInitialDirection(1),fNumIt(3),fBlowUpFactor(50.),fMomLow(-100.0),fMomHigh(100.0),fMaxUpdate(1.0),fErrScaleSTh(1.0),fErrScaleMTh(1.0), fGENfPRINT(false)
40 {
42 }
43 
45 
48  if ((direction != 1) && (direction != -1))
49  throw GFException(std::string(__func__) + ": wrong direction", __LINE__, __FILE__).setFatal();
50  // trk->clearGFBookkeeping();
51  trk->clearRepAtHit();
52  /*
53  int nreps=trk->getNumReps();
54  for(int i=0; i<nreps; ++i){
55  trk->getBK(i)->setNhits(trk->getNumHits());
56  trk->getBK(i)->bookMatrices("fPredSt");
57  trk->getBK(i)->bookMatrices("fPredCov");
58  trk->getBK(i)->bookMatrices("bPredSt");
59  trk->getBK(i)->bookMatrices("bPredCov");
60  trk->getBK(i)->bookGFDetPlanes("f");
61  trk->getBK(i)->bookGFDetPlanes("b");
62  trk->getBK(i)->bookNumbers("fChi2");
63  trk->getBK(i)->bookNumbers("bChi2");
64  }
65  */
66 
67  /*why is there a factor of two here (in the for statement)?
68  Because we consider one full iteration to be one back and
69  one forth fitting pass */
70  for(int ipass=0; ipass<2*fNumIt; ipass++){
71  // std::cout << "GFKalman:: numreps, numhits, ipass, direction"<< trk->getNumReps()<<", "<< trk->getNumHits()<<", "<< ipass<<", " <<direction << std::endl;
72  // if(floor(ipass/2)==fNumIt && direction==1) blowUpCovsDiag(trk);
73  if(ipass>0) blowUpCovs(trk); // Back to >0, not >=0 EC, 9-11-Feb-2013
74  if(direction==1){
75  trk->setNextHitToFit(0);
76  }
77  else {
78  trk->setNextHitToFit(trk->getNumHits()-1);
79  }
80 
81 
82  try
83  {
84  fittingPass(trk,direction);
85  }
86  catch(GFException &)
87  {
88  throw cet::exception("GFKalman.cxx: ") << " Line " << __LINE__ << ", " << __FILE__ << " ...Rethrow. \n";
89  }
90 
91  //save first and last plane,state&cov after the fitting pass
92  if(direction==1){//forward at last hit
93  int nreps=trk->getNumReps();
94  for(int i=0; i<nreps; ++i){
95  trk->getTrackRep(i)->
96  setLastPlane( trk->getTrackRep(i)->getReferencePlane() );
97  trk->getTrackRep(i)->
98  setLastState( trk->getTrackRep(i)->getState() );
99  trk->getTrackRep(i)->
100  setLastCov( trk->getTrackRep(i)->getCov() );
101  }
102  }
103  else{//backward at first hit
104  int nreps=trk->getNumReps();
105  for(int i=0; i<nreps; ++i){
106  trk->getTrackRep(i)->
107  setFirstPlane( trk->getTrackRep(i)->getReferencePlane() );
108  trk->getTrackRep(i)->
109  setFirstState( trk->getTrackRep(i)->getState() );
110  trk->getTrackRep(i)->
111  setFirstCov( trk->getTrackRep(i)->getCov() );
112  }
113  }
114 
115  //switch direction of fitting and also inside all the reps
116  if(direction==1) direction=-1;
117  else direction=1;
118  switchDirection(trk);
119 
120  }
121  return;
122 }
123 
124 void
126  int nreps=trk->getNumReps();
127  for(int i=0; i<nreps; ++i){
128  trk->getTrackRep(i)->switchDirection();
129  }
130 }
131 
133  int nreps=trk->getNumReps();
134  for(int irep=0; irep<nreps; ++irep){
135  GFAbsTrackRep* arep=trk->getTrackRep(irep);
136  //dont do it for already compromsied reps, since they wont be fitted anyway
137  if(arep->getStatusFlag()==0) {
138  TMatrixT<Double_t> cov = arep->getCov();
139  for(int i=0;i<cov.GetNrows();++i){
140  for(int j=0;j<cov.GetNcols();++j){
141  if(i!=j){//off diagonal
142  cov[i][j]=0.;
143  }
144  else{//diagonal
145  cov[i][j] = cov[i][j] * fBlowUpFactor;
146  // cov[0][0] = 0.1;
147  }
148  }
149  }
150  arep->setCov(cov);
151  }
152  }
153 }
155  int nreps=trk->getNumReps();
156  for(int irep=0; irep<nreps; ++irep){
157  GFAbsTrackRep* arep=trk->getTrackRep(irep);
158  //dont do it for already compromsied reps, since they wont be fitted anyway
159  if(arep->getStatusFlag()==0) {
160  TMatrixT<Double_t> cov = arep->getCov();
161  for(int i=0;i<cov.GetNrows();++i){
162  for(int j=0;j<cov.GetNcols();++j){
163  //if(i!=j){//off diagonal
164  // cov[i][j]=0.;
165  //}
166  //else{//diagonal
167  cov[i][j] = cov[i][j] * fBlowUpFactor;
168  //}
169  }
170  }
171  arep->setCov(cov);
172  }
173  }
174 }
175 
176 void
178  //loop over hits
179 
180  unsigned int nhits=trk->getNumHits();
181  unsigned int starthit=trk->getNextHitToFit();
182 
183  int nreps=trk->getNumReps();
184  int ihit=(int)starthit;
185 
186  //clear chi2 sum and ndf sum in track reps
187  for(int i=0;i<nreps;++i){
188  GFAbsTrackRep* arep=trk->getTrackRep(i);
189  arep->setChiSqu(0.);
190  arep->setNDF(0);
191  }
192 
193  //clear failedHits and outliers
194  trk->getBK()->clearFailedHits();
195 
196  while((ihit<(int)nhits && direction==1) || (ihit>-1 && direction==-1)){
197  // GFAbsRecoHit* ahit=trk->getHit(ihit);
198  // loop over reps
199  for(int irep=0; irep<nreps; ++irep){
200  GFAbsTrackRep* arep=trk->getTrackRep(irep);
201  if(arep->getStatusFlag()==0) {
202  try {
203  processHit(trk,ihit,irep,direction);
204  }
205  catch(GFException& e) {
206  trk->addFailedHit(irep,ihit);
207  std::cerr << e.what();
208  e.info();
209  if(e.isFatal()) {
210  arep->setStatusFlag(1);
211  continue; // go to next rep immediately
212  }
213  }
214  }
215  }// end loop over reps
216  ihit+=direction;
217  }// end loop over hits
218  trk->setNextHitToFit(ihit-2*direction);
219  //trk->printGFBookkeeping();
220 }
221 
222 double genf::GFKalman::chi2Increment(const TMatrixT<Double_t>& r,const TMatrixT<Double_t>& H,
223  const TMatrixT<Double_t>& cov,const TMatrixT<Double_t>& V){
224 
225  // residuals covariances:R=(V - HCH^T)
226  TMatrixT<Double_t> R(V);
227  TMatrixT<Double_t> covsum1(cov,TMatrixT<Double_t>::kMultTranspose,H);
228  TMatrixT<Double_t> covsum(H,TMatrixT<Double_t>::kMult,covsum1);
229 
230  R-=covsum;
231 
232  // chisq= r^TR^(-1)r
233  double det=0.;
234  TMatrixT<Double_t> Rsave(R);
235  R.SetTol(1.0e-30); // to avoid inversion problem, EC, 8-Aug-2011. Was23, 9-Jan-2012.
236 
237  try
238  {
239  R.Invert(&det);
240  }
241  catch (cet::exception &)
242  {
243  GFException e("Kalman Chi2Increment: R is not even invertible. But keep plowing on ... ",__LINE__,__FILE__);
244  //e.setFatal();
245  //throw e;
246  }
247  if(TMath::IsNaN(det)) {
248  throw GFException("Kalman Chi2Increment: det of covsum is nan",__LINE__,__FILE__).setFatal();
249  }
250  TMatrixT<Double_t> residTranspose(r);
251  residTranspose.T();
252  TMatrixT<Double_t> chisq=residTranspose*(R*r);
253  assert(chisq.GetNoElements()==1);
254 
255  if(TMath::IsNaN(chisq[0][0])){
256  GFException exc("chi2 is nan",__LINE__,__FILE__);
257  exc.setFatal();
258  std::vector<double> numbers;
259  numbers.push_back(det);
260  exc.setNumbers("det",numbers);
261  std::vector< TMatrixT<Double_t> > matrices;
262  matrices.push_back(r);
263  matrices.push_back(V);
264  matrices.push_back(Rsave);
265  matrices.push_back(R);
266  matrices.push_back(cov);
267  exc.setMatrices("r, V, Rsave, R, cov",matrices);
268  throw exc;
269  }
270 
271  return chisq[0][0];
272 }
273 
274 
275 double
277 {
278  // get prototypes for matrices
279  int repDim=rep->getDim();
280  TMatrixT<Double_t> state(repDim,1);
281  TMatrixT<Double_t> cov(repDim,repDim);;
282  GFDetPlane pl=hit->getDetPlane(rep);
283  rep->extrapolate(pl,state,cov);
284 
285 
286  TMatrixT<Double_t> H = hit->getHMatrix(rep);
287 
288  // get hit covariances
289  TMatrixT<Double_t> V=hit->getHitCov(pl);
290  V[0][0] *= fErrScaleMTh;
291  TMatrixT<Double_t> r=hit->residualVector(rep,state,pl);
292  assert(r.GetNrows()>0);
293 
294  r[0][0] = fabs(r[0][0]);
295  //this is where chi2 is calculated
296  double chi2 = chi2Increment(r,H,cov,V);
297 
298  return chi2/r.GetNrows();
299 }
300 
301 
302 void
303 genf::GFKalman::processHit(GFTrack* tr, int ihit, int irep,int direction){
304  GFAbsRecoHit* hit = tr->getHit(ihit);
305  GFAbsTrackRep* rep = tr->getTrackRep(irep);
306 
307  // get prototypes for matrices
308  int repDim = rep->getDim();
309  TMatrixT<Double_t> state(repDim,1);
310  TMatrixT<Double_t> cov(repDim,repDim);
311  static TMatrixT<Double_t> covFilt(cov);
312  const double pi2(10.0);
313  GFDetPlane pl, plPrev;
314  unsigned int nhits=tr->getNumHits();
315  int phit=ihit;
316  static TMatrixT<Double_t> oldState(5,1);
317  static std::vector<TVector3> pointsPrev;
318 
319  const double eps(1.0e-6);
320  if (direction>0 && ihit>0)
321  {
322  phit = ihit - 1;
323  }
324  if (direction<0 && ihit<((int)nhits-1))
325  {
326  phit = ihit + 1;
327  }
328  GFAbsRecoHit* hitPrev = tr->getHit(phit);
329 
330  /* do an extrapolation, if the trackrep irep is not given
331  * at this ihit position. This will usually be the case, but
332  * not if the fit turns around
333  */
334  // std::cout << "GFKalman::ProcessHit(): ihit is " << ihit << std::endl;
335  // std::cout << "GFKalman::ProcessHit(): direction is " << direction << std::endl;
336  //rep->Print();
337 
338 
339  Double_t thetaPlanes(0.0);
340  if(ihit!=tr->getRepAtHit(irep)){
341  //std::cout << "not same" << std::endl;
342  // get the (virtual) detector plane. This call itself calls
343  // extrapolateToPoint(), which calls Extrap() with just 2 args and
344  // default 3rd arg, which propagates track through
345  // material to find the next plane. But it in fact seems
346  // to just return this pl and plPrev, as desired. Something in Extrap()
347  // kicks it out... even though I can see it walking through material ...
348  // A mystery.
349 
350  pl=hit->getDetPlane(rep);
351  plPrev=hitPrev->getDetPlane(rep);
352 
353  // std::cout << "GFKalman::ProcessHit(): hit is ... " << ihit << std::endl;
354  //hit->Print();
355  //std::cout << "GFKalman::ProcessHit(): plane is ... " << std::endl;
356  //pl.Print();
357 
358  //do the extrapolation. This calls Extrap(), which wraps up RKutta and
359  // GFMaterials calls. state is intialized, pl is unchanged. Latter behavior
360  // I will alter below.
361  try{
362  rep->extrapolate(pl,state,cov);
363  /*
364  if ( std::isnan(cov[0][0]) )
365  {
366  cov = covFilt;
367  }
368  else
369  {
370  covFilt = cov;
371  }
372  */
373  }
374  catch (cet::exception &)
375  {
376  throw cet::exception("GFKalman.cxx: ") << " Line " << __LINE__ << ", " << __FILE__ << " ...Rethrow. \n";
377  }
378 
379  }
380  else{
381  pl = rep->getReferencePlane();
382  plPrev = hitPrev->getDetPlane(rep);
383  state = rep->getState();
384  cov = rep->getCov();
385  }
386 
387  // Update plane. This is a code change from Genfit out of box.
388  // state has accumulated the material/magnetic field changes since last
389  // step. Here, we make the *plane* at this step know about that bend.
390  // We will not, in the end, save this plane, cuz Genfit knows to properly
391  // calculate it later.
392  // Actually, I do *not* change the plane. EC, 11-May-2012.
393  TVector3 u(pl.getU());
394  TVector3 v(pl.getV());
395  TVector3 wold(u.Cross(v));
396  //Double_t sign(1.0);
397  //if ((direction==-1) && ihit==((int)nhits-1)) sign = -1.0;
398 
399  TVector3 pTilde = direction * (wold + state[1][0] * u + state[2][0] * v);
400  TVector3 w(pTilde.Unit());
401  // Find angle/vector through which we rotate. Use it subsequently.
402  TVector3 rot(wold.Cross(w));
403  Double_t ang(TMath::ACos(w*wold));
404  ang = TMath::Min(ang,0.95*TMath::Pi()/2.0);
405  /*
406  {
407  u.Rotate(ang,rot);
408  v.Rotate(ang,rot);
409 
410  pl.setNormal(w);
411  pl.setU(u.Unit());
412  pl.setV(v.Unit());
413  }
414  */
415  if(cov[0][0]<1.E-50 || TMath::IsNaN(cov[0][0])){
416  //std::cout<<"GFKalman::processHit() 0. Calling Exception."<<std::endl;
417  GFException exc(COVEXC,__LINE__,__FILE__);
418  if (fGENfPRINT) cov.Print();
419  if (fGENfPRINT) std::cout<<"GFKalman::processHit() 1. No longer throw exception. Force cov[0][0] to 0.01."<<std::endl;
420  // std::cout<<"GFKalman::processHit() 1. About to throw GFException."<<std::endl;
421  cov = covFilt;
422  // throw exc;
423  // std::cout<<"GFKalman::processHit() 2. Back from GFException."<<std::endl;
424 
425  }
426 
427  /*
428  if(direction==1){
429  tr->getBK(irep)->setMatrix("fPredSt",ihit,state);
430  tr->getBK(irep)->setMatrix("fPredCov",ihit,cov);
431  tr->getBK(irep)->setGFDetPlane("f",ihit,pl);
432  }
433  else{
434  tr->getBK(irep)->setMatrix("bPredSt",ihit,state);
435  tr->getBK(irep)->setMatrix("bPredCov",ihit,cov);
436  tr->getBK(irep)->setGFDetPlane("b",ihit,pl);
437  }
438  */
439 
440  // TMatrixT<Double_t> origcov=rep->getCov();
441 
442  int pdg = tr->getPDG();
443  TParticlePDG * part = TDatabasePDG::Instance()->GetParticle(pdg);
444  Double_t mass = part->Mass();
445 
446  Double_t dist = (pl.getO()-plPrev.getO()).Mag();
447  Double_t mom = fabs(1.0/state[0][0]);
448  Double_t beta = mom/sqrt(mass*mass+mom*mom);
449  const Double_t lowerLim(0.01);
450  if (std::isnan(dist) || dist<=0.0) dist=lowerLim; // don't allow 0s here.
451  if (std::isnan(beta) || beta<0.01) beta=0.01;
452  TMatrixT<Double_t> H=hit->getHMatrix(rep,beta,dist);
453 
454 
455  // Can force huge error here on p and du,v/dw, and thus insensitivity to p,
456  // by setting V[0][0]->inf.
457  TMatrixT<Double_t> V=hit->getHitCov(pl,plPrev,state,mass);
458  V[0][0] *= fErrScaleMTh;
459  tr->setHitMeasuredCov(V);
460 
461  // calculate kalman gain ------------------------------
462  TMatrixT<Double_t> Gain(calcGain(cov,V,H));
463 
464  //std::cout << "GFKalman:: processHits(), state is " << std::endl;
465  //rep->getState().Print();
466 
467  TMatrixT<Double_t> res=hit->residualVector(rep,state,pl,plPrev,mass);
468  // calculate update -----------------------------------
469 
470 
471  // TVector3 pointer((pl.getO()-plPrev.getO()).Unit());
472 
473  TMatrixT<Double_t> rawcoord = hit->getRawHitCoord();
474  TVector3 point(rawcoord[0][0],rawcoord[1][0],rawcoord[2][0]);
475  TMatrixT<Double_t> prevrawcoord = hitPrev->getRawHitCoord();
476  TVector3 pointPrev(prevrawcoord[0][0],prevrawcoord[1][0],prevrawcoord[2][0]);
477  pointsPrev.push_back(pointPrev);
478  TMatrixT<Double_t> Hnew(H);
479  if ((ihit==((int)nhits-1)&&direction==-1) || (ihit==0&&direction==1))
480  pointsPrev.clear();
481 
482  TVector3 pointer((point-pointPrev).Unit());
483  static TVector3 pointerPrev(pointer);
484  if (ihit==0&&direction==1 )
485  {pointer[0] = 0.0;pointer[1] = 0.0;pointer[2] = 1.0;}
486  if (ihit==((int)nhits-1)&&direction==-1)
487  {pointer[0] = 0.0;pointer[1] = 0.0;pointer[2] = -1.0;}
488  double thetaMeas = TMath::Min(fabs(pointer.Angle(pointerPrev)),0.95*TMath::Pi()/2.0);
489  // Below line introduced because it's not true we predict the angle to be
490  // precisely ang. If we'd taken this quantity from our transport result
491  // (with measurement errors in it) we'd expect a smear like this on ang.
492  // EC, 22-Feb-2012.
493  // Error is taken on th^2
494  ang = (Hnew*state)[0][0]; // H->Hnew
495 
496  // fErrScale is extra-fun bonus factor!
497  //thetaPlanes = ang*ang; // + fErrScaleSTh*gRandom->Gaus(0.0,V[0][0]);
498  thetaPlanes = gRandom->Gaus(0.0,ang*ang);
499  thetaPlanes = TMath::Min(sqrt(fabs(thetaPlanes)),0.95*TMath::Pi()/2.0);
500 
501  Double_t dtheta = thetaMeas - thetaPlanes; // was fabs(res[0][0]). EC, 26-Jan-2012
502  if (((ihit==((int)nhits-1)||ihit==((int)nhits-2))&&direction==-1) || ((ihit==0||ihit==1)&&direction==1))
503  {
504  dtheta = 0.0; // at the bare minimum (Jan-2013). Now (Feb-2013), ....
505  // Do not let these 4 points*direction influence the state or chi2.
506  rep->setData(state,pl,&cov); // This is where fState,
507  // fRefPlane and fCov are updated.
508  pointerPrev = pointer;
509  return;
510  }
511 
512  oldState = state;
513 
514  res[0][0] = dtheta/Hnew[0][0];
515  // The particle was propagated through material until its poca
516  // to this current hit was found. At that point its plane is defined,
517  // in which du/dw=dv/dw=0 by construction. But, there's a deviation
518  // in those two quantities measured in the *previous* plane which we want.
519  TVector3 uPrev(plPrev.getU());
520  TVector3 vPrev(plPrev.getV());
521  TVector3 wPrev(u.Cross(v));
522  // We want the newest momentum vector's d{u,v}/dw defined in the prev plane.
523  // That is the predicted du,v/dw. We will subtract from the actual.
524  // But the u's and v's need to come from the State not the Planes, cuz I
525  // haven't updated pl,plPrev per latest State. plFilt, the updated plPrev,
526  // is what we want. w is from updated new plane pl, per latest State.
527  // e.g. plFilt.
528  static GFDetPlane plFilt;//(pl);
529  if (plFilt.getO().Mag()>eps)
530  {
531  uPrev = plFilt.getU();
532  vPrev = plFilt.getV();
533  wPrev = plFilt.getNormal();
534  }
535  //res[1][0] = (pointer*uPrev)/(pointer*wPrev) - (w*uPrev)/(w*wPrev);
536  //res[2][0] = (pointer*vPrev)/(pointer*wPrev) - (w*vPrev)/(w*wPrev);
537  // res[1][0] = 0.0;
538  // res[2][0] = 0.0;
539 
540  TMatrixT<Double_t> update=Gain*res;
541 
542  // Don't let crazy ass spacepoints which pull the fit all over the place
543  // be allowed to update the state. Use __ xRMS for fMaxUpdate.
544  // Use a larger limit (~0.1) when starting with too-high seed momentum.
545  // Smaller, when starting with an underestimate.
546  //
547  // Don't allow fits above 20 GeV, or so. Else, 1/p will likely
548  // migrate across zero. Similarly, don't allow tiny momentum.
549  //
550  tr->setHitUpdate(update);
551  if (fabs(update[0][0])>fMaxUpdate )
552  { // zero the gain so as to not update cov, state
553  // zero the updates
554  update[0][0] = 0.0;/* */
555  // update[0][0] = fMaxUpdate*update[0][0]/fabs(update[0][0]);
556  // either of below leads to craziness, somehow.
557  // update.Zero();
558  // Gain.Zero();
559  }
560  state+=update; // prediction overwritten!
561 
562  // Debugging purposes
563  TMatrixT<Double_t> GH(Gain*Hnew);
564  if (GH[0][0] > 1.)
565  {
566  // std::cout << "GFKalman:: Beginnings of a problem." << std::endl;
567  Hnew[0][0] = Hnew[0][0] - eps/Gain[0][0];
568  }
569  TMatrixT<Double_t> GHc(GH*cov);
570 
571  cov-=Gain*(Hnew*cov);
572 
573  // Below is protection required at end of contained track when
574  // momentum is tiny and cov[0][0] gets huge.
575 
576  // if (cov[0][0]>pi2/Hnew[0][0] || cov[0][0] <= 0.0 || TMath::IsNaN(cov[0][0]))
577  if (cov[0][0] <= 0.0 || TMath::IsNaN(cov[0][0]))
578  {
579  cov[0][0]=pi2/Hnew[0][0]; // some big value.a
580  //cov = covFilt;
581  }
582  else // store away a good cov matrix
583  {
584  covFilt = cov;
585  }
586 
587 
588  TMatrixT<double> cov7x7(calcCov7x7(cov,pl));
589  tr->setHitCov(cov);
590  tr->setHitCov7x7(cov7x7);
591  tr->setHitState(state);
592 
593  if (fabs(1.0/state[0][0])<fMomLow)
594  state[0][0] = 1.0/fMomLow*fabs(state[0][0])/state[0][0];
595  if (fabs(1.0/state[0][0])>fMomHigh)
596  state[0][0] = 1.0/fMomHigh*fabs(state[0][0])/state[0][0];
597 
598 
599  // Let's also calculate the "filtered" plane, and pointer here.
600  // Use those to calculate filtered chisq and pass to setHitPlane()
601  // for plane-by-plane correct "filtered" track pointing that gets
602  // stuffed into the Track().
603  TVector3 uf(pl.getU());
604  TVector3 vf(pl.getV());
605  TVector3 wf(uf.Cross(vf));
606  TVector3 Of(pl.getO());
607  // direction *
608  //Double_t sign2(1.0);
609  //if (direction==-1 && ihit==(int)nhits-1) sign2 = -1.0;
610 
611  TVector3 pf = direction*(wf + state[1][0] * uf + state[2][0] * vf);
612  TVector3 pposf = Of + state[3][0] * uf + state[4][0] * vf;
613  Double_t angf = TMath::Min(fabs(pf.Angle(wf)),0.95*TMath::Pi()/2.0);
614  TVector3 rotf(wf.Cross(pf.Unit()));
615 
616  uf.Rotate(angf,rotf);
617  vf.Rotate(angf,rotf);
618  wf = uf.Cross(vf);
619  plFilt.setU(uf.Unit());
620  plFilt.setV(vf.Unit());
621  plFilt.setO(pposf);
622  plFilt.setNormal(pf.Unit());
623 
624  tr->setHitPlaneXYZ(plFilt.getO());
625  tr->setHitPlaneUxUyUz(plFilt.getNormal());
626  tr->setHitPlaneU(plFilt.getU());
627  tr->setHitPlaneV(plFilt.getV());
628 
629  // calculate filtered chisq from filtered residuals
630  TMatrixT<Double_t> r=hit->residualVector(rep,state,plFilt,plPrev,mass);
631  if (direction==-1) wold.Rotate(TMath::Pi(),wold.Orthogonal());
632  dtheta = thetaMeas - TMath::Min(fabs(wold.Angle(plFilt.getNormal())),0.95*TMath::Pi()/2.);
633  r[0][0] = dtheta/Hnew[0][0];
634  //r[1][0] = (pointer*uPrev)/(pointer*wPrev) - (wf*uPrev)/(wf*wPrev);
635  //r[2][0] = (pointer*vPrev)/(pointer*wPrev) - (wf*vPrev)/(wf*wPrev);
636  //r[1][0] = 0.;
637  //r[2][0] = 0.;
638 
639  double chi2 = chi2Increment(r,Hnew,cov,V);
640  int ndf = r.GetNrows();
641  if (update[0][0]==0.0) {chi2=0.0; ndf=0;};
642  rep->addChiSqu( chi2 );
643  rep->addNDF( ndf );
644  pointerPrev = pointer;
645  tr->setHitChi2(chi2);
646  /*
647  if(direction==1){
648  tr->getBK(irep)->setNumber("fChi2",ihit,chi2/ndf);
649  }
650  else{
651  tr->getBK(irep)->setNumber("bChi2",ihit,chi2/ndf);
652  }
653  */
654 
655  // if we survive until here: update TrackRep
656  //rep->setState(state);
657  //rep->setCov(cov);
658  //rep->setReferencePlane(pl);
659 
660  // Since I've tilted the planes, d(u,v)/dw=0 by construction.
661  // But, I *haven't* tilted the planes! EC, 11-May-2012.
662  //state[1][0] = 0.0;
663  //state[2][0] = 0.0;
664  rep->setData(state,pl,&cov); // This is where fState,
665  // fRefPlane and fCov are updated.
666  tr->setRepAtHit(irep,ihit);
667 
668 }
669 
670 TMatrixT<Double_t>
671 genf::GFKalman::calcCov7x7(const TMatrixT<Double_t>& cov, const GFDetPlane& plane)
672 {
673  // This ends up, confusingly, as: 7 columns, 5 rows!
674  TMatrixT<Double_t> jac(7,5); // X,Y,Z,UX,UY,UZ,Theta in detector coords
675 
676  TVector3 u=plane.getU();
677  TVector3 v=plane.getV();
678  TVector3 w=u.Cross(v);
679 
680  // Below line has been done, with code change in GFKalman that has updated
681  // the plane orientation by now.
682  // TVector3 pTilde = 1.0 * (w + state[1][0] * u + state[2][0] * v);
683  TVector3 pTilde = w;
684  double pTildeMag = pTilde.Mag();
685 
686  jac.Zero();
687  jac[6][0] = 1.; // Should be C as in GFSpacepointHitPolicy. 16-Feb-2013.
688 
689  jac[0][3] = u[0];
690  jac[1][3] = u[1];
691  jac[2][3] = u[2];
692 
693  jac[0][4] = v[0];
694  jac[1][4] = v[1];
695  jac[2][4] = v[2];
696 
697  // cnp'd from RKTrackRep.cxx, line ~496
698  // da{x,y,z}/du'
699  jac[3][1] = 1.0/pTildeMag*(u.X()-pTilde.X()/(pTildeMag*pTildeMag)*u*pTilde);
700  jac[4][1] = 1.0/pTildeMag*(u.Y()-pTilde.Y()/(pTildeMag*pTildeMag)*u*pTilde);
701  jac[5][1] = 1.0/pTildeMag*(u.Z()-pTilde.Z()/(pTildeMag*pTildeMag)*u*pTilde);
702  // da{x,y,z}/dv'
703  jac[3][2] = 1.0/pTildeMag*(v.X()-pTilde.X()/(pTildeMag*pTildeMag)*v*pTilde);
704  jac[4][2] = 1.0/pTildeMag*(v.Y()-pTilde.Y()/(pTildeMag*pTildeMag)*v*pTilde);
705  jac[5][2] = 1.0/pTildeMag*(v.Z()-pTilde.Z()/(pTildeMag*pTildeMag)*v*pTilde);
706 
707  // y = A.x => x = A^T.A.A^T.y
708  // Thus, y's Jacobians Jac become for x (Jac^T.Jac)^(-1) Jac^T
709 
710  TMatrixT<Double_t> jac_t(TMatrixT<Double_t>::kTransposed,jac);
711  TMatrixT<Double_t> jjInv(jac_t * jac);
712  //jInv.Zero();
713 
714  double det(0.0);
715  try {
716  jjInv.Invert(&det); // this is all 1s on the diagonal, perhaps
717  // to no one's surprise.
718  }
719  catch (cet::exception &)
720  {
721  throw GFException("GFKalman: Jac.T*Jac is not invertible. But keep plowing on ... ",__LINE__,__FILE__).setFatal();
722  }
723  if(TMath::IsNaN(det)) {
724  throw GFException("GFKalman: det of Jac.T*Jac is nan",__LINE__,__FILE__).setFatal();
725  }
726 
727  TMatrixT<Double_t> j5x7 = jjInv*jac_t;
728  TMatrixT<Double_t> c7x7 = j5x7.T() * (cov * j5x7);
729  return c7x7;
730 }
731 
732 TMatrixT<Double_t>
733 genf::GFKalman::calcGain(const TMatrixT<Double_t>& cov,
734  const TMatrixT<Double_t>& HitCov,
735  const TMatrixT<Double_t>& H){
736 
737  // calculate covsum (V + HCH^T)
738 
739  // Comment next 3 out for normal running.
740  //cov.Print();
741  //H.Print();
742  //HitCov.Print();
743  TMatrixT<Double_t> covsum1(cov,TMatrixT<Double_t>::kMultTranspose,H);
744  // covsum1.Print();
745  TMatrixT<Double_t> covsum(H,TMatrixT<Double_t>::kMult,covsum1);
746 
747  //covsum.Print();
748 
749  covsum+=HitCov;
750  //covsum = (covsum,TMatrixT<Double_t>::kPlus,HitCov);
751  // covsum.Print();
752 
753  // invert
754  double det=0;
755  covsum.SetTol(1.0e-23); // to avoid inversion problem, EC, 8-Aug-2011.
756  covsum.Invert(&det);
757  // std::cout << "GFKalman:: calGain(), det is "<< det << std::endl;
758  if(TMath::IsNaN(det)) {
759  throw GFException("Kalman Gain: det of covsum is nan",__LINE__,__FILE__).setFatal();
760  }
761 
762  if(det==0){
763  GFException exc("cannot invert covsum in Kalman Gain - det=0",
764  __LINE__,__FILE__);
765  exc.setFatal();
766  std::vector< TMatrixT<Double_t> > matrices;
767  matrices.push_back(cov);
768  matrices.push_back(HitCov);
769  matrices.push_back(covsum1);
770  matrices.push_back(covsum);
771  exc.setMatrices("cov, HitCov, covsum1 and covsum",matrices);
772  throw exc;
773 
774  }
775 
776  // gain is CH^T/(V + HCH^T)
777  // calculate gain
778  TMatrixT<Double_t> gain1(H,TMatrixT<Double_t>::kTransposeMult,covsum);
779  TMatrixT<Double_t> gain(cov,TMatrixT<Double_t>::kMult,gain1);
780 
781  return gain;
782 }
void clearRepAtHit()
clear the hit indices at which plane,state&cov of reps are defined
Definition: GFTrack.h:405
void setNDF(unsigned int n)
GFException & setNumbers(std::string, const std::vector< double > &)
set list of numbers with description
Definition: GFException.cxx:34
double getChi2Hit(GFAbsRecoHit *, GFAbsTrackRep *)
Calculates chi2 of a given hit with respect to a given track representation.
Definition: GFKalman.cxx:276
double beta(double KE, const simb::MCParticle *part)
const GFDetPlane & getReferencePlane() const
Double_t fErrScaleMTh
Definition: GFKalman.h:163
std::string string
Definition: nybbler.cc:12
void setNextHitToFit(unsigned int i)
Set next hit to be used in a fit.
Definition: GFTrack.h:192
void setRepAtHit(unsigned int irep, int ihit)
set the hit index at which plane,state&cov of rep irep is defined
Definition: GFTrack.h:389
void processHit(GFTrack *, int, int, int)
One Kalman step.
Definition: GFKalman.cxx:303
virtual const char * what() const
standard error message handling for exceptions. use like "std::cerr << e.what();" ...
Definition: GFException.cxx:48
TMatrixT< Double_t > calcGain(const TMatrixT< Double_t > &cov, const TMatrixT< Double_t > &HitCov, const TMatrixT< Double_t > &H)
Calculate Kalman Gain.
Definition: GFKalman.cxx:733
GFAbsTrackRep * getTrackRep(int id) const
Accessor for track representations.
Definition: GFTrack.h:196
void fittingPass(GFTrack *, int dir)
Performs fit on a GFTrack beginning with the current hit.
Definition: GFKalman.cxx:177
virtual TMatrixT< Double_t > getHMatrix(const GFAbsTrackRep *stateVector)=0
Get transformation matrix. Transformation between hit coordinates and track representation coordinate...
virtual const GFDetPlane & getDetPlane(GFAbsTrackRep *)=0
Get detector plane for a given track representation.
TMatrixT< Double_t > getRawHitCoord() const
Get raw hit coordinates.
Definition: GFAbsRecoHit.h:183
const TMatrixT< Double_t > & getState() const
void blowUpCovs(GFTrack *trk)
this is needed to blow up the covariance matrix before a fitting pass drops off-diagonal elements and...
Definition: GFKalman.cxx:154
Base Class for genfit track representations. Defines interface for track parameterizations.
Definition: GFAbsTrackRep.h:83
void setCov(const TMatrixT< Double_t > &aCov)
int fInitialDirection
Definition: GFKalman.h:155
void setChiSqu(double aChiSqu)
Double_t fMomLow
Definition: GFKalman.h:159
void setHitState(TMatrixT< Double_t > mat)
Definition: GFTrack.h:354
Double_t fMaxUpdate
Definition: GFKalman.h:161
void addNDF(unsigned int n)
TVector3 getO() const
Definition: GFDetPlane.h:81
const double e
void switchDirection(GFTrack *trk)
Used to switch between forward and backward filtering.
Definition: GFKalman.cxx:125
unsigned int getNumReps() const
Get number of track represenatations.
Definition: GFTrack.h:202
void setHitMeasuredCov(TMatrixT< Double_t > mat)
Definition: GFTrack.h:351
virtual void setData(const TMatrixT< Double_t > &st, const GFDetPlane &pl, const TMatrixT< Double_t > *cov=NULL)
double chi2Increment(const TMatrixT< Double_t > &r, const TMatrixT< Double_t > &H, const TMatrixT< Double_t > &cov, const TMatrixT< Double_t > &V)
this returns the reduced chi2 increment for a hit
Definition: GFKalman.cxx:222
void addFailedHit(unsigned int irep, unsigned int id)
Definition: GFTrack.h:294
virtual TMatrixT< Double_t > getHitCov(const GFDetPlane &)=0
Get hit covariances in a specific detector plane.
bool isFatal()
get fatal flag.
Definition: GFException.h:80
void addChiSqu(double aChiSqu)
GFAbsRecoHit * getHit(int id) const
Definition: GFTrack.h:160
int getRepAtHit(unsigned int irep)
get the hit index at which plane,state&cov of rep irep is defined
Definition: GFTrack.h:397
void setHitPlaneV(TVector3 pl)
Definition: GFTrack.h:360
void setHitPlaneUxUyUz(TVector3 pl)
Definition: GFTrack.h:358
virtual double extrapolate(const GFDetPlane &plane, TMatrixT< Double_t > &statePred)
returns the tracklength spanned in this extrapolation
void processTrack(GFTrack *)
Performs fit on a GFTrack.
Definition: GFKalman.cxx:46
void setHitChi2(Double_t mat)
Definition: GFTrack.h:353
const TMatrixT< Double_t > & getCov() const
GFBookkeeping * getBK(int index=-1)
get GFBookKeeping object for particular track rep (default is cardinal rep)
Definition: GFTrack.h:335
GFException & setMatrices(std::string, const std::vector< TMatrixT< Double_t > > &)
set list of matrices with description
Definition: GFException.cxx:41
void setHitCov7x7(TMatrixT< Double_t > mat)
Definition: GFTrack.h:356
void blowUpCovsDiag(GFTrack *trk)
Definition: GFKalman.cxx:132
Detector simulation of raw signals on wires.
void info()
print information in the exception object
Definition: GFException.cxx:58
double fBlowUpFactor
Definition: GFKalman.h:157
unsigned int getNumHits() const
Definition: GFTrack.h:164
void setHitCov(TMatrixT< Double_t > mat)
Definition: GFTrack.h:355
bool fGENfPRINT
Definition: GFKalman.h:164
Double_t fMomHigh
Definition: GFKalman.h:160
Exception class for error handling in GENFIT (provides storage for diagnostic information) ...
Definition: GFException.h:48
TVector3 getU() const
Definition: GFDetPlane.h:82
constexpr double dist(const TReal *x, const TReal *y, const unsigned int dimension)
void setStatusFlag(int _val)
void setHitUpdate(TMatrixT< Double_t > mat)
Definition: GFTrack.h:352
TMatrixT< Double_t > calcCov7x7(const TMatrixT< Double_t > &cov, const genf::GFDetPlane &plane)
Definition: GFKalman.cxx:671
GFException & setFatal(bool b=true)
set fatal flag. if this is true, the fit stops for this current track repr.
Definition: GFException.h:78
#define COVEXC
Definition: GFKalman.cxx:36
void setHitPlaneXYZ(TVector3 pl)
Definition: GFTrack.h:357
int getPDG()
Definition: GFTrack.h:373
TVector3 getV() const
Definition: GFDetPlane.h:83
unsigned int getDim() const
returns dimension of state vector
unsigned int getNextHitToFit() const
Accessor for fNextHitToFit.
Definition: GFTrack.h:188
virtual TMatrixT< Double_t > residualVector(const GFAbsTrackRep *stateVector, const TMatrixT< Double_t > &state, const GFDetPlane &d)
Calculate residual with respect to a track representation.
Definition: GFAbsRecoHit.h:148
virtual void switchDirection()=0
void setHitPlaneU(TVector3 pl)
Definition: GFTrack.h:359
cet::coded_exception< error, detail::translate > exception
Definition: exception.h:33
QTextStream & endl(QTextStream &s)
Int_t fNumIt
Definition: GFKalman.h:156