]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/STEER/AliMillePede2.cxx
Changes for #92834: Implement trigger aliases in AliRoot, RAW data reconstruction
[u/mrichter/AliRoot.git] / STEER / STEER / AliMillePede2.cxx
1 /**********************************************************************************************/
2 /* General class for alignment with large number of degrees of freedom                        */
3 /* Based on the original milliped2 by Volker Blobel                                           */
4 /* and AliMillepede class by Javier                                                           */ 
5 /* Allows operations with large sparse matrices                                               */ 
6 /* http://www.desy.de/~blobel/mptalks.html                                                    */
7 /*                                                                                            */ 
8 /* Author: ruben.shahoyan@cern.ch                                                             */
9 /*                                                                                            */ 
10 /**********************************************************************************************/
11
12 #include "AliMillePede2.h"
13 #include "AliLog.h"
14 #include <TStopwatch.h>
15 #include <TFile.h>
16 #include <TChain.h>
17 #include <TMath.h>
18 #include <TVectorD.h>
19 #include <TArrayL.h>
20 #include <TSystem.h>
21 #include "AliMatrixSq.h"
22 #include "AliSymMatrix.h"
23 #include "AliRectMatrix.h"
24 #include "AliMatrixSparse.h"
25 #include <stdio.h>
26 #include <stdlib.h>
27 #include <unistd.h> 
28 #include <sys/types.h>
29 #include <sys/stat.h>
30 #include <fcntl.h>
31 #include <fstream>
32
33 using namespace std;
34
35
36 ClassImp(AliMillePede2)
37
38 Bool_t   AliMillePede2::fgInvChol        = kTRUE;     // Invert global matrix with Cholesky solver
39 Bool_t   AliMillePede2::fgWeightSigma    = kTRUE;     // weight local constraint by module statistics
40 Bool_t   AliMillePede2::fgIsMatGloSparse = kFALSE;    // use faster dense matrix by default
41 Int_t    AliMillePede2::fgMinResCondType = 1;         // Jacoby preconditioner by default
42 Double_t AliMillePede2::fgMinResTol      = 1.e-11;    // default tolerance
43 Int_t    AliMillePede2::fgMinResMaxIter  = 10000;     // default max number of iterations 
44 Int_t    AliMillePede2::fgIterSol        = AliMinResSolve::kSolMinRes; // default iterative solver
45 Int_t    AliMillePede2::fgNKrylovV       = 240;        // default number of Krylov vectors to keep
46
47 //_____________________________________________________________________________________________
48 AliMillePede2::AliMillePede2() 
49 : fNLocPar(0),
50   fNGloPar(0),
51   fNGloSize(0),
52 //
53   fNLocEquations(0),
54   fIter(0),
55   fMaxIter(10),
56   fNStdDev(3),
57   fNGloConstraints(0),
58   fNLagrangeConstraints(0),
59   fNLocFits(0),
60   fNLocFitsRejected(0),
61   fNGloFix(0),
62   fGloSolveStatus(kFailed),
63 //
64   fChi2CutFactor(1.),
65   fChi2CutRef(1.),
66   fResCutInit(100.),
67   fResCut(100.),
68   fMinPntValid(1),
69 //
70   fNGroupsSet(0),
71   fParamGrID(0),
72   fProcPnt(0),
73   fVecBLoc(0),
74   fDiagCGlo(0),
75   fVecBGlo(0),
76   fInitPar(0),
77   fDeltaPar(0),
78   fSigmaPar(0),
79   fIsLinear(0),
80   fConstrUsed(0),
81 //
82   fGlo2CGlo(0),
83   fCGlo2Glo(0),
84 //
85   fMatCLoc(0),
86   fMatCGlo(0),
87   fMatCGloLoc(0),
88   //
89   fFillIndex(0),
90   fFillValue(0),
91   //
92   fRecDataTreeName("AliMillePedeRecords_Data"),
93   fRecConsTreeName("AliMillePedeRecords_Consaints"),
94   fRecDataBranchName("Record_Data"),
95   fRecConsBranchName("Record_Consaints"),
96   //
97   fDataRecFName("/tmp/mp2_data_records.root"),
98   fRecord(0),
99   fDataRecFile(0),  
100   fTreeData(0),
101   fRecFileStatus(0),
102   //
103   fConstrRecFName("/tmp/mp2_constraints_records.root"),
104   fTreeConstr(0),
105   fConsRecFile(0),
106   fCurrRecDataID(0),
107   fCurrRecConstrID(0),
108   fLocFitAdd(kTRUE),
109   fSelFirst(1),
110   fSelLast(-1),
111   fRejRunList(0),
112   fAccRunList(0)
113 {}
114
115 //_____________________________________________________________________________________________
116 AliMillePede2::AliMillePede2(const AliMillePede2& src) : 
117   TObject(src),fNLocPar(0),fNGloPar(0),fNGloSize(0),fNLocEquations(0),fIter(0),
118   fMaxIter(10),fNStdDev(3),fNGloConstraints(0),fNLagrangeConstraints(0),
119   fNLocFits(0),fNLocFitsRejected(0),
120   fNGloFix(0),fGloSolveStatus(0),fChi2CutFactor(0),fChi2CutRef(0),fResCutInit(0),
121   fResCut(0),fMinPntValid(1),fNGroupsSet(0),fParamGrID(0),fProcPnt(0),fVecBLoc(0),fDiagCGlo(0),fVecBGlo(0),
122   fInitPar(0),fDeltaPar(0),fSigmaPar(0),fIsLinear(0),fConstrUsed(0),fGlo2CGlo(0),fCGlo2Glo(0),
123   fMatCLoc(0),fMatCGlo(0),fMatCGloLoc(0),fFillIndex(0),fFillValue(0),
124   fRecDataTreeName(0),fRecConsTreeName(0),fRecDataBranchName(0),fRecConsBranchName(0),
125   fDataRecFName(0),fRecord(0),fDataRecFile(0),
126   fTreeData(0),fRecFileStatus(0),fConstrRecFName(0),fTreeConstr(0),fConsRecFile(0),fCurrRecDataID(0),
127   fCurrRecConstrID(0),fLocFitAdd(kTRUE),
128   fSelFirst(1),
129   fSelLast(-1),
130   fRejRunList(0),
131   fAccRunList(0)
132 {printf("Dummy\n");}
133
134 //_____________________________________________________________________________________________
135 AliMillePede2::~AliMillePede2() 
136 {
137   // destructor
138   CloseDataRecStorage();
139   CloseConsRecStorage();
140   //
141   delete[] fParamGrID;
142   delete[] fProcPnt;
143   delete[] fVecBLoc;
144   delete[] fDiagCGlo;
145   delete[] fVecBGlo;
146   delete[] fInitPar;
147   delete[] fDeltaPar;
148   delete[] fSigmaPar;
149   delete[] fGlo2CGlo;
150   delete[] fCGlo2Glo;
151   delete[] fIsLinear;
152   delete[] fConstrUsed;
153   delete[] fFillIndex;
154   delete[] fFillValue;
155   //
156   delete fRecord;
157   delete fMatCLoc;
158   delete fMatCGlo;
159   delete fMatCGloLoc;
160   delete fRejRunList;
161   delete fAccRunList;
162
163
164 //_____________________________________________________________________________________________
165 Int_t AliMillePede2::InitMille(int nGlo, int nLoc, int lNStdDev,double lResCut, double lResCutInit)
166 {
167   // init all
168   if (nLoc>0)        fNLocPar = nLoc;
169   if (nGlo>0)        fNGloPar = nGlo;
170   if (lResCutInit>0) fResCutInit = lResCutInit; 
171   if (lResCut>0)     fResCut     = lResCut;
172   if (lNStdDev>0)    fNStdDev    = lNStdDev;
173   //
174   AliInfo(Form("NLoc: %d NGlo: %d",fNLocPar,fNGloPar));
175
176   fNGloSize = fNGloPar;
177   //
178   if (fgIsMatGloSparse) {fMatCGlo = new AliMatrixSparse(fNGloPar); fMatCGlo->SetSymmetric(kTRUE);}
179   else                   fMatCGlo = new AliSymMatrix(fNGloPar);
180   //
181   fFillIndex    = new Int_t[fNGloPar];
182   fFillValue    = new Double_t[fNGloPar];
183   //
184   fMatCLoc      = new AliSymMatrix(fNLocPar);
185   fMatCGloLoc   = new AliRectMatrix(fNGloPar,fNLocPar);
186   //
187   fParamGrID    = new Int_t[fNGloPar];
188   fProcPnt      = new Int_t[fNGloPar];
189   fVecBLoc      = new Double_t[fNLocPar];
190   fDiagCGlo     = new Double_t[fNGloPar];
191   //
192   fInitPar      = new Double_t[fNGloPar];
193   fDeltaPar     = new Double_t[fNGloPar];
194   fSigmaPar     = new Double_t[fNGloPar];
195   fIsLinear     = new Bool_t[fNGloPar];
196   //
197   fGlo2CGlo     = new Int_t[fNGloPar];
198   fCGlo2Glo     = new Int_t[fNGloPar];
199   //
200   memset(fVecBLoc   ,0,fNLocPar*sizeof(Double_t));
201   memset(fDiagCGlo  ,0,fNGloPar*sizeof(Double_t));
202   memset(fInitPar   ,0,fNGloPar*sizeof(Double_t));
203   memset(fDeltaPar  ,0,fNGloPar*sizeof(Double_t));
204   memset(fSigmaPar  ,0,fNGloPar*sizeof(Double_t));
205   memset(fProcPnt   ,0,fNGloPar*sizeof(Int_t));
206   //
207   for (int i=fNGloPar;i--;) {
208     fGlo2CGlo[i] = fCGlo2Glo[i] = -1;
209     fIsLinear[i] = kTRUE;
210     fParamGrID[i] = -1;
211   }
212   //
213   return 1;
214 }
215
216 //_____________________________________________________________________________________________
217 Bool_t AliMillePede2::ImposeDataRecFile(const char* fname)
218 {
219   // set filename for records
220   CloseDataRecStorage();
221   SetDataRecFName(fname);
222   return InitDataRecStorage(kTRUE); // open in read mode
223 }
224
225 //_____________________________________________________________________________________________
226 Bool_t AliMillePede2::ImposeConsRecFile(const char* fname)
227 {
228   // set filename for constraints
229   CloseConsRecStorage();
230   SetConsRecFName(fname);
231   return InitConsRecStorage(kTRUE); // open in read mode
232 }
233
234 //_____________________________________________________________________________________________
235 Bool_t AliMillePede2::InitDataRecStorage(Bool_t read)
236 {
237   // initialize the buffer for processed measurements records
238   //
239   if (fTreeData) {AliInfo("Data Records File is already initialized"); return kFALSE;} 
240   //
241   if (!fRecord) fRecord = new AliMillePedeRecord();
242   //
243   if (!read) { // write mode: cannot use chain
244     fDataRecFile = TFile::Open(GetDataRecFName(),"recreate");
245     if (!fDataRecFile) {AliFatal(Form("Failed to initialize data records file %s",GetDataRecFName())); return kFALSE;}
246     AliInfo(Form("File %s used for derivatives records",GetDataRecFName()));
247     fTreeData = new TTree(GetRecDataTreeName(),"Data Records for AliMillePede2");
248     fTreeData->Branch(GetRecDataBranchName(),"AliMillePedeRecord",&fRecord,32000,99);
249   }
250   else { // use chain
251     TChain* ch = new TChain(GetRecDataTreeName());
252     //
253     if (fDataRecFName.EndsWith(".root")) ch->AddFile(fDataRecFName);
254     else { // assume text file with list of filenames
255       //
256       ifstream inpf(fDataRecFName.Data());
257       if (!inpf.good()) {AliInfo(Form("Failed on input records list %s\n",fDataRecFName.Data())); return kFALSE;}
258       //
259       TString recfName;
260       while ( !(recfName.ReadLine(inpf)).eof() ) {
261         recfName = recfName.Strip(TString::kBoth,' ');
262         if (recfName.BeginsWith("//") || recfName.BeginsWith("#") || !recfName.EndsWith(".root")) {  // comment
263           AliInfo(Form("Skip %s\n",recfName.Data()));
264           continue;
265         }
266         //
267         recfName = recfName.Strip(TString::kBoth,',');
268         recfName = recfName.Strip(TString::kBoth,'"');
269         gSystem->ExpandPathName(recfName);
270         printf("Adding %s\n",recfName.Data());
271         ch->AddFile(recfName.Data());
272       }      
273     }
274     //
275     Long64_t nent = ch->GetEntries();
276     if (nent<1) { AliInfo("Obtained chain is empty"); return kFALSE;}
277     fTreeData = ch;
278     fTreeData->SetBranchAddress(GetRecDataBranchName(),&fRecord);
279     AliInfo(Form("Found %lld derivatives records",nent));
280   }
281   fCurrRecDataID = -1;
282   fRecFileStatus = read ? 1:2;
283   //
284   return kTRUE;
285 }
286
287 //_____________________________________________________________________________________________
288 Bool_t AliMillePede2::InitConsRecStorage(Bool_t read)
289 {
290   // initialize the buffer for processed measurements records
291   //
292   if (fConsRecFile) {AliInfo("Constraints Records File is already initialized"); return kFALSE;} 
293   //
294   if (!fRecord) fRecord = new AliMillePedeRecord();
295   //
296   fConsRecFile = TFile::Open(GetConsRecFName(),read ? "":"recreate");
297   if (!fConsRecFile) {AliInfo(Form("Failed to initialize constraints records file %s",GetConsRecFName())); return kFALSE;}
298   //
299   AliInfo(Form("File %s used for constraints records",GetConsRecFName()));
300   if (read) {
301     fTreeConstr = (TTree*)fConsRecFile->Get(GetRecConsTreeName());
302     if (!fTreeConstr) {AliInfo(Form("Did not find constraints records tree in %s",GetConsRecFName())); return kFALSE;}
303     fTreeConstr->SetBranchAddress(GetRecConsBranchName(),&fRecord);
304     AliInfo(Form("Found %lld constraints records",fTreeConstr->GetEntries()));
305     //
306   }
307   else {
308     //
309     fTreeConstr = new TTree(GetRecConsTreeName(),"Constraints Records for AliMillePede2");
310     fTreeConstr->Branch(GetRecConsBranchName(),"AliMillePedeRecord",&fRecord,32000,99);
311   }
312   fCurrRecConstrID = -1;
313   //
314   return kTRUE;
315 }
316
317 //_____________________________________________________________________________________________
318 void AliMillePede2::CloseDataRecStorage()
319 {
320   // close records file
321   if (fTreeData) {
322     if (fDataRecFile && fDataRecFile->IsWritable()) {
323       fDataRecFile->cd();
324       fTreeData->Write();
325     }
326     delete fTreeData;  
327     fTreeData = 0;
328     if (fDataRecFile) {
329       fDataRecFile->Close();
330       delete fDataRecFile;
331       fDataRecFile = 0;
332     }
333   }
334   fRecFileStatus = 0;
335   //
336 }
337
338 //_____________________________________________________________________________________________
339 void AliMillePede2::CloseConsRecStorage()
340 {
341   // close constraints file
342   if (fTreeConstr) {
343     if (fConsRecFile->IsWritable()) {
344       fConsRecFile->cd();
345       fTreeConstr->Write();
346     }
347     delete fTreeConstr;  
348     fTreeConstr = 0;
349     fConsRecFile->Close();
350     delete fConsRecFile;
351     fConsRecFile = 0;
352   }
353   //
354 }
355
356 //_____________________________________________________________________________________________
357 Bool_t AliMillePede2::ReadNextRecordData()
358 {
359   // read next data record (if any)
360   if (!fTreeData || ++fCurrRecDataID >= fTreeData->GetEntries()) { fCurrRecDataID--; return kFALSE;}
361   fTreeData->GetEntry(fCurrRecDataID);
362   return kTRUE;
363 }
364
365 //_____________________________________________________________________________________________
366 Bool_t AliMillePede2::ReadNextRecordConstraint()
367 {
368   // read next constraint record (if any)
369   if (!fTreeConstr || ++fCurrRecConstrID >= fTreeConstr->GetEntries()) { fCurrRecConstrID--; return kFALSE;}
370   fTreeConstr->GetEntry(fCurrRecConstrID);
371   return kTRUE;
372 }
373
374 //_____________________________________________________________________________________________
375 void AliMillePede2::SetRecordWeight(double wgh)
376 {
377   // assign weight
378   if (fRecFileStatus<2) InitDataRecStorage(); // create a buffer to store the data
379   fRecord->SetWeight(wgh);
380 }
381
382 //_____________________________________________________________________________________________
383 void AliMillePede2::SetRecordRun(Int_t run)
384 {
385   // assign run
386   if (fRecFileStatus<2) InitDataRecStorage(); // create a buffer to store the data
387   fRecord->SetRunID(run);
388 }
389
390 //_____________________________________________________________________________________________
391 void AliMillePede2::SetLocalEquation(double *dergb, double *derlc, double lMeas, double lSigma)
392 {
393   // assing derivs of loc.eq.
394   if (fRecFileStatus<2) InitDataRecStorage(); // create a buffer to store the data
395   //
396   // write data of single measurement
397   if (lSigma<=0.0) { // If parameter is fixed, then no equation
398     for (int i=fNLocPar; i--;) derlc[i] = 0.0;
399     for (int i=fNGloPar; i--;) dergb[i] = 0.0;
400     return;
401   }
402   //
403   fRecord->AddResidual(lMeas);
404   //
405   // Retrieve local param interesting indices
406   for (int i=0;i<fNLocPar;i++) if (!IsZero(derlc[i])) {fRecord->AddIndexValue(i,derlc[i]); derlc[i] = 0.0;}
407   //
408   fRecord->AddWeight( 1.0/lSigma/lSigma );
409   //
410   // Idem for global parameters
411   for (int i=0;i<fNGloPar;i++) if (!IsZero(dergb[i])) {
412     fRecord->AddIndexValue(i,dergb[i]); dergb[i] = 0.0;
413     fRecord->MarkGroup(fParamGrID[i]);
414   }
415   //  fRecord->Print();
416   //
417 }
418
419 //_____________________________________________________________________________________________
420 void AliMillePede2::SetLocalEquation(int *indgb, double *dergb, int ngb, int *indlc,
421                                      double *derlc,int nlc,double lMeas,double lSigma)
422 {       
423   // write data of single measurement
424   if (lSigma<=0.0) { // If parameter is fixed, then no equation
425     for (int i=nlc;i--;)  derlc[i] = 0.0;
426     for (int i=ngb;i--;)  dergb[i] = 0.0;
427     return;
428   }
429   //
430   if (fRecFileStatus<2) InitDataRecStorage(); // create a buffer to store the data
431   //
432   fRecord->AddResidual(lMeas);
433   //
434   // Retrieve local param interesting indices
435   for (int i=0;i<nlc;i++) if (!IsZero(derlc[i])) {fRecord->AddIndexValue(indlc[i],derlc[i]); derlc[i]=0.; indlc[i]=0;}
436   //
437   fRecord->AddWeight( 1./lSigma/lSigma );
438   //
439   // Idem for global parameters
440   for (int i=0;i<ngb;i++) if (!IsZero(dergb[i])) {fRecord->AddIndexValue(indgb[i],dergb[i]); dergb[i]=0.; indgb[i]=0;} 
441   //
442 }
443
444
445 //_____________________________________________________________________________________________
446 void AliMillePede2::SetGlobalConstraint(const double *dergb, double val, double sigma)
447 {       
448   // Define a constraint equation.
449   if (!fConsRecFile || !fConsRecFile->IsWritable()) InitConsRecStorage(); // create a buffer to store the data
450   //
451   fRecord->Reset();
452   fRecord->AddResidual(val);
453   fRecord->AddWeight(sigma);
454   for (int i=0; i<fNGloPar; i++) if (!IsZero(dergb[i]))  fRecord->AddIndexValue(i,dergb[i]);
455   fNGloConstraints++;
456   if (IsZero(sigma)) fNLagrangeConstraints++;
457   //  printf("NewConstraint:\n"); fRecord->Print(); //RRR
458   SaveRecordConstraint();
459   //
460 }
461
462 //_____________________________________________________________________________________________
463 void AliMillePede2::SetGlobalConstraint(const int *indgb, const double *dergb, int ngb, double val,double sigma)
464 {       
465   // Define a constraint equation.
466   if (!fConsRecFile || !fConsRecFile->IsWritable()) InitConsRecStorage(); // create a buffer to store the data
467   fRecord->Reset();
468   fRecord->AddResidual(val);
469   fRecord->AddWeight(sigma);   // dummy
470   for (int i=0; i<ngb; i++) if (!IsZero(dergb[i]))  fRecord->AddIndexValue(indgb[i],dergb[i]);
471   fNGloConstraints++;
472   if (IsZero(sigma)) fNLagrangeConstraints++;
473   SaveRecordConstraint();
474   //
475 }
476
477 //_____________________________________________________________________________________________
478 Int_t AliMillePede2::LocalFit(double *localParams)
479 {
480   /*
481     Perform local parameters fit once all the local equations have been set
482     -----------------------------------------------------------
483     localParams = (if !=0) will contain the fitted track parameters and
484     related errors
485   */
486   static int     nrefSize = 0;
487   //  static TArrayI refLoc,refGlo,nrefLoc,nrefGlo;
488   static Int_t  *refLoc=0,*refGlo=0,*nrefLoc=0,*nrefGlo=0;
489   int nPoints = 0;
490   //
491   AliSymMatrix    &matCLoc    = *fMatCLoc;
492   AliMatrixSq     &matCGlo    = *fMatCGlo;
493   AliRectMatrix   &matCGloLoc = *fMatCGloLoc;
494   //
495   memset(fVecBLoc,0,fNLocPar*sizeof(double));
496   matCLoc.Reset();      
497   //
498   int cnt=0;
499   int recSz = fRecord->GetSize();
500   //
501   while(cnt<recSz) {  // Transfer the measurement records to matrices
502     //
503     // extract addresses of residual, weight and pointers on local and global derivatives for each point
504     if (nrefSize<=nPoints) {
505       int *tmpA = 0;
506       nrefSize = 2*(nPoints+1); 
507       tmpA = refLoc;  refLoc  = new Int_t[nrefSize]; if (tmpA) memcpy(refLoc,tmpA,nPoints*sizeof(int));
508       tmpA = refGlo;  refGlo  = new Int_t[nrefSize]; if (tmpA) memcpy(refGlo,tmpA,nPoints*sizeof(int));
509       tmpA = nrefLoc; nrefLoc = new Int_t[nrefSize]; if (tmpA) memcpy(nrefLoc,tmpA,nPoints*sizeof(int));
510       tmpA = nrefGlo; nrefGlo = new Int_t[nrefSize]; if (tmpA) memcpy(nrefGlo,tmpA,nPoints*sizeof(int));
511     }
512     //
513     refLoc[nPoints]  = ++cnt;
514     int nLoc = 0;
515     while(!fRecord->IsWeight(cnt)) {nLoc++; cnt++;}
516     nrefLoc[nPoints] = nLoc;
517     //
518     refGlo[nPoints]  = ++cnt;
519     int nGlo = 0;
520     while(!fRecord->IsResidual(cnt) && cnt<recSz) {nGlo++; cnt++;} 
521     nrefGlo[nPoints] = nGlo;
522     //
523     nPoints++;
524   }
525   double vl;
526   //
527   double gloWgh = fRecord->GetWeight(); // global weight for this set
528   Int_t maxLocUsed = 0;
529   //
530   for (int ip=nPoints;ip--;) {  // Transfer the measurement records to matrices
531     double  resid  = fRecord->GetValue( refLoc[ip]-1 );
532     double  weight = fRecord->GetValue( refGlo[ip]-1 )*gloWgh;    
533     double *derLoc = fRecord->GetValue()+refLoc[ip];
534     double *derGlo = fRecord->GetValue()+refGlo[ip];
535     int    *indLoc = fRecord->GetIndex()+refLoc[ip];
536     int    *indGlo = fRecord->GetIndex()+refGlo[ip];
537     //
538     for (int i=nrefGlo[ip];i--;) {      // suppress the global part (only relevant with iterations)
539       int iID = indGlo[i];              // Global param indice
540       if (fSigmaPar[iID]<=0.) continue;                                    // fixed parameter RRRCheck
541       if (fIsLinear[iID]) resid -= derGlo[i]*(fInitPar[iID]+fDeltaPar[iID]);  // linear parameter
542       else                resid -= derGlo[i]*fDeltaPar[iID];                  // nonlinear parameter
543     }
544     //
545     // Symmetric matrix, don't bother j>i coeffs
546     for (int i=nrefLoc[ip];i--;) {         // Fill local matrix and vector
547       fVecBLoc[ indLoc[i] ] += weight*resid*derLoc[i];
548       if (indLoc[i]>maxLocUsed) maxLocUsed = indLoc[i];  
549       for (int j=i+1;j--;) matCLoc(indLoc[i] ,indLoc[j]) += weight*derLoc[i]*derLoc[j];
550     }
551     //
552   } // end of the transfer of the measurement record to matrices
553   //
554   matCLoc.SetSizeUsed(++maxLocUsed);   // data with B=0 may use less than declared nLocals 
555   //
556   /* //RRR
557   fRecord->Print("l");
558   printf("\nBefore\nLocalMatrix: "); matCLoc.Print("l");
559   printf("RHSLoc: "); for (int i=0;i<fNLocPar;i++) printf("%+e |",fVecBLoc[i]); printf("\n");
560   */
561   // first try to solve by faster Cholesky decomposition, then by Gaussian elimination
562   if (!matCLoc.SolveChol(fVecBLoc,kTRUE)) {
563     AliInfo("Failed to solve locals by Cholesky, trying Gaussian Elimination");
564     if (!matCLoc.SolveSpmInv(fVecBLoc,kTRUE)) { 
565       AliInfo("Failed to solve locals by Gaussian Elimination, skip...");
566       matCLoc.Print("d");
567       return 0; // failed to solve
568     }
569   }
570   //
571   // If requested, store the track params and errors
572   //RRR  printf("locfit: "); for (int i=0;i<fNLocPar;i++) printf("%+e |",fVecBLoc[i]); printf("\n");
573   
574   if (localParams) for (int i=maxLocUsed; i--;) {
575       localParams[2*i]   = fVecBLoc[i];
576       localParams[2*i+1] = TMath::Sqrt(TMath::Abs(matCLoc.QueryDiag(i)));
577     }
578   //
579   float lChi2 = 0;
580   int   nEq   = 0;
581   //
582   for (int ip=nPoints;ip--;) {  // Calculate residuals
583     double  resid  = fRecord->GetValue( refLoc[ip]-1 );
584     double  weight = fRecord->GetValue( refGlo[ip]-1 )*gloWgh;    
585     double *derLoc = fRecord->GetValue()+refLoc[ip];
586     double *derGlo = fRecord->GetValue()+refGlo[ip];
587     int    *indLoc = fRecord->GetIndex()+refLoc[ip];
588     int    *indGlo = fRecord->GetIndex()+refGlo[ip];
589     //
590     // Suppress local and global contribution in residuals;
591     for (int i=nrefLoc[ip];i--;) resid -= derLoc[i]*fVecBLoc[ indLoc[i] ];     // local part 
592     //
593     for (int i=nrefGlo[ip];i--;) {                                             // global part
594       int iID = indGlo[i];
595       if ( fSigmaPar[iID] <= 0.) continue;                                     // fixed parameter RRRCheck
596       if (fIsLinear[iID]) resid -= derGlo[i]*(fInitPar[iID]+fDeltaPar[iID]);   // linear parameter
597       else                resid -= derGlo[i]*fDeltaPar[iID];                   // nonlinear parameter
598     }
599     //
600     // reject the track if the residual is too large (outlier)
601     double absres = TMath::Abs(resid);
602     if ( (absres >= fResCutInit && fIter ==1 ) ||
603          (absres >= fResCut     && fIter > 1)) {
604       if (fLocFitAdd) fNLocFitsRejected++;      
605       //      printf("reject res %5ld %+e\n",fCurrRecDataID,resid);
606       return 0;
607     }
608     // 
609     lChi2 += weight*resid*resid ; // total chi^2
610     nEq++;                        // number of equations                        
611   } // end of Calculate residuals
612   //
613   lChi2 /= gloWgh;  
614   int nDoF = nEq-maxLocUsed;
615   lChi2 = (nDoF>0) ? lChi2/nDoF : 0;  // Chi^2/dof  
616   //
617   if (fNStdDev != 0 && nDoF>0 && lChi2 > Chi2DoFLim(fNStdDev,nDoF)*fChi2CutFactor) { // check final chi2
618     if (fLocFitAdd) fNLocFitsRejected++;      
619     //    printf("reject chi2 %5ld: %+e\n",fCurrRecDataID, lChi2);
620     return 0;
621   }
622   //
623   if (fLocFitAdd) {
624     fNLocFits++;
625     fNLocEquations += nEq;
626   }
627   else {
628     fNLocFits--;
629     fNLocEquations -= nEq;
630   }
631   //
632   //  local operations are finished, track is accepted 
633   //  We now update the global parameters (other matrices)
634   //
635   int nGloInFit = 0;
636   //
637   for (int ip=nPoints;ip--;) {  // Update matrices
638     double  resid  = fRecord->GetValue( refLoc[ip]-1 );
639     double  weight = fRecord->GetValue( refGlo[ip]-1 )*gloWgh;    
640     double *derLoc = fRecord->GetValue()+refLoc[ip];
641     double *derGlo = fRecord->GetValue()+refGlo[ip];
642     int    *indLoc = fRecord->GetIndex()+refLoc[ip];
643     int    *indGlo = fRecord->GetIndex()+refGlo[ip];
644     //
645     for (int i=nrefGlo[ip];i--;) {    // suppress the global part
646       int iID = indGlo[i];           // Global param indice
647       if ( fSigmaPar[iID] <= 0.) continue;                                         // fixed parameter RRRCheck
648       if (fIsLinear[iID])  resid -= derGlo[i]*(fInitPar[iID]+fDeltaPar[iID]);      // linear parameter
649       else                 resid -= derGlo[i]*fDeltaPar[iID];                      // nonlinear parameter
650     }
651     //
652     for (int ig=nrefGlo[ip];ig--;) {
653       int iIDg = indGlo[ig];   // Global param indice (the matrix line)          
654       if ( fSigmaPar[iIDg] <= 0.) continue;                                    // fixed parameter RRRCheck
655       if (fLocFitAdd) fVecBGlo[ iIDg ] += weight*resid*derGlo[ig]; //!!!
656       else            fVecBGlo[ iIDg ] -= weight*resid*derGlo[ig]; //!!!
657       //
658       // First of all, the global/global terms (exactly like local matrix)
659       int nfill = 0;
660       for (int jg=ig+1;jg--;) { // matCGlo is symmetric by construction  
661         int jIDg = indGlo[jg];
662         if ( fSigmaPar[jIDg] <= 0.) continue;                                    // fixed parameter RRRCheck
663         if ( !IsZero(vl = weight*derGlo[ig]*derGlo[jg]) ) {
664           fFillIndex[nfill]   = jIDg;
665           fFillValue[nfill++] = fLocFitAdd ? vl:-vl;
666         }
667       }
668       if (nfill) matCGlo.AddToRow(iIDg,fFillValue,fFillIndex,nfill);
669       //
670       // Now we have also rectangular matrices containing global/local terms.
671       int iCIDg = fGlo2CGlo[iIDg];  // compressed Index of index          
672       if (iCIDg == -1) {
673         Double_t *rowGL = matCGloLoc(nGloInFit);
674         for (int k=maxLocUsed;k--;) rowGL[k] = 0.0;  // reset the row
675         iCIDg = fGlo2CGlo[iIDg] = nGloInFit;
676         fCGlo2Glo[nGloInFit++] = iIDg;
677       }
678       //
679       Double_t *rowGLIDg = matCGloLoc(iCIDg);
680       for (int il=nrefLoc[ip];il--;) rowGLIDg[ indLoc[il] ] += weight*derGlo[ig]*derLoc[il];
681       fProcPnt[iIDg] += fLocFitAdd ? 1:-1;       // update counter
682       //
683     }
684   } // end of Update matrices
685   //
686   /*//RRR
687   printf("After GLO\n");
688   printf("MatCLoc: "); fMatCLoc->Print("l");
689   printf("MatCGlo: "); fMatCGlo->Print("l");
690   printf("MatCGlLc:"); fMatCGloLoc->Print("l");
691   printf("BGlo: "); for (int i=0; i<fNGloPar; i++) printf("%+e |",fVecBGlo[i]); printf("\n");
692   */
693   // calculate fMatCGlo -= fMatCGloLoc * fMatCLoc * fMatCGloLoc^T
694   // and       fVecBGlo -= fMatCGloLoc * fVecBLoc
695   //
696   //-------------------------------------------------------------- >>>
697   double vll;
698   for (int iCIDg=0; iCIDg<nGloInFit; iCIDg++) {
699     int iIDg = fCGlo2Glo[iCIDg];
700     //
701     vl = 0;
702     Double_t *rowGLIDg =  matCGloLoc(iCIDg);
703     for (int kl=0;kl<maxLocUsed;kl++) if (rowGLIDg[kl]) vl += rowGLIDg[kl]*fVecBLoc[kl];
704     if  (!IsZero(vl)) fVecBGlo[iIDg] -= fLocFitAdd ? vl : -vl;
705     //
706     int nfill = 0;
707     for (int jCIDg=0;jCIDg<=iCIDg; jCIDg++) {  
708       int jIDg = fCGlo2Glo[jCIDg];
709       //
710       vl = 0;
711       Double_t *rowGLJDg =  matCGloLoc(jCIDg);
712       for (int kl=0;kl<maxLocUsed;kl++) {
713         // diag terms
714         if ( (!IsZero(vll=rowGLIDg[kl]*rowGLJDg[kl])) ) vl += matCLoc.QueryDiag(kl)*vll;
715         //
716         // off-diag terms
717         for (int ll=0;ll<kl;ll++) {
718           if ( !IsZero(vll=rowGLIDg[kl]*rowGLJDg[ll]) ) vl += matCLoc(kl,ll)*vll;
719           if ( !IsZero(vll=rowGLIDg[ll]*rowGLJDg[kl]) ) vl += matCLoc(kl,ll)*vll;
720         }
721       }
722       if (!IsZero(vl)) {
723         fFillIndex[nfill]   = jIDg;
724         fFillValue[nfill++] = fLocFitAdd ? -vl : vl;
725       }
726     }
727     if (nfill)  matCGlo.AddToRow(iIDg,fFillValue,fFillIndex,nfill);
728   }
729   //
730   // reset compressed index array
731   //
732   /*//RRR
733   printf("After GLOLoc\n");
734   printf("MatCGlo: "); fMatCGlo->Print("");
735   printf("BGlo: "); for (int i=0; i<fNGloPar; i++) printf("%+e |",fVecBGlo[i]); printf("\n");
736   */
737   for (int i=nGloInFit;i--;) { 
738     fGlo2CGlo[ fCGlo2Glo[i] ] = -1;
739     fCGlo2Glo[i] = -1;
740   }
741   //
742   //---------------------------------------------------- <<<
743   return 1;
744 }
745
746 //_____________________________________________________________________________________________
747 Int_t AliMillePede2::GlobalFit(Double_t *par, Double_t *error, Double_t *pull)
748 {
749   // performs a requested number of global iterations
750   fIter = 1;
751   //
752   TStopwatch sw; sw.Start();
753   //
754   int res = 0;
755   AliInfo("Starting Global fit.");
756   while (fIter<=fMaxIter) {
757     //
758     res = GlobalFitIteration();
759     if (!res) break;
760     //
761     if (!IsZero(fChi2CutFactor-fChi2CutRef)) {    
762       fChi2CutFactor = TMath::Sqrt(fChi2CutFactor);
763       if (fChi2CutFactor < 1.2*fChi2CutRef) {
764         fChi2CutFactor = fChi2CutRef;
765         //RRR   fIter = fMaxIter - 1;     // Last iteration
766       }
767     }
768     fIter++;
769   }
770   //
771   sw.Stop();
772   AliInfo(Form("Global fit %s, CPU time: %.1f",res ? "Converged":"Failed",sw.CpuTime()));  
773   if (!res) return 0;
774   //
775   if (par) for (int i=fNGloPar;i--;) par[i] = fInitPar[i]+fDeltaPar[i]; 
776   //
777   if (fGloSolveStatus==kInvert) { // errors on params are available
778     if (error) for (int i=fNGloPar;i--;) error[i] = fProcPnt[i]>0 ? TMath::Sqrt(TMath::Abs(fMatCGlo->QueryDiag(i))) : 0.;
779     if (pull)  for (int i=fNGloPar;i--;) pull[i]  = fProcPnt[i]>0 && (fSigmaPar[i]*fSigmaPar[i]-fMatCGlo->QueryDiag(i))>0. && fSigmaPar[i]>0 
780                                            ? fDeltaPar[i]/TMath::Sqrt(fSigmaPar[i]*fSigmaPar[i]-fMatCGlo->QueryDiag(i)) : 0;
781   }
782   //
783   return 1;
784 }
785
786 //_____________________________________________________________________________________________
787 Int_t AliMillePede2::GlobalFitIteration()
788 {
789   // perform global parameters fit once all the local equations have been fitted
790   //
791   AliInfo(Form("Global Fit Iteration#%2d (Local Fit Chi^2 cut factor: %.2f)",fIter,fChi2CutFactor));
792   //
793   if (!fNGloPar || !fTreeData) {
794     AliInfo("No data was stored, stopping iteration");
795     return 0;
796   }
797   TStopwatch sw,sws; 
798   sw.Start();
799   sws.Stop();
800   //
801   if (!fConstrUsed) {
802     fConstrUsed = new Bool_t[fNGloConstraints];
803     memset(fConstrUsed,0,fNGloConstraints*sizeof(Bool_t));
804   }
805   // Reset all info specific for this step
806   AliMatrixSq& matCGlo = *fMatCGlo;
807   matCGlo.Reset();
808   memset(fProcPnt,0,fNGloPar*sizeof(Int_t));
809   //
810   fNGloConstraints = fTreeConstr ? fTreeConstr->GetEntries() : 0;
811   //
812   // count number of Lagrange constraints: they need new row/cols to be added
813   fNLagrangeConstraints = 0;
814   for (int i=0; i<fNGloConstraints; i++) {
815     ReadRecordConstraint(i);
816     if ( IsZero(fRecord->GetValue(1)) ) fNLagrangeConstraints++; // exact constraint (no error) -> Lagrange multiplier 
817   }
818   //
819   // if needed, readjust the size of the global vector (for matrices this is done automatically)
820   if (!fVecBGlo || fNGloSize!=fNGloPar+fNLagrangeConstraints) {
821     delete[] fVecBGlo;   // in case some constraint was added between the two manual iterations
822     fNGloSize = fNGloPar+fNLagrangeConstraints;
823     fVecBGlo = new Double_t[fNGloSize];
824   }
825   memset(fVecBGlo,0,fNGloSize*sizeof(double));
826   //
827   fNLocFits         = 0;
828   fNLocFitsRejected = 0;
829   fNLocEquations    = 0;
830   //
831   //  Process data records and build the matrices
832   Long_t ndr = fTreeData->GetEntries();
833   Long_t first = fSelFirst>0  ? fSelFirst : 0;
834   Long_t last  = fSelLast<1   ? ndr : (fSelLast>=ndr ? ndr : fSelLast+Long_t(1));
835   ndr = last - first;
836   //
837   AliInfo(Form("Building the Global matrix from data records %ld : %ld",first,last));
838   if (ndr<1) return 0;
839   //
840   TStopwatch swt; swt.Start();
841   fLocFitAdd = kTRUE;  // add contributions of matching tracks
842   for (Long_t i=0;i<ndr;i++) {
843     Long_t iev = i+first;
844     ReadRecordData(iev);
845     if (!IsRecordAcceptable()) continue;
846     LocalFit();
847     if ( (i%int(0.2*ndr)) == 0) printf("%.1f%% of local fits done\n", double(100.*i)/ndr);
848   }
849   swt.Stop();
850   printf("%ld local fits done: ", ndr);
851   /*
852   printf("MatCGlo: "); fMatCGlo->Print("l");
853   printf("BGlo: "); for (int i=0; i<fNGloPar; i++) printf("%+e |",fVecBGlo[i]); printf("\n");
854   swt.Print();
855   */
856   sw.Start(kFALSE);
857   //
858   //
859   // ---------------------- Reject parameters with low statistics ------------>>
860   fNGloFix = 0;
861   if (fMinPntValid>1 && fNGroupsSet) {
862     //
863     printf("Checking parameters with statistics < %d\n",fMinPntValid);
864     TStopwatch swsup;
865     swsup.Start();
866     // 1) build the list of parameters to fix
867     Int_t fixArrSize = 10;
868     Int_t nFixedGroups = 0;
869     TArrayI fixGroups(fixArrSize);
870     //
871     int grIDold = -2;
872     int oldStart = -1;
873     double oldMin = 1.e20;
874     double oldMax =-1.e20;
875     //
876     for (int i=fNGloPar;i--;) { // // Reset row and column of fixed params and add 1/sig^2 to free ones
877       int grID = fParamGrID[i];
878       if (grID<0) continue; // not in the group
879       //
880       if (grID!=grIDold) { // starting new group
881         if (grIDold>=0) { // decide if the group has enough statistics
882           if (oldMin<fMinPntValid && oldMax<2*fMinPntValid) { // suppress group
883             for (int iold=oldStart;iold>i;iold--) fProcPnt[iold] = 0;
884             Bool_t fnd = kFALSE;    // check if the group is already accounted
885             for (int j=nFixedGroups;j--;) if (fixGroups[j]==grIDold) {fnd=kTRUE; break;}
886             if (!fnd) {
887               if (nFixedGroups>=fixArrSize) {fixArrSize*=2; fixGroups.Set(fixArrSize);}
888               fixGroups[nFixedGroups++] = grIDold; // add group to fix
889             }
890           }       
891         }
892         grIDold = grID; // mark the start of the new group
893         oldStart = i;
894         oldMin =  1.e20;
895         oldMax = -1.e20;
896       }
897       if (oldMin>fProcPnt[i]) oldMin = fProcPnt[i];
898       if (oldMax<fProcPnt[i]) oldMax = fProcPnt[i];
899       //
900     }
901     // extra check for the last group
902     if (grIDold>=0 && oldMin<fMinPntValid && oldMax<2*fMinPntValid) { // suppress group
903       for (int iold=oldStart;iold--;) fProcPnt[iold] = 0;
904       Bool_t fnd = kFALSE;    // check if the group is already accounted
905       for (int j=nFixedGroups;j--;) if (fixGroups[j]==grIDold) {fnd=kTRUE; break;}
906       if (!fnd) {
907         if (nFixedGroups>=fixArrSize) {fixArrSize*=2; fixGroups.Set(fixArrSize);}
908         fixGroups[nFixedGroups++] = grIDold; // add group to fix
909       }
910     }
911     //
912     // 2) loop over records and add contributions of fixed groups with negative sign
913     fLocFitAdd = kFALSE;
914     //
915     for (Long_t i=0;i<ndr;i++) {
916       Long_t iev = i+first;
917       ReadRecordData(iev);
918       if (!IsRecordAcceptable()) continue;
919       Bool_t suppr = kFALSE;
920       for (int ifx=nFixedGroups;ifx--;)if (fRecord->IsGroupPresent(fixGroups[ifx])) suppr = kTRUE;
921       if (suppr) LocalFit();
922     }
923     fLocFitAdd = kTRUE;
924     //
925     if (nFixedGroups) {
926       printf("Suppressed contributions of groups with NPoints<%d :\n",fMinPntValid);
927       for (int i=0;i<nFixedGroups;i++) printf("%d ",fixGroups[i]); printf("\n");
928     }
929     swsup.Stop();
930     swsup.Print();
931   }
932   // ---------------------- Reject parameters with low statistics ------------<<
933   //
934   // add large number to diagonal of fixed params  
935   //
936   for (int i=fNGloPar;i--;) { // // Reset row and column of fixed params and add 1/sig^2 to free ones
937     //    printf("#%3d : Nproc : %5d   grp: %d\n",i,fProcPnt[i],fParamGrID[i]);
938     if (fProcPnt[i]<1) {
939       fNGloFix++; 
940       fVecBGlo[i] = 0.;
941       matCGlo.DiagElem(i) = 1.;//float(fNLocEquations*fNLocEquations);
942       //      matCGlo.DiagElem(i) = float(fNLocEquations*fNLocEquations);
943     }
944     else matCGlo.DiagElem(i) += (fgWeightSigma ? fProcPnt[i] : 1.)/(fSigmaPar[i]*fSigmaPar[i]);
945   }
946   //
947   for (int i=fNGloPar;i--;) fDiagCGlo[i] = matCGlo.QueryDiag(i); // save the diagonal elements  
948   //
949   // add constraint equations
950   int nVar = fNGloPar;                    // Current size of global matrix      
951   for (int i=0; i<fNGloConstraints; i++) {
952     ReadRecordConstraint(i);
953     double val   = fRecord->GetValue(0);  
954     double sig   = fRecord->GetValue(1);  
955     int    *indV = fRecord->GetIndex()+2;
956     double *der  = fRecord->GetValue()+2;
957     int    csize = fRecord->GetSize()-2;
958     //
959     // check if after suppression of fixed variables there are non-0 derivatives
960     // and determine the max statistics of involved params
961     int nSuppressed = 0;
962     int maxStat = 1;
963     for (int j=csize;j--;) {
964       if (fProcPnt[indV[j]]<1) nSuppressed++; 
965       else {
966         maxStat = TMath::Max(maxStat,fProcPnt[indV[j]]);
967       }
968     }
969     //
970     if (nSuppressed==csize) {
971       //      AliInfo(Form("Neglecting constraint %d of %d derivatives since no free parameters left",i,csize));
972       //
973       // was this constraint ever created ? 
974       if ( sig==0 && fConstrUsed[i] ) { // this is needed only for constraints with Lagrange multiplier
975         // to avoid empty row impose dummy constraint on "Lagrange multiplier"
976         matCGlo.DiagElem(nVar) = 1.;
977         fVecBGlo[nVar++] = 0;
978       }
979       continue;
980     }
981     //
982     // account for already accumulated corrections
983     for (int j=csize; j--;) val -= der[j]*(fInitPar[ indV[j] ]+fDeltaPar[ indV[j] ]);
984     //
985     if (sig > 0) {  // this is a gaussian constriant: no Lagrange multipliers are added
986       //
987       double sig2i = (fgWeightSigma ? TMath::Sqrt(maxStat) : 1.)/sig/sig;
988       for (int ir=0;ir<csize;ir++) {
989         int iID = indV[ir];
990         for (int ic=0;ic<=ir;ic++) { // matrix is symmetric
991           int jID = indV[ic];
992           double vl = der[ir]*der[ic]*sig2i;
993           if (!IsZero(vl)) matCGlo(iID,jID) += vl;
994         }
995         fVecBGlo[iID] += val*der[ir]*sig2i;
996       }
997     }
998     else {   // this is exact constriant:  Lagrange multipliers must be added
999       for (int j=csize; j--;) {
1000         int jID = indV[j];
1001         if (fProcPnt[jID]<1) continue;                      // this parameter was fixed, don't put it into constraint
1002         matCGlo(nVar,jID) = float(fNLocEquations)*der[j];   // fMatCGlo is symmetric, only lower triangle is filled  
1003       }
1004       //
1005       if (matCGlo.QueryDiag(nVar)) matCGlo.DiagElem(nVar) = 0.0;
1006       fVecBGlo[nVar++] = float(fNLocEquations)*val; //RS ? should we use here fNLocFits ? 
1007       fConstrUsed[i] = kTRUE;
1008     }
1009   }
1010   //
1011   AliInfo(Form("Obtained %-7ld equations from %-7ld records (%-7ld rejected). Fixed %-4d globals",
1012                fNLocEquations,fNLocFits,fNLocFitsRejected,fNGloFix));
1013
1014   //
1015   sws.Start();
1016
1017   /*
1018   printf("Solving:\n");
1019   matCGlo.Print();
1020   for (int i=0;i<fNGloSize;i++) printf("b%2d : %+e\n",i,fVecBGlo[i]);
1021   */
1022   fGloSolveStatus = SolveGlobalMatEq();                     // obtain solution for this step
1023   sws.Stop();
1024   printf("Solve %d |",fIter); sws.Print();
1025   //
1026   sw.Stop();
1027   AliInfo(Form("Iteration#%2d %s. CPU time: %.1f",fIter,fGloSolveStatus==kFailed ? "Failed":"Converged",sw.CpuTime()));
1028   if (fGloSolveStatus==kFailed) return 0;
1029   //
1030   for (int i=fNGloPar;i--;) fDeltaPar[i] += fVecBGlo[i];    // Update global parameters values (for iterations)
1031   //
1032   //  PrintGlobalParameters();
1033   return 1;
1034 }
1035
1036 //_____________________________________________________________________________________________
1037 Int_t AliMillePede2::SolveGlobalMatEq()
1038 {
1039   //
1040   // solve global matrix equation MatCGlob*X=VecBGlo and store the result in the VecBGlo
1041   //
1042   /*
1043   printf("GlobalMatrix\n");
1044   fMatCGlo->Print();
1045   printf("RHS\n");
1046   for (int i=0;i<fNGloPar;i++) printf("%d %+e\n",i,fVecBGlo[i]);
1047   */
1048   //
1049   if (!fgIsMatGloSparse) {
1050     //
1051     if (fNLagrangeConstraints==0) { // pos-def systems are faster to solve by Cholesky
1052       if ( ((AliSymMatrix*)fMatCGlo)->SolveChol(fVecBGlo, fgInvChol) ) return fgInvChol ? kInvert:kNoInversion;
1053       else AliInfo("Solution of Global Dense System by Cholesky failed, trying Gaussian Elimiation");
1054     }
1055     //
1056     if (((AliSymMatrix*)fMatCGlo)->SolveSpmInv(fVecBGlo, kTRUE)) return kInvert;
1057     else AliInfo("Solution of Global Dense System by Gaussian Elimination failed, trying iterative methods");
1058   }
1059   // try to solve by minres
1060   TVectorD sol(fNGloSize);
1061   //
1062   AliMinResSolve *slv = new AliMinResSolve(fMatCGlo,fVecBGlo);
1063   if (!slv) return kFailed;
1064   //
1065   Bool_t res = kFALSE;
1066   if      (fgIterSol == AliMinResSolve::kSolMinRes) 
1067     res =  slv->SolveMinRes(sol,fgMinResCondType,fgMinResMaxIter,fgMinResTol);
1068   else if (fgIterSol == AliMinResSolve::kSolFGMRes) 
1069     res =  slv->SolveFGMRES(sol,fgMinResCondType,fgMinResMaxIter,fgMinResTol,fgNKrylovV);
1070   else 
1071     AliInfo(Form("Undefined Iteritive Solver ID=%d, only %d are defined",fgIterSol,AliMinResSolve::kNSolvers));
1072   //
1073   if (!res) {
1074     const char* faildump = "fgmr_failed.dat";
1075     int defout = dup(1);
1076     if (defout<0) {
1077       AliInfo("Failed on dup");
1078       return kFailed;
1079     }
1080     int slvDump = open(faildump, O_RDWR|O_CREAT, 0666);
1081     if (slvDump>=0) {
1082       dup2(slvDump,1);
1083       //
1084       printf("#Failed to solve using solver %d with PreCond: %d MaxIter: %d Tol: %e NKrylov: %d\n",
1085              fgIterSol,fgMinResCondType,fgMinResMaxIter,fgMinResTol,fgNKrylovV);
1086       printf("#Dump of matrix:\n");
1087       fMatCGlo->Print("10");
1088       printf("#Dump of RHS:\n");
1089       for (int i=0;i<fNGloSize;i++) printf("%d %+.10f\n",i,fVecBGlo[i]);
1090       //
1091       dup2(defout,1);
1092       close(slvDump);
1093       printf("#Dumped failed matrix and RHS to %s\n",faildump);
1094     }
1095     else AliInfo("Failed on file open for matrix dumping");
1096     close(defout);
1097     return kFailed;
1098   }
1099   for (int i=fNGloSize;i--;) fVecBGlo[i] = sol[i];
1100   return kNoInversion;
1101   //
1102 }
1103
1104 //_____________________________________________________________________________________________
1105 Float_t AliMillePede2::Chi2DoFLim(int nSig, int nDoF) const
1106 {
1107   /// return the limit in chi^2/nd for n sigmas stdev authorized
1108   // Only n=1, 2, and 3 are expected in input
1109   int lNSig;
1110   float sn[3]        =  {0.47523, 1.690140, 2.782170};
1111   float table[3][30] = {{1.0000, 1.1479, 1.1753, 1.1798, 1.1775, 1.1730, 1.1680, 1.1630,
1112                          1.1581, 1.1536, 1.1493, 1.1454, 1.1417, 1.1383, 1.1351, 1.1321,
1113                          1.1293, 1.1266, 1.1242, 1.1218, 1.1196, 1.1175, 1.1155, 1.1136,
1114                          1.1119, 1.1101, 1.1085, 1.1070, 1.1055, 1.1040},
1115                         {4.0000, 3.0900, 2.6750, 2.4290, 2.2628, 2.1415, 2.0481, 1.9736,
1116                          1.9124, 1.8610, 1.8171, 1.7791, 1.7457, 1.7161, 1.6897, 1.6658,
1117                          1.6442, 1.6246, 1.6065, 1.5899, 1.5745, 1.5603, 1.5470, 1.5346,
1118                          1.5230, 1.5120, 1.5017, 1.4920, 1.4829, 1.4742},
1119                         {9.0000, 5.9146, 4.7184, 4.0628, 3.6410, 3.3436, 3.1209, 2.9468,
1120                          2.8063, 2.6902, 2.5922, 2.5082, 2.4352, 2.3711, 2.3143, 2.2635,
1121                          2.2178, 2.1764, 2.1386, 2.1040, 2.0722, 2.0428, 2.0155, 1.9901,
1122                          1.9665, 1.9443, 1.9235, 1.9040, 1.8855, 1.8681}};
1123   
1124   if (nDoF < 1) {
1125     return 0.0;
1126   }
1127   else {  
1128     lNSig = TMath::Max(1,TMath::Min(nSig,3));
1129     
1130     if (nDoF <= 30) {    
1131       return table[lNSig-1][nDoF-1];
1132     }
1133     else { // approximation
1134       return ((sn[lNSig-1]+TMath::Sqrt(float(2*nDoF-3)))*
1135               (sn[lNSig-1]+TMath::Sqrt(float(2*nDoF-3))))/float(2*nDoF-2);
1136     }
1137   }
1138 }
1139
1140 //_____________________________________________________________________________________________
1141 Int_t AliMillePede2::SetIterations(double lChi2CutFac)
1142 {
1143   // Number of iterations is calculated from lChi2CutFac 
1144   fChi2CutFactor = TMath::Max(1.0, lChi2CutFac);
1145   //
1146   AliInfo(Form("Initial cut factor is %f",fChi2CutFactor));
1147   fIter = 1; // Initializes the iteration process
1148   return 1;
1149 }
1150
1151 //_____________________________________________________________________________________________
1152 Double_t AliMillePede2::GetParError(int iPar) const
1153 {
1154   // return error for parameter iPar
1155   if (fGloSolveStatus==kInvert) {
1156     double res = fMatCGlo->QueryDiag(iPar);
1157     if (res>=0) return TMath::Sqrt(res);
1158   } 
1159   return -1.;
1160 }
1161
1162
1163 //_____________________________________________________________________________________________
1164 Int_t AliMillePede2::PrintGlobalParameters() const
1165 {
1166   ///  Print the final results into the logfile
1167   double lError = 0.;
1168   double lGlobalCor =0.;
1169         
1170   AliInfo("");
1171   AliInfo("   Result of fit for global parameters");
1172   AliInfo("   ===================================");
1173   AliInfo("    I       initial       final       differ        lastcor        error       gcor       Npnt");
1174   AliInfo("----------------------------------------------------------------------------------------------");
1175   //
1176   for (int i=0; i<fNGloPar; i++) {
1177     lError = GetParError(i);
1178     lGlobalCor = 0.0;
1179     //          
1180     double dg;
1181     if (fGloSolveStatus==kInvert && TMath::Abs( (dg=fMatCGlo->QueryDiag(i)) *fDiagCGlo[i]) > 0) {    
1182       lGlobalCor = TMath::Sqrt(TMath::Abs(1.0-1.0/(dg*fDiagCGlo[i])));
1183       AliInfo(Form("%d\t %.6f\t %.6f\t %.6f\t %.6f\t %.6f\t %.6f\t %6d",
1184                    i,fInitPar[i],fInitPar[i]+fDeltaPar[i],fDeltaPar[i],fVecBGlo[i],lError,lGlobalCor,fProcPnt[i]));
1185     }
1186     else {    
1187       AliInfo(Form("%d\t %.6f\t %.6f\t %.6f\t %.6f\t OFF\t OFF\t %6d",i,fInitPar[i],fInitPar[i]+fDeltaPar[i],
1188                    fDeltaPar[i],fVecBGlo[i],fProcPnt[i]));
1189     }
1190   }
1191   return 1;
1192 }
1193
1194 //_____________________________________________________________________________________________
1195 Bool_t AliMillePede2::IsRecordAcceptable() const
1196 {
1197   // validate record according run lists set by the user
1198   static Long_t prevRunID = kMaxInt;
1199   static Bool_t prevAns   = kTRUE;
1200   Long_t runID = fRecord->GetRunID();
1201   if (runID!=prevRunID) {
1202     int n = 0;
1203     prevRunID = runID;
1204     // is run to be rejected?
1205     if (fRejRunList && (n=fRejRunList->GetSize())) {
1206       prevAns = kTRUE;
1207       for (int i=n;i--;) if (runID == (*fRejRunList)[i]) {
1208           prevAns = kFALSE; 
1209           printf("New Run to reject: %ld -> %d\n",runID,prevAns);
1210           break;
1211         }
1212     }
1213     else if (fAccRunList && (n=fAccRunList->GetSize())) {     // is run specifically selected
1214       prevAns = kFALSE;
1215       for (int i=n;i--;) if (runID == (*fAccRunList)[i]) {prevAns = kTRUE; break;}
1216     }
1217   }
1218   //
1219   return prevAns;
1220   //
1221 }
1222
1223 //_____________________________________________________________________________________________
1224 void AliMillePede2::SetRejRunList(const UInt_t *runs, Int_t nruns)
1225 {
1226   // set the list of runs to be rejected
1227   if (fRejRunList) delete fRejRunList; 
1228   fRejRunList = 0;
1229   if (nruns<1 || !runs) return;
1230   fRejRunList = new TArrayL(nruns);
1231   for (int i=0;i<nruns;i++) (*fRejRunList)[i] = runs[i];
1232 }
1233
1234 //_____________________________________________________________________________________________
1235 void AliMillePede2::SetAccRunList(const UInt_t *runs, Int_t nruns)
1236 {
1237   // set the list of runs to be selected
1238   if (fAccRunList) delete fAccRunList;
1239   fAccRunList = 0;
1240   if (nruns<1 || !runs) return;
1241   fAccRunList = new TArrayL(nruns);
1242   for (int i=0;i<nruns;i++) (*fAccRunList)[i] = runs[i];
1243 }