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