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