FairRoot/PandaRoot
RhoKinHyperonFitter.cxx
Go to the documentation of this file.
1 #include <iostream>
2 #include "RhoKinHyperonFitter.h"
4 #include "RhoBase/RhoFactory.h"
5 #include "TDecompLU.h"
6 #include "TMatrixD.h"
7 #include "TMatrixDSym.h"
9 using namespace std;
10 
12 
13 TBuffer& operator>>(TBuffer& buf, RhoKinHyperonFitter *&obj)
14 {
15  obj = (RhoKinHyperonFitter*) buf.ReadObject(RhoKinHyperonFitter::Class());
16  return buf;
17 }
18 
20  RhoFitterBase( b )
21 {
22  f4MomConstraint=-1;
23  fMomConstraint=-1;
24  fMassConstraint =-1;
26  fTotEConstraint=-1;
27 }
28 
30 {
31 }
32 
33 
35 {
36  f4MomConstraint = 1;
37  flmm=lv;
38 }
39 
41 {
42  fMomConstraint = 1;
43  fmm=v;
44 }
45 
47 {
48  fTotEConstraint = 1;
49  fEc=energy;
50 }
52 {
54  fMom=momentum;
55 }
57 {
58  fMassConstraint = 1;
59  fMass=mass;
60 }
61 
63 {
64  fDaughters.clear();
65  //FindAndAddFinalStateDaughters(fHeadOfTree);
67  //int nd=fDaughters.size(); //unused?
68 
69  fNumCon=0;
70  if(f4MomConstraint >0) {
71  int n4Mom = 4;
72  fNumCon = fNumCon + n4Mom;
73  }
74  if(fMomConstraint >0) {
75  int nMom = 3;
76  fNumCon = fNumCon + nMom;
77  }
78  if(fTotEConstraint >0) {
79  int nE = 1;
80  fNumCon = fNumCon + nE;
81  }
82  if(fMassConstraint >0) {
83  int nMass = 1;
84  fNumCon = fNumCon + nMass;
85  }
86  if(fTotMomConstraint >0) {
87  int nTotMom = 1;
88  fNumCon = fNumCon + nTotMom;
89  }
90 
92 {
93  SetMatrices();
94  ZeroMatrices();
95  ReadMatrix();
96 }
97  Bool_t check = Solve();
98  SetOutput();
99  return check;
100 }
101 
103 {
104  int nd;
105 
106  if(fMassConstraint>0 || fMomConstraint>0) {nd=1;}
107  else{
108  fDaughters.clear();
109  //FindAndAddFinalStateDaughters(fHeadOfTree);
111  nd=fDaughters.size();
112  //cout<<"SetMatrix nd="<<nd<<endl;
113  }
114 
115  fNvar=7;
116  fNpar = nd*fNvar;
117  fNcon= fNumCon;
118  fNc=0;
119  fNiter=0;
120 
121 
122  fAl0.ResizeTo(7*nd,1);
123  fV_al0.ResizeTo(fNpar,fNpar);
124  fAl1.ResizeTo(7*nd,1);
125  fV_al1.ResizeTo(fNpar,fNpar);
126  fmD.ResizeTo(fNcon,fNpar);
127  fmE.ResizeTo(fNcon,3);
128  fmd.ResizeTo(fNcon,1);
129  fmPull.ResizeTo(7*nd,1);
130 }
131 
132 
134 {
135  fAl0.Zero();
136  fV_al0.Zero();
137  fAl1.Zero();
138  fV_al1.Zero();
139  fmPull.Zero();
140  fmD.Zero();
141  fmd.Zero();
142  fmE.Zero();
143 }
144 
145 
147 {
148  double ierr; // used to check inversions
149  fAl1=fAl0;
150  TMatrixD chi2_0(1,1);
151  chi2_0[0][0]=10000000.;// Xinying Song added a very large initial chi2
152  int j1Max=40;
153 
154  for(Int_t j1=0;j1<j1Max;++j1)
155  {
157  if(fMomConstraint >0) { ReadMomKinMatrix();}
161 
162 
163  TMatrixD mD_t=fmD;
164  mD_t.T();
165 
166  TMatrixD Vd_inv = fmD*fV_al0*mD_t;
167  TMatrixD Vd = Vd_inv.Invert(&ierr);
168  TMatrixD lam = Vd* ( fmD*(fAl0 - fAl1) + fmd);
169  TMatrixD lam_t=lam;
170  lam_t.T();
171  TMatrixD al_new=fAl0-fV_al0*mD_t*lam;
172  TMatrixD V_al_new=fV_al0-fV_al0*mD_t*Vd*fmD*fV_al0;
173 // TMatrixD chi2_new = lam_t*(fmD*(fAl0-fAl1)+fmd );
174  TMatrixD chi2_new = lam_t*(fmD*(fAl0-al_new)+fmd );
175 
176  if(TMath::IsNaN(chi2_new[0][0])) continue;
177  double delchi2 = chi2_0[0][0]- chi2_new[0][0];
178 
179 
180  if( chi2_new[0][0] < chi2_0[0][0]|| fabs(delchi2)<1.){
181  fAl1= al_new;
182  // TMatrixD V_al_new=fV_al0-fV_al0*mD_t*Vd*fmD*fV_al0;
183  // fV_al1= V_al_new;
184  chi2_0[0][0]= chi2_new[0][0];
185  fChiSquare=chi2_new[0][0];
186  if(fabs(delchi2)<0.01||j1==j1Max-1){ fV_al1= V_al_new; break;}
187  }
188 }
189 
190 fAl0=fAl1;
191 fV_al0=fV_al1;
192 return kTRUE;
193 
194 
195 }
196 
197 //Write output
199 {
200 
201 
202  fNDegreesOfFreedom = 0;
203  if(f4MomConstraint >0)
204  //fNDegreesOfFreedom += 4+fDaughters.size();
205  fNDegreesOfFreedom += 4;
206  else {
207  if(fMomConstraint >0)
208  fNDegreesOfFreedom += 3;
210  fNDegreesOfFreedom += 1;
211  }
212 
213 
214 
216  int nd=1;
217  TMatrixD m(nd,1);
218  TVector3 pos(fAl0[4][0],fAl0[5][0],fAl0[6][0]);
219  TLorentzVector mom4(fAl0[0][0],fAl0[1][0],fAl0[2][0],fAl0[3][0]);
220  fHeadOfTree->SetP7(pos,mom4);
221 
222 
223  TMatrixD p1Cov(7,7);
224  TMatrixD HeadCov(7,7);
225  for(int i=0; i<7; i++) {
226  for (int j=0; j<7; j++) {
227 
228  p1Cov[i][j]= fV_al0[i][j];
229  //fHeadOfTree->SetCov7(p1Cov); //New covariance matrix without correlations
230  }
231  }
232 
233 //Change from px,py,pz,E,x,y,z
234  // to x,y,z,px,py,pz,E
235  for(int i=0; i<7; i++) {
236  for(int j=0; j<7; j++) {
237  if(i>=4) {
238  if(j>=4) {
239  HeadCov[i-4][j-4] = p1Cov[i][j];
240  } else { HeadCov[i-4][j+3] = p1Cov[i][j]; }
241  } else {
242  if(j>=4) {
243  HeadCov[i+3][j-4] = p1Cov[i][j];
244  } else { HeadCov[i+3][j+3] = p1Cov[i][j]; }
245  }
246  }
247  }
248 
249 
250 fHeadOfTree->SetCov7(HeadCov);
251  return;
252 
253 }else{//in the 4C case, we need to output the p7 and cov7 of the daughters
254  int nd=fDaughters.size();
255  TMatrixD m(nd,1);
256  TLorentzVector tot(0,0.0,0.0,0.0);
257  for (int k=0; k<nd; k++) {
258  TLorentzVector p1=fDaughters[k]->P4();
259  m[k][0]=p1.M();
260  TVector3 pos(fAl0[k*7+4][0],fAl0[k*7+5][0],fAl0[k*7+6][0]);
261  double E =sqrt(fAl0[k*7+0][0]*fAl0[k*7+0][0]+
262  fAl0[k*7+1][0]*fAl0[k*7+1][0]+
263  fAl0[k*7+2][0]*fAl0[k*7+2][0]+
264  m[k][0]*m[k][0]);
265 
266  TLorentzVector mom4(fAl0[k*7+0][0],fAl0[k*7+1][0],fAl0[k*7+2][0],E);
267  fDaughters[k]->SetP7(pos,mom4);
268 
269 
270  TMatrixD p1Cov(7,7);
271  TMatrixD p2Cov(7,7);
272  for(int i=0; i<7; i++) {
273  for (int j=0; j<7; j++) {
274 
275  p1Cov[i][j]= fV_al0[k*7+i][k*7+j];
276  // fDaughters[k]->SetCov7(p1Cov); //New covariance matrix without correlations
277  }
278  }
279 
280 //Change from px,py,pz,E,x,y,z
281 // to x,y,z,px,py,pz,E
282 
283  for(int i=0; i<7; i++) {
284  for(int j=0; j<7; j++) {
285  if(i>=4) {
286  if(j>=4) {
287  p2Cov[i-4][j-4] = p1Cov[i][j];
288  } else { p2Cov[i-4][j+3] = p1Cov[i][j]; }
289  } else {
290  if(j>=4) {
291  p2Cov[i+3][j-4] = p1Cov[i][j];
292  } else { p2Cov[i+3][j+3] = p1Cov[i][j]; }
293  }
294  }
295  }
296 tot+=mom4;
297 fDaughters[k]->SetCov7(p2Cov);
298  }
299 fHeadOfTree->SetP4(tot);
300 return;
301  }
302 
303 
304 }
305 
306 
307 
308 //Read the input vector
310 {
311 
312 //Here the m.c.f is treated differently with other fit, e.g. 4c
313 
315 
316 //px,py,pz,E,x,y,z
317  TLorentzVector p1=fHeadOfTree->P4();
318  TVector3 p2=fHeadOfTree->Pos();
319  fAl0[0][0]=p1.X();
320  fAl0[1][0]=p1.Y();
321  fAl0[2][0]=p1.Z();
322  fAl0[3][0]=p1.E();
323  fAl0[4][0]=p2.X();
324  fAl0[5][0]=p2.Y();
325  fAl0[6][0]=p2.Z();
326 
327 
328 // Read Covariance Matrix .... Can read 6x6 matrices..................
329  TMatrixD p1Cov(7,7);
330  TMatrixD HeadCov(7,7);
331  TMatrixD p3Cov(7,7);
332  TMatrixD p4Cov(7,7);
333  p1Cov=fHeadOfTree->Cov7(); //Cov Matrix x,y,z,px,py,pz,E
334 
335  for (int ii=0; ii<6; ii++) {for(int jj=0; jj<6; jj++) {HeadCov[ii][jj]=p1Cov[ii][jj];}} //test
336 
337  //Extend matrix for energy for each candidates .....6x6 to 7x7
338 
339  p3Cov=p1Cov;
340 //cout<<"in readMatrix before trans"<<endl;
341 // p3Cov.Print();
342 //Change to px,py,pz,E,x,y,z
343  for(int i=0; i<7; i++) {
344  for(int j=0; j<7; j++) {
345  if(i>=3) {
346  if(j>=3) {
347  p4Cov[i-3][j-3] = p3Cov[i][j];
348  } else { p4Cov[i-3][j+4] = p3Cov[i][j]; }
349  } else {
350  if(j>=3) {
351  p4Cov[i+4][j-3] = p3Cov[i][j];
352  } else { p4Cov[i+4][j+4] = p3Cov[i][j]; }
353  }
354  }
355  }
356 
357  //cout<<"after trans"<<endl;
358  // p4Cov.Print();
359  for(int i=0; i<7; i++) {
360  for (int j=0; j<7; j++) {
361  fV_al0[i][j] = p4Cov[i][j];
362  }
363  }
364 
365  }
366 
367  }
368 
369 
370 
372 {
373 
374 // Here these Matrix must be read in from the vertex fit
375  fNumCon=0;
376  if(f4MomConstraint >0) {
377  int n4Mom = 4;
378  fNumCon = fNumCon + n4Mom;
379  }
380  if(fMomConstraint >0) {
381  int nMom = 3;
382  fNumCon = fNumCon + nMom;
383  }
384  if(fTotEConstraint >0) {
385  int nE = 1;
386  fNumCon = fNumCon + nE;
387  }
388  if(fMassConstraint >0) {
389  int nMass = 1;
390  fNumCon = fNumCon + nMass;
391  }
392  if(fTotMomConstraint >0) {
393  int nTotMom = 1;
394  fNumCon = fNumCon + nTotMom;
395  }
396 
397 
398 SetMatrices();
399 ZeroMatrices();
400 fAl0=M1;
401 fV_al0=V1;
402 
403 }
404 
405 
406 
407 
408 
409 
411 {
412 
413  fmD.Zero();
414  //double Etot = 0.;
415  //double Px = 0.;
416  //double Py = 0.;
417  //double Pz = 0.;
418  TMatrixD al1p(fAl1);
419  //double a;
420 
421 
422  double px = al1p[0][0];
423  double py = al1p[1][0];
424  double pz = al1p[2][0];
425  double E = al1p[3][0];
426 
427  fmd[0][0] = E*E - px*px - py*py - pz*pz - fMass*fMass ;
428 
429 // V.J. - force head mass to constraint mass
430 //....................................................
431  fmD[0][0] = -2.*px;
432  fmD[0][1] = -2.*py;;
433  fmD[0][2] = -2.*pz;
434  fmD[0][3] = 2.*E;
435  fmD[0][4] = 0.0;
436  fmD[0][5] = 0.0;
437  fmD[0][6] = 0.0;
438 //.....................................
439 
440 }
441 
443 {
444  int nd=fDaughters.size();
445  TMatrixD alp(fAl1);
446  TMatrixD m(nd,1);
447  fmD.Zero();
448  fmd.Zero();
449  int k,i;
450  for (k=0; k<nd; k++) {
451 
452  TLorentzVector p1=fDaughters[k]->P4();
453  m[k][0]=p1.M();
454 
455  // double E =alp[k*7+3][0];
456  double E =sqrt(alp[k*7+0][0]*alp[k*7+0][0]+
457  alp[k*7+1][0]*alp[k*7+1][0]+
458  alp[k*7+2][0]*alp[k*7+2][0]+
459  m[k][0]*m[k][0]);
460 
461  for(i=0; i<3; i++) {
462  fmD[fNc+i][k*7+i] = 1;
463  fmD[fNc+3][k*7+i] = alp[k*7+i][0]/E;
464  //if(i==3){
465  //fmD[fNc][k*7+i]=E/alp[k*7][0];
466  //fmD[fNc+1][k*7+i]=E/alp[k*7+1][0];
467  //fmD[fNc+2][k*7+i]=E/alp[k*7+2][0];
468  //fmD[fNc+3][k*7+i]=1;
469  //}
470 
471  }
472 
473 // cout << "value of D" << mD[3][k*7+i] << endl;
474 
475  for (i=0; i<3; i++) {
476  fmd[fNc+i][0] += alp[k*7+i][0];
477  }
478  fmd[fNc+3][0] += E;
479  }
480 
481  fmd[fNc+0][0] -= flmm.X();
482  fmd[fNc+1][0] -= flmm.Y();
483  fmd[fNc+2][0] -= flmm.Z();
484  fmd[fNc+3][0] -= flmm.T();
485 
486 // fNc += 4;
487 }
488 
489 
491 {
492  int nd=fDaughters.size();
493  TMatrixD alp(fAl1);
494  TMatrixD m(nd,1);
495  int k,i;
496  for (k=0; k<nd; k++) {
497  for (i=0; i<3; i++) {
498  fmD[fNc+i][k*7+i] = 1;
499  }
500  for (i=0; i<3; i++) {
501  fmd[fNc+i][0] += alp[k*7+i][0];
502  }
503  }
504  fmd[fNc+0][0]-= fmm.X();
505  fmd[fNc+1][0]-= fmm.Y();
506  fmd[fNc+2][0]-= fmm.Z();
507 
508  // fNc += 3;
509 }
510 
511 
513 {
514  int nd=fDaughters.size();
515  TMatrixD alp(fAl1);
516  TMatrixD m(nd,1);
517  int k,i;
518  for (k=0; k<nd; k++) {
519  TLorentzVector p1=fDaughters[k]->P4();
520  m[k][0]=p1.M();
521  double E =sqrt(alp[k*7+0][0]*alp[k*7+0][0]+
522  alp[k*7+1][0]*alp[k*7+1][0]+
523  alp[k*7+2][0]*alp[k*7+2][0]+
524  m[k][0]*m[k][0]);
525  for (i=0; i<3; i++) {
526  fmD[fNc+0][k*7+i] = 0.;
527  // mD[0][k*7+i] = alp[k*7+i][0]/E;
528  }
529  fmD[fNc+0][k*7+3] = 1.;
530  fmd[fNc+0][0] += E;
531  }
532  fmd[fNc+0][0] -= fEc;
533  //fNc +=1;
534 }
535 
537 {
538  int nd=fDaughters.size();
539  TMatrixD alp(fAl1);
540  TMatrixD m(nd,1);
541  int k,i;
542  for (k=0; k<nd; k++) {
543  double Ptot =sqrt(alp[k*7+0][0]*alp[k*7+0][0]+
544  alp[k*7+1][0]*alp[k*7+1][0]+
545  alp[k*7+2][0]*alp[k*7+2][0]);
546  for (i=0; i<3; i++) {
547  // mD[0][k*7+i] = 1.;
548  fmD[fNc+0][k*7+i] = alp[k*7+i][0]/Ptot;
549  }
550  fmd[fNc+0][0] += Ptot;
551  }
552  fmd[fNc+0][0] -= fMom;
553  // fNc +=1;
554 }
555 
556 /*
557 void RhoKinHyperonFitter::ReadEqMassKinMatrix()
558 {int nd=fDaughters.size();
559  TMatrixD alp(al1);
560  TMatrixD m(nd,1);
561  int k,i,j;
562  for (k=0;k<nd;k++)
563 {
564  TLorentzVector p1=fDaughters[k]->P4();
565  double Mtot =sqrt(alp[k*7+3][0]*alp[k*7+3][0]
566  -(alp[k*7+0][0]*alp[k*7+0][0]+
567  alp[k*7+1][0]*alp[k*7+1][0]+
568  alp[k*7+2][0]*alp[k*7+2][0])
569  for (i=0;i<3;i++)
570  {
571  // mD[0][k*7+i] = 1.;
572  mD[fNc+0][k*7+i] -= 2.*alp[k*7+i][0];
573  mD[fNc+0][k*7+3] -= -2.*alp[k*7+3][0];
574  }
575  md[fNc+0][0] -= Mtot;
576  }
577  md[fNc+0][0] -= 0.0;
578  //fNc +=1;
579 }
580 */
581 
TVector3 pos
RhoKinHyperonFitter(RhoCandidate *b)
void SetP7(const TVector3 &pos, const TLorentzVector &p4)
void SetInputMatrix(TMatrixD M1, TMatrixD V1)
Int_t i
Definition: run_full.C:25
__m128 m
Definition: P4_F32vec4.h:28
TTree * b
friend F32vec4 sqrt(const F32vec4 &a)
Definition: P4_F32vec4.h:29
TVector3 Pos() const
Definition: RhoCandidate.h:186
void Add4MomConstraint(TLorentzVector lv)
__m128 v
Definition: P4_F32vec4.h:4
void SetP4(Double_t mass, const TVector3 &p3)
TLorentzVector P4() const
Definition: RhoCandidate.h:195
friend F32vec4 fabs(const F32vec4 &a)
Definition: P4_F32vec4.h:47
std::vector< RhoCandidate * > fDaughters
Definition: RhoFitterBase.h:69
TPad * p2
Definition: hist-t7.C:117
TMatrixD Cov7() const
void AddTotMomConstraint(double momentum)
RhoCandidate * fHeadOfTree
Definition: RhoFitterBase.h:62
void AddTotEConstraint(double energy)
double fChiSquare
Definition: RhoFitterBase.h:74
ClassImp(PndAnaContFact)
TPad * p1
Definition: hist-t7.C:116
void AddMomConstraint(TVector3 v)
void SetDaugthersFromComposite(RhoCandidate *cand)
void SetCov7(const TMatrixD &cov7)
void AddMassConstraint(double mass)
TMatrixT< double > TMatrixD
Definition: PndLmdDim.h:52
double pz[39]
Definition: pipisigmas.h:14
Double_t energy
Definition: plot_dirc.C:15
int fNDegreesOfFreedom
Definition: RhoFitterBase.h:75