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