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