1 /**********************************************************************************************/
2 /* General class for alignment with large number of degrees of freedom */
3 /* Based on the original milliped2 by Volker Blobel */
4 /* http://www.desy.de/~blobel/mptalks.html */
6 /* Author: ruben.shahoyan@cern.ch */
8 /**********************************************************************************************/
10 #include "AliMillePede2.h"
12 #include <TStopwatch.h>
17 #include "AliMatrixSq.h"
18 #include "AliSymMatrix.h"
19 #include "AliRectMatrix.h"
20 #include "AliMatrixSparse.h"
24 #include <sys/types.h>
31 ClassImp(AliMillePede2)
33 Bool_t AliMillePede2::fgInvChol = kTRUE; // Invert global matrix with Cholesky solver
34 Bool_t AliMillePede2::fgWeightSigma = kTRUE; // weight local constraint by module statistics
35 Bool_t AliMillePede2::fgIsMatGloSparse = kFALSE; // use faster dense matrix by default
36 Int_t AliMillePede2::fgMinResCondType = 1; // Jacoby preconditioner by default
37 Double_t AliMillePede2::fgMinResTol = 1.e-11; // default tolerance
38 Int_t AliMillePede2::fgMinResMaxIter = 10000; // default max number of iterations
39 Int_t AliMillePede2::fgIterSol = AliMinResSolve::kSolMinRes; // default iterative solver
40 Int_t AliMillePede2::fgNKrylovV = 240; // default number of Krylov vectors to keep
42 //_____________________________________________________________________________________________
43 AliMillePede2::AliMillePede2()
53 fNLagrangeConstraints(0),
57 fGloSolveStatus(gkFailed),
87 fDataRecFName("/tmp/mp2_data_records.root"),
93 fConstrRecFName("/tmp/mp2_constraints_records.root"),
105 //_____________________________________________________________________________________________
106 AliMillePede2::AliMillePede2(const AliMillePede2& src) :
107 TObject(src),fNLocPar(0),fNGloPar(0),fNGloSize(0),fNLocEquations(0),fIter(0),
108 fMaxIter(10),fNStdDev(3),fNGloConstraints(0),fNLagrangeConstraints(0),
109 fNLocFits(0),fNLocFitsRejected(0),
110 fNGloFix(0),fGloSolveStatus(0),fChi2CutFactor(0),fChi2CutRef(0),fResCutInit(0),
111 fResCut(0),fMinPntValid(1),fNGroupsSet(0),fParamGrID(0),fProcPnt(0),fVecBLoc(0),fDiagCGlo(0),fVecBGlo(0),
112 fInitPar(0),fDeltaPar(0),fSigmaPar(0),fIsLinear(0),fConstrUsed(0),fGlo2CGlo(0),fCGlo2Glo(0),
113 fMatCLoc(0),fMatCGlo(0),fMatCGloLoc(0),fFillIndex(0),fFillValue(0),
114 fDataRecFName(0),fRecord(0),fDataRecFile(0),
115 fTreeData(0),fRecFileStatus(0),fConstrRecFName(0),fTreeConstr(0),fConsRecFile(0),fCurrRecDataID(0),
116 fCurrRecConstrID(0),fLocFitAdd(kTRUE),
123 //_____________________________________________________________________________________________
124 AliMillePede2::~AliMillePede2()
126 CloseDataRecStorage();
127 CloseConsRecStorage();
140 delete[] fConstrUsed;
152 //_____________________________________________________________________________________________
153 Int_t AliMillePede2::InitMille(int nGlo, int nLoc, int lNStdDev,double lResCut, double lResCutInit)
156 if (nLoc>0) fNLocPar = nLoc;
157 if (nGlo>0) fNGloPar = nGlo;
158 if (lResCutInit>0) fResCutInit = lResCutInit;
159 if (lResCut>0) fResCut = lResCut;
160 if (lNStdDev>0) fNStdDev = lNStdDev;
162 fNGloSize = fNGloPar;
166 if (fgIsMatGloSparse) {fMatCGlo = new AliMatrixSparse(fNGloPar); fMatCGlo->SetSymmetric(kTRUE);}
167 else fMatCGlo = new AliSymMatrix(fNGloPar);
169 fFillIndex = new Int_t[fNGloPar];
170 fFillValue = new Double_t[fNGloPar];
172 fMatCLoc = new AliSymMatrix(fNLocPar);
173 fMatCGloLoc = new AliRectMatrix(fNGloPar,fNLocPar);
175 fParamGrID = new Int_t[fNGloPar];
176 fProcPnt = new Int_t[fNGloPar];
177 fVecBLoc = new Double_t[fNLocPar];
178 fDiagCGlo = new Double_t[fNGloPar];
180 fInitPar = new Double_t[fNGloPar];
181 fDeltaPar = new Double_t[fNGloPar];
182 fSigmaPar = new Double_t[fNGloPar];
183 fIsLinear = new Bool_t[fNGloPar];
185 fGlo2CGlo = new Int_t[fNGloPar];
186 fCGlo2Glo = new Int_t[fNGloPar];
189 AliInfo(Form("Failed to allocate the memory for %d global and %d local parameters",fNGloPar,fNLocPar));
193 memset(fVecBLoc ,0,fNLocPar*sizeof(Double_t));
194 memset(fDiagCGlo ,0,fNGloPar*sizeof(Double_t));
195 memset(fInitPar ,0,fNGloPar*sizeof(Double_t));
196 memset(fDeltaPar ,0,fNGloPar*sizeof(Double_t));
197 memset(fSigmaPar ,0,fNGloPar*sizeof(Double_t));
198 memset(fProcPnt ,0,fNGloPar*sizeof(Int_t));
200 for (int i=fNGloPar;i--;) {
201 fGlo2CGlo[i] = fCGlo2Glo[i] = -1;
202 fIsLinear[i] = kTRUE;
209 //_____________________________________________________________________________________________
210 Bool_t AliMillePede2::ImposeDataRecFile(const char* fname)
212 CloseDataRecStorage();
213 SetDataRecFName(fname);
214 return InitDataRecStorage(kTRUE); // open in read mode
217 //_____________________________________________________________________________________________
218 Bool_t AliMillePede2::ImposeConsRecFile(const char* fname)
220 CloseConsRecStorage();
221 SetConsRecFName(fname);
222 return InitConsRecStorage(kTRUE); // open in read mode
225 //_____________________________________________________________________________________________
226 Bool_t AliMillePede2::InitDataRecStorage(Bool_t read)
228 // initialize the buffer for processed measurements records
230 if (fDataRecFile) {AliInfo("Data Records File is already initialized"); return kFALSE;}
232 if (!fRecord) fRecord = new AliMillePedeRecord();
234 fDataRecFile = TFile::Open(GetDataRecFName(),read ? "":"recreate");
235 if (!fDataRecFile) {AliFatal(Form("Failed to initialize data records file %s",GetDataRecFName())); return kFALSE;}
237 AliInfo(Form("File %s used for derivatives records",GetDataRecFName()));
239 fTreeData = (TTree*)fDataRecFile->Get("AliMillePedeRecords_Data");
240 if (!fTreeData) {AliFatal(Form("Did not find data records tree in %s",GetDataRecFName())); return kFALSE;}
241 fTreeData->SetBranchAddress("Record_Data",&fRecord);
242 AliInfo(Form("Found %lld derivatives records",fTreeData->GetEntries()));
245 fTreeData = new TTree("AliMillePedeRecords_Data","Data Records for AliMillePede2");
246 fTreeData->Branch("Record_Data","AliMillePedeRecord",&fRecord,32000,99);
249 fRecFileStatus = read ? 1:2;
254 //_____________________________________________________________________________________________
255 Bool_t AliMillePede2::InitConsRecStorage(Bool_t read)
257 // initialize the buffer for processed measurements records
259 if (fConsRecFile) {AliInfo("Constraints Records File is already initialized"); return kFALSE;}
261 if (!fRecord) fRecord = new AliMillePedeRecord();
263 fConsRecFile = TFile::Open(GetConsRecFName(),read ? "":"recreate");
264 if (!fConsRecFile) {AliInfo(Form("Failed to initialize constraints records file %s",GetConsRecFName())); return kFALSE;}
266 AliInfo(Form("File %s used for constraints records",GetConsRecFName()));
268 fTreeConstr = (TTree*)fConsRecFile->Get("AliMillePedeRecords_Constraints");
269 if (!fTreeConstr) {AliInfo(Form("Did not find constraints records tree in %s",GetConsRecFName())); return kFALSE;}
270 fTreeConstr->SetBranchAddress("Record_Constraints",&fRecord);
271 AliInfo(Form("Found %lld constraints records",fTreeConstr->GetEntries()));
276 fTreeConstr = new TTree("AliMillePedeRecords_Constraints","Constraints Records for AliMillePede2");
277 fTreeConstr->Branch("Record_Constraints","AliMillePedeRecord",&fRecord,32000,99);
279 fCurrRecConstrID = -1;
284 //_____________________________________________________________________________________________
285 void AliMillePede2::CloseDataRecStorage()
288 if (fDataRecFile->IsWritable()) {
294 fDataRecFile->Close();
302 //_____________________________________________________________________________________________
303 void AliMillePede2::CloseConsRecStorage()
306 if (fConsRecFile->IsWritable()) {
308 fTreeConstr->Write();
312 fConsRecFile->Close();
319 //_____________________________________________________________________________________________
320 Bool_t AliMillePede2::ReadNextRecordData()
322 // read next data record (if any)
323 if (!fTreeData || ++fCurrRecDataID >= fTreeData->GetEntries()) { fCurrRecDataID--; return kFALSE;}
324 fTreeData->GetEntry(fCurrRecDataID);
328 //_____________________________________________________________________________________________
329 Bool_t AliMillePede2::ReadNextRecordConstraint()
331 // read next constraint record (if any)
332 if (!fTreeConstr || ++fCurrRecConstrID >= fTreeConstr->GetEntries()) { fCurrRecConstrID--; return kFALSE;}
333 fTreeConstr->GetEntry(fCurrRecConstrID);
337 //_____________________________________________________________________________________________
338 void AliMillePede2::SetRecordWeight(double wgh)
340 if (fRecFileStatus<2) InitDataRecStorage(); // create a buffer to store the data
341 fRecord->SetWeight(wgh);
344 //_____________________________________________________________________________________________
345 void AliMillePede2::SetRecordRun(Int_t run)
347 if (fRecFileStatus<2) InitDataRecStorage(); // create a buffer to store the data
348 fRecord->SetRunID(run);
351 //_____________________________________________________________________________________________
352 void AliMillePede2::SetLocalEquation(double *dergb, double *derlc, double lMeas, double lSigma)
354 if (fRecFileStatus<2) InitDataRecStorage(); // create a buffer to store the data
356 // write data of single measurement
357 if (lSigma<=0.0) { // If parameter is fixed, then no equation
358 for (int i=fNLocPar; i--;) derlc[i] = 0.0;
359 for (int i=fNGloPar; i--;) dergb[i] = 0.0;
363 fRecord->AddResidual(lMeas);
365 // Retrieve local param interesting indices
366 for (int i=0;i<fNLocPar;i++) if (!IsZero(derlc[i])) {fRecord->AddIndexValue(i,derlc[i]); derlc[i] = 0.0;}
368 fRecord->AddWeight( 1.0/lSigma/lSigma );
370 // Idem for global parameters
371 for (int i=0;i<fNGloPar;i++) if (!IsZero(dergb[i])) {
372 fRecord->AddIndexValue(i,dergb[i]); dergb[i] = 0.0;
373 fRecord->MarkGroup(fParamGrID[i]);
378 //_____________________________________________________________________________________________
379 void AliMillePede2::SetLocalEquation(int *indgb, double *dergb, int ngb, int *indlc,
380 double *derlc,int nlc,double lMeas,double lSigma)
382 // write data of single measurement
383 if (lSigma<=0.0) { // If parameter is fixed, then no equation
384 for (int i=nlc;i--;) derlc[i] = 0.0;
385 for (int i=ngb;i--;) dergb[i] = 0.0;
389 if (fRecFileStatus<2) InitDataRecStorage(); // create a buffer to store the data
391 fRecord->AddResidual(lMeas);
393 // Retrieve local param interesting indices
394 for (int i=0;i<nlc;i++) if (!IsZero(derlc[i])) {fRecord->AddIndexValue(indlc[i],derlc[i]); derlc[i]=0.; indlc[i]=0;}
396 fRecord->AddWeight( 1./lSigma/lSigma );
398 // Idem for global parameters
399 for (int i=0;i<ngb;i++) if (!IsZero(dergb[i])) {fRecord->AddIndexValue(indgb[i],dergb[i]); dergb[i]=0.; indgb[i]=0;}
404 //_____________________________________________________________________________________________
405 void AliMillePede2::SetGlobalConstraint(double *dergb, double val, double sigma)
407 // Define a constraint equation.
408 if (!fConsRecFile || !fConsRecFile->IsWritable()) InitConsRecStorage(); // create a buffer to store the data
411 fRecord->AddResidual(val);
412 fRecord->AddWeight(sigma);
413 for (int i=0; i<fNGloPar; i++) if (!IsZero(dergb[i])) fRecord->AddIndexValue(i,dergb[i]);
415 if (IsZero(sigma)) fNLagrangeConstraints++;
416 // printf("NewConstraint:\n"); fRecord->Print(); //RRR
417 SaveRecordConstraint();
421 //_____________________________________________________________________________________________
422 void AliMillePede2::SetGlobalConstraint(const int *indgb, double *dergb, int ngb, double val,double sigma)
424 // Define a constraint equation.
425 if (!fConsRecFile || !fConsRecFile->IsWritable()) InitConsRecStorage(); // create a buffer to store the data
427 fRecord->AddResidual(val);
428 fRecord->AddWeight(sigma); // dummy
429 for (int i=0; i<ngb; i++) if (!IsZero(dergb[i])) fRecord->AddIndexValue(indgb[i],dergb[i]);
431 if (IsZero(sigma)) fNLagrangeConstraints++;
432 SaveRecordConstraint();
436 //_____________________________________________________________________________________________
437 Int_t AliMillePede2::LocalFit(double *localParams)
440 Perform local parameters fit once all the local equations have been set
441 -----------------------------------------------------------
442 localParams = (if !=0) will contain the fitted track parameters and
445 static int nrefSize = 0;
446 // static TArrayI refLoc,refGlo,nrefLoc,nrefGlo;
447 static Int_t *refLoc=0,*refGlo=0,*nrefLoc=0,*nrefGlo=0;
450 AliSymMatrix &matCLoc = *fMatCLoc;
451 AliMatrixSq &matCGlo = *fMatCGlo;
452 AliRectMatrix &matCGloLoc = *fMatCGloLoc;
454 memset(fVecBLoc,0,fNLocPar*sizeof(double));
458 int recSz = fRecord->GetSize();
460 while(cnt<recSz) { // Transfer the measurement records to matrices
462 // extract addresses of residual, weight and pointers on local and global derivatives for each point
463 if (nrefSize<=nPoints) {
465 nrefSize = 2*(nPoints+1);
466 tmpA = refLoc; refLoc = new Int_t[nrefSize]; if (tmpA) memcpy(refLoc,tmpA,nPoints*sizeof(int));
467 tmpA = refGlo; refGlo = new Int_t[nrefSize]; if (tmpA) memcpy(refGlo,tmpA,nPoints*sizeof(int));
468 tmpA = nrefLoc; nrefLoc = new Int_t[nrefSize]; if (tmpA) memcpy(nrefLoc,tmpA,nPoints*sizeof(int));
469 tmpA = nrefGlo; nrefGlo = new Int_t[nrefSize]; if (tmpA) memcpy(nrefGlo,tmpA,nPoints*sizeof(int));
472 refLoc[nPoints] = ++cnt;
474 while(!fRecord->IsWeight(cnt)) {nLoc++; cnt++;}
475 nrefLoc[nPoints] = nLoc;
477 refGlo[nPoints] = ++cnt;
479 while(!fRecord->IsResidual(cnt) && cnt<recSz) {nGlo++; cnt++;}
480 nrefGlo[nPoints] = nGlo;
486 double gloWgh = fRecord->GetWeight(); // global weight for this set
487 Int_t maxLocUsed = 0;
489 for (int ip=nPoints;ip--;) { // Transfer the measurement records to matrices
490 double resid = fRecord->GetValue( refLoc[ip]-1 );
491 double weight = fRecord->GetValue( refGlo[ip]-1 )*gloWgh;
492 double *derLoc = fRecord->GetValue()+refLoc[ip];
493 double *derGlo = fRecord->GetValue()+refGlo[ip];
494 int *indLoc = fRecord->GetIndex()+refLoc[ip];
495 int *indGlo = fRecord->GetIndex()+refGlo[ip];
497 for (int i=nrefGlo[ip];i--;) { // suppress the global part (only relevant with iterations)
498 int iID = indGlo[i]; // Global param indice
499 if (fSigmaPar[iID]<=0.) continue; // fixed parameter RRRCheck
500 if (fIsLinear[iID]) resid -= derGlo[i]*(fInitPar[iID]+fDeltaPar[iID]); // linear parameter
501 else resid -= derGlo[i]*fDeltaPar[iID]; // nonlinear parameter
504 // Symmetric matrix, don't bother j>i coeffs
505 for (int i=nrefLoc[ip];i--;) { // Fill local matrix and vector
506 fVecBLoc[ indLoc[i] ] += weight*resid*derLoc[i];
507 if (indLoc[i]>maxLocUsed) maxLocUsed = indLoc[i];
508 for (int j=i+1;j--;) matCLoc(indLoc[i] ,indLoc[j]) += weight*derLoc[i]*derLoc[j];
511 } // end of the transfer of the measurement record to matrices
513 matCLoc.SetSizeUsed(++maxLocUsed); // data with B=0 may use less than declared nLocals
515 // first try to solve by faster Cholesky decomposition, then by Gaussian elimination
516 if (!matCLoc.SolveChol(fVecBLoc,kTRUE)) {
517 AliInfo("Failed to solve locals by Cholesky, trying Gaussian Elimination");
518 if (!matCLoc.SolveSpmInv(fVecBLoc,kTRUE)) {
519 AliInfo("Failed to solve locals by Gaussian Elimination, skip...");
521 return 0; // failed to solve
525 // If requested, store the track params and errors
526 // printf("locfit: "); for (int i=0;i<fNLocPar;i++) printf("%+e |",fVecBLoc[i]); printf("\n");
527 if (localParams) for (int i=maxLocUsed; i--;) {
528 localParams[2*i] = fVecBLoc[i];
529 localParams[2*i+1] = TMath::Sqrt(TMath::Abs(matCLoc.QueryDiag(i)));
535 for (int ip=nPoints;ip--;) { // Calculate residuals
536 double resid = fRecord->GetValue( refLoc[ip]-1 );
537 double weight = fRecord->GetValue( refGlo[ip]-1 )*gloWgh;
538 double *derLoc = fRecord->GetValue()+refLoc[ip];
539 double *derGlo = fRecord->GetValue()+refGlo[ip];
540 int *indLoc = fRecord->GetIndex()+refLoc[ip];
541 int *indGlo = fRecord->GetIndex()+refGlo[ip];
543 // Suppress local and global contribution in residuals;
544 for (int i=nrefLoc[ip];i--;) resid -= derLoc[i]*fVecBLoc[ indLoc[i] ]; // local part
546 for (int i=nrefGlo[ip];i--;) { // global part
548 if ( fSigmaPar[iID] <= 0.) continue; // fixed parameter RRRCheck
549 if (fIsLinear[iID]) resid -= derGlo[i]*(fInitPar[iID]+fDeltaPar[iID]); // linear parameter
550 else resid -= derGlo[i]*fDeltaPar[iID]; // nonlinear parameter
553 // reject the track if the residual is too large (outlier)
554 double absres = TMath::Abs(resid);
555 if ( (absres >= fResCutInit && fIter ==1 ) ||
556 (absres >= fResCut && fIter > 1)) {
557 if (fLocFitAdd) fNLocFitsRejected++;
558 // printf("reject res %5ld %+e\n",fCurrRecDataID,resid);
562 lChi2 += weight*resid*resid ; // total chi^2
563 nEq++; // number of equations
564 } // end of Calculate residuals
567 int nDoF = nEq-maxLocUsed;
568 lChi2 = (nDoF>0) ? lChi2/nDoF : 0; // Chi^2/dof
570 if (fNStdDev != 0 && nDoF>0 && lChi2 > Chi2DoFLim(fNStdDev,nDoF)*fChi2CutFactor) { // check final chi2
571 if (fLocFitAdd) fNLocFitsRejected++;
572 // printf("reject chi2 %5ld: %+e\n",fCurrRecDataID, lChi2);
578 fNLocEquations += nEq;
582 fNLocEquations -= nEq;
585 // local operations are finished, track is accepted
586 // We now update the global parameters (other matrices)
590 for (int ip=nPoints;ip--;) { // Update matrices
591 double resid = fRecord->GetValue( refLoc[ip]-1 );
592 double weight = fRecord->GetValue( refGlo[ip]-1 )*gloWgh;
593 double *derLoc = fRecord->GetValue()+refLoc[ip];
594 double *derGlo = fRecord->GetValue()+refGlo[ip];
595 int *indLoc = fRecord->GetIndex()+refLoc[ip];
596 int *indGlo = fRecord->GetIndex()+refGlo[ip];
598 for (int i=nrefGlo[ip];i--;) { // suppress the global part
599 int iID = indGlo[i]; // Global param indice
600 if ( fSigmaPar[iID] <= 0.) continue; // fixed parameter RRRCheck
601 if (fIsLinear[iID]) resid -= derGlo[i]*(fInitPar[iID]+fDeltaPar[iID]); // linear parameter
602 else resid -= derGlo[i]*fDeltaPar[iID]; // nonlinear parameter
605 for (int ig=nrefGlo[ip];ig--;) {
606 int iIDg = indGlo[ig]; // Global param indice (the matrix line)
607 if ( fSigmaPar[iIDg] <= 0.) continue; // fixed parameter RRRCheck
608 if (fLocFitAdd) fVecBGlo[ iIDg ] += weight*resid*derGlo[ig]; //!!!
609 else fVecBGlo[ iIDg ] -= weight*resid*derGlo[ig]; //!!!
611 // First of all, the global/global terms (exactly like local matrix)
613 for (int jg=ig+1;jg--;) { // matCGlo is symmetric by construction
614 int jIDg = indGlo[jg];
615 if ( fSigmaPar[jIDg] <= 0.) continue; // fixed parameter RRRCheck
616 if ( !IsZero(vl = weight*derGlo[ig]*derGlo[jg]) ) {
617 fFillIndex[nfill] = jIDg;
618 fFillValue[nfill++] = fLocFitAdd ? vl:-vl;
621 if (nfill) matCGlo.AddToRow(iIDg,fFillValue,fFillIndex,nfill);
623 // Now we have also rectangular matrices containing global/local terms.
624 int iCIDg = fGlo2CGlo[iIDg]; // compressed Index of index
626 Double_t *rowGL = matCGloLoc(nGloInFit);
627 for (int k=maxLocUsed;k--;) rowGL[k] = 0.0; // reset the row
628 iCIDg = fGlo2CGlo[iIDg] = nGloInFit;
629 fCGlo2Glo[nGloInFit++] = iIDg;
632 Double_t *rowGLIDg = matCGloLoc(iCIDg);
633 for (int il=nrefLoc[ip];il--;) rowGLIDg[ indLoc[il] ] += weight*derGlo[ig]*derLoc[il];
634 fProcPnt[iIDg] += fLocFitAdd ? 1:-1; // update counter
637 } // end of Update matrices
639 // calculate fMatCGlo -= fMatCGloLoc * fMatCLoc * fMatCGloLoc^T
640 // and fVecBGlo -= fMatCGloLoc * fVecBLoc
642 //-------------------------------------------------------------- >>>
644 for (int iCIDg=0; iCIDg<nGloInFit; iCIDg++) {
645 int iIDg = fCGlo2Glo[iCIDg];
648 Double_t *rowGLIDg = matCGloLoc(iCIDg);
649 for (int kl=0;kl<maxLocUsed;kl++) if (rowGLIDg[kl]) vl += rowGLIDg[kl]*fVecBLoc[kl];
650 if (!IsZero(vl)) fVecBGlo[iIDg] -= fLocFitAdd ? vl : -vl;
653 for (int jCIDg=0;jCIDg<=iCIDg; jCIDg++) {
654 int jIDg = fCGlo2Glo[jCIDg];
657 Double_t *rowGLJDg = matCGloLoc(jCIDg);
658 for (int kl=0;kl<maxLocUsed;kl++) {
660 if ( (!IsZero(vll=rowGLIDg[kl]*rowGLJDg[kl])) ) vl += matCLoc.QueryDiag(kl)*vll;
663 for (int ll=0;ll<kl;ll++) {
664 if ( !IsZero(vll=rowGLIDg[kl]*rowGLJDg[ll]) ) vl += matCLoc(kl,ll)*vll;
665 if ( !IsZero(vll=rowGLIDg[ll]*rowGLJDg[kl]) ) vl += matCLoc(kl,ll)*vll;
669 fFillIndex[nfill] = jIDg;
670 fFillValue[nfill++] = fLocFitAdd ? -vl : vl;
673 if (nfill) matCGlo.AddToRow(iIDg,fFillValue,fFillIndex,nfill);
676 // reset compressed index array
678 for (int i=nGloInFit;i--;) {
679 fGlo2CGlo[ fCGlo2Glo[i] ] = -1;
683 //---------------------------------------------------- <<<
687 //_____________________________________________________________________________________________
688 Int_t AliMillePede2::GlobalFit(Double_t *par, Double_t *error, Double_t *pull)
690 // performs a requested number of global iterations
693 TStopwatch sw; sw.Start();
696 AliInfo("Starting Global fit.");
697 while (fIter<=fMaxIter) {
699 res = GlobalFitIteration();
702 if (!IsZero(fChi2CutFactor-fChi2CutRef)) {
703 fChi2CutFactor = TMath::Sqrt(fChi2CutFactor);
704 if (fChi2CutFactor < 1.2*fChi2CutRef) {
705 fChi2CutFactor = fChi2CutRef;
706 //RRR fIter = fMaxIter - 1; // Last iteration
713 AliInfo(Form("Global fit %s, CPU time: %.1f",res ? "Converged":"Failed",sw.CpuTime()));
716 if (par) for (int i=fNGloPar;i--;) par[i] = fInitPar[i]+fDeltaPar[i];
718 if (fGloSolveStatus==gkInvert) { // errors on params are available
719 if (error) for (int i=fNGloPar;i--;) error[i] = fProcPnt[i]>0 ? TMath::Sqrt(TMath::Abs(fMatCGlo->QueryDiag(i))) : 0.;
720 if (pull) for (int i=fNGloPar;i--;) pull[i] = fProcPnt[i]>0 && (fSigmaPar[i]*fSigmaPar[i]-fMatCGlo->QueryDiag(i))>0. && fSigmaPar[i]>0
721 ? fDeltaPar[i]/TMath::Sqrt(fSigmaPar[i]*fSigmaPar[i]-fMatCGlo->QueryDiag(i)) : 0;
727 //_____________________________________________________________________________________________
728 Int_t AliMillePede2::GlobalFitIteration()
730 // perform global parameters fit once all the local equations have been fitted
732 AliInfo(Form("Global Fit Iteration#%2d (Local Fit Chi^2 cut factor: %.2f)",fIter,fChi2CutFactor));
734 if (!fNGloPar || !fTreeData) {
735 AliInfo("No data was stored, aborting iteration");
743 fConstrUsed = new Bool_t[fNGloConstraints];
744 memset(fConstrUsed,0,fNGloConstraints*sizeof(Bool_t));
746 // Reset all info specific for this step
747 AliMatrixSq& matCGlo = *fMatCGlo;
749 memset(fProcPnt,0,fNGloPar*sizeof(Int_t));
751 fNGloConstraints = fTreeConstr ? fTreeConstr->GetEntries() : 0;
753 // count number of Lagrange constraints: they need new row/cols to be added
754 fNLagrangeConstraints = 0;
755 for (int i=0; i<fNGloConstraints; i++) {
756 ReadRecordConstraint(i);
757 if ( IsZero(fRecord->GetValue(1)) ) fNLagrangeConstraints++; // exact constraint (no error) -> Lagrange multiplier
760 // if needed, readjust the size of the global vector (for matrices this is done automatically)
761 if (!fVecBGlo || fNGloSize!=fNGloPar+fNLagrangeConstraints) {
762 delete[] fVecBGlo; // in case some constraint was added between the two manual iterations
763 fNGloSize = fNGloPar+fNLagrangeConstraints;
764 fVecBGlo = new Double_t[fNGloSize];
766 memset(fVecBGlo,0,fNGloSize*sizeof(double));
769 fNLocFitsRejected = 0;
772 // Process data records and build the matrices
773 Long_t ndr = fTreeData->GetEntries();
774 Long_t first = fSelFirst>0 ? fSelFirst : 0;
775 Long_t last = fSelLast<1 ? ndr : (fSelLast>=ndr ? ndr : fSelLast+Long_t(1));
778 AliInfo(Form("Building the Global matrix from data records %ld : %ld",first,last));
781 TStopwatch swt; swt.Start();
782 fLocFitAdd = kTRUE; // add contributions of matching tracks
783 for (Long_t i=0;i<ndr;i++) {
784 Long_t iev = i+first;
786 if (!IsRecordAcceptable()) continue;
788 if ( (i%int(0.2*ndr)) == 0) printf("%.1f%% of local fits done\n", double(100.*i)/ndr);
791 printf("%ld local fits done: ", ndr);
796 // ---------------------- Reject parameters with low statistics ------------>>
798 if (fMinPntValid>1 && fNGroupsSet) {
800 printf("Checking parameters with statistics < %d\n",fMinPntValid);
803 // 1) build the list of parameters to fix
804 Int_t fixArrSize = 10;
805 Int_t nFixedGroups = 0;
806 TArrayI fixGroups(fixArrSize);
810 double oldMin = 1.e20;
811 double oldMax =-1.e20;
813 for (int i=fNGloPar;i--;) { // // Reset row and column of fixed params and add 1/sig^2 to free ones
814 int grID = fParamGrID[i];
815 if (grID<0) continue; // not in the group
817 if (grID!=grIDold) { // starting new group
818 if (grIDold>=0) { // decide if the group has enough statistics
819 if (oldMin<fMinPntValid && oldMax<2*fMinPntValid) { // suppress group
820 for (int iold=oldStart;iold>i;iold--) fProcPnt[iold] = 0;
821 Bool_t fnd = kFALSE; // check if the group is already accounted
822 for (int j=nFixedGroups;j--;) if (fixGroups[j]==grIDold) {fnd=kTRUE; break;}
824 if (nFixedGroups>=fixArrSize) {fixArrSize*=2; fixGroups.Set(fixArrSize);}
825 fixGroups[nFixedGroups++] = grIDold; // add group to fix
829 grIDold = grID; // mark the start of the new group
834 if (oldMin>fProcPnt[i]) oldMin = fProcPnt[i];
835 if (oldMax<fProcPnt[i]) oldMax = fProcPnt[i];
838 // extra check for the last group
839 if (grIDold>=0 && oldMin<fMinPntValid && oldMax<2*fMinPntValid) { // suppress group
840 for (int iold=oldStart;iold--;) fProcPnt[iold] = 0;
841 Bool_t fnd = kFALSE; // check if the group is already accounted
842 for (int j=nFixedGroups;j--;) if (fixGroups[j]==grIDold) {fnd=kTRUE; break;}
844 if (nFixedGroups>=fixArrSize) {fixArrSize*=2; fixGroups.Set(fixArrSize);}
845 fixGroups[nFixedGroups++] = grIDold; // add group to fix
849 // 2) loop over records and add contributions of fixed groups with negative sign
852 for (Long_t i=0;i<ndr;i++) {
853 Long_t iev = i+first;
855 if (!IsRecordAcceptable()) continue;
856 Bool_t suppr = kFALSE;
857 for (int ifx=nFixedGroups;ifx--;)if (fRecord->IsGroupPresent(fixGroups[ifx])) suppr = kTRUE;
858 if (suppr) LocalFit();
863 printf("Suppressed contributions of groups with NPoints<%d :\n",fMinPntValid);
864 for (int i=0;i<nFixedGroups;i++) printf("%d ",fixGroups[i]); printf("\n");
869 // ---------------------- Reject parameters with low statistics ------------<<
871 // add large number to diagonal of fixed params
873 for (int i=fNGloPar;i--;) { // // Reset row and column of fixed params and add 1/sig^2 to free ones
874 // printf("#%3d : Nproc : %5d grp: %d\n",i,fProcPnt[i],fParamGrID[i]);
878 matCGlo.DiagElem(i) = 1.;//float(fNLocEquations*fNLocEquations);
879 // matCGlo.DiagElem(i) = float(fNLocEquations*fNLocEquations);
881 else matCGlo.DiagElem(i) += (fgWeightSigma ? fProcPnt[i] : 1.)/(fSigmaPar[i]*fSigmaPar[i]);
884 for (int i=fNGloPar;i--;) fDiagCGlo[i] = matCGlo.QueryDiag(i); // save the diagonal elements
886 // add constraint equations
887 int nVar = fNGloPar; // Current size of global matrix
888 for (int i=0; i<fNGloConstraints; i++) {
889 ReadRecordConstraint(i);
890 double val = fRecord->GetValue(0);
891 double sig = fRecord->GetValue(1);
892 int *indV = fRecord->GetIndex()+2;
893 double *der = fRecord->GetValue()+2;
894 int csize = fRecord->GetSize()-2;
896 // check if after suppression of fixed variables there are non-0 derivatives
897 // and determine the max statistics of involved params
900 for (int j=csize;j--;) {
901 if (fProcPnt[indV[j]]<1) nSuppressed++;
903 maxStat = TMath::Max(maxStat,fProcPnt[indV[j]]);
907 if (nSuppressed==csize) {
908 // AliInfo(Form("Neglecting constraint %d of %d derivatives since no free parameters left",i,csize));
910 // was this constraint ever created ?
911 if ( sig==0 && fConstrUsed[i] ) { // this is needed only for constraints with Lagrange multiplier
912 // to avoid empty row impose dummy constraint on "Lagrange multiplier"
913 matCGlo.DiagElem(nVar) = 1.;
914 fVecBGlo[nVar++] = 0;
919 // account for already accumulated corrections
920 for (int j=csize; j--;) val -= der[j]*(fInitPar[ indV[j] ]+fDeltaPar[ indV[j] ]);
922 if (sig > 0) { // this is a gaussian constriant: no Lagrange multipliers are added
924 double sig2i = (fgWeightSigma ? TMath::Sqrt(maxStat) : 1.)/sig/sig;
925 for (int ir=0;ir<csize;ir++) {
927 for (int ic=0;ic<=ir;ic++) { // matrix is symmetric
929 double vl = der[ir]*der[ic]*sig2i;
930 if (!IsZero(vl)) matCGlo(iID,jID) += vl;
932 fVecBGlo[iID] += val*der[ir]*sig2i;
935 else { // this is exact constriant: Lagrange multipliers must be added
936 for (int j=csize; j--;) {
938 if (fProcPnt[jID]<1) continue; // this parameter was fixed, don't put it into constraint
939 matCGlo(nVar,jID) = float(fNLocEquations)*der[j]; // fMatCGlo is symmetric, only lower triangle is filled
942 if (matCGlo.QueryDiag(nVar)) matCGlo.DiagElem(nVar) = 0.0;
943 fVecBGlo[nVar++] = float(fNLocEquations)*val; //RS ? should we use here fNLocFits ?
944 fConstrUsed[i] = kTRUE;
948 AliInfo(Form("Obtained %-7ld equations from %-7ld records (%-7ld rejected). Fixed %-4d globals",
949 fNLocEquations,fNLocFits,fNLocFitsRejected,fNGloFix));
953 fGloSolveStatus = SolveGlobalMatEq(); // obtain solution for this step
955 printf("Solve %d |",fIter); sws.Print();
958 AliInfo(Form("Iteration#%2d %s. CPU time: %.1f",fIter,fGloSolveStatus==gkFailed ? "Failed":"Converged",sw.CpuTime()));
959 if (fGloSolveStatus==gkFailed) return 0;
961 for (int i=fNGloPar;i--;) fDeltaPar[i] += fVecBGlo[i]; // Update global parameters values (for iterations)
963 // PrintGlobalParameters();
967 //_____________________________________________________________________________________________
968 Int_t AliMillePede2::SolveGlobalMatEq()
971 // solve global matrix equation MatCGlob*X=VecBGlo and store the result in the VecBGlo
974 printf("GlobalMatrix\n");
977 for (int i=0;i<fNGloPar;i++) printf("%d %+e\n",i,fVecBGlo[i]);
980 if (!fgIsMatGloSparse) {
982 if (fNLagrangeConstraints==0) { // pos-def systems are faster to solve by Cholesky
983 if ( ((AliSymMatrix*)fMatCGlo)->SolveChol(fVecBGlo, fgInvChol) ) return fgInvChol ? gkInvert:gkNoInversion;
984 else AliInfo("Solution of Global Dense System by Cholesky failed, trying Gaussian Elimiation");
987 if (((AliSymMatrix*)fMatCGlo)->SolveSpmInv(fVecBGlo, kTRUE)) return gkInvert;
988 else AliInfo("Solution of Global Dense System by Gaussian Elimination failed, trying iterative methods");
990 // try to solve by minres
991 TVectorD sol(fNGloSize);
993 AliMinResSolve *slv = new AliMinResSolve(fMatCGlo,fVecBGlo);
994 if (!slv) return gkFailed;
997 if (fgIterSol == AliMinResSolve::kSolMinRes)
998 res = slv->SolveMinRes(sol,fgMinResCondType,fgMinResMaxIter,fgMinResTol);
999 else if (fgIterSol == AliMinResSolve::kSolFGMRes)
1000 res = slv->SolveFGMRES(sol,fgMinResCondType,fgMinResMaxIter,fgMinResTol,fgNKrylovV);
1002 AliInfo(Form("Undefined Iteritive Solver ID=%d, only %d are defined",fgIterSol,AliMinResSolve::kNSolvers));
1005 const char* faildump = "fgmr_failed.dat";
1006 int defout = dup(1);
1008 AliInfo("Failed on dup");
1011 int slvDump = open(faildump, O_RDWR|O_CREAT, 0666);
1014 printf("#Failed to solve using solver %d with PreCond: %d MaxIter: %d Tol: %e NKrylov: %d\n",
1015 fgIterSol,fgMinResCondType,fgMinResMaxIter,fgMinResTol,fgNKrylovV);
1016 printf("#Dump of matrix:\n");
1017 fMatCGlo->Print("10");
1018 printf("#Dump of RHS:\n");
1019 for (int i=0;i<fNGloSize;i++) printf("%d %+.10f\n",i,fVecBGlo[i]);
1024 printf("#Dumped failed matrix and RHS to %s\n",faildump);
1027 for (int i=fNGloSize;i--;) fVecBGlo[i] = sol[i];
1028 return gkNoInversion;
1032 //_____________________________________________________________________________________________
1033 Float_t AliMillePede2::Chi2DoFLim(int nSig, int nDoF) const
1035 /// return the limit in chi^2/nd for n sigmas stdev authorized
1036 // Only n=1, 2, and 3 are expected in input
1038 float sn[3] = {0.47523, 1.690140, 2.782170};
1039 float table[3][30] = {{1.0000, 1.1479, 1.1753, 1.1798, 1.1775, 1.1730, 1.1680, 1.1630,
1040 1.1581, 1.1536, 1.1493, 1.1454, 1.1417, 1.1383, 1.1351, 1.1321,
1041 1.1293, 1.1266, 1.1242, 1.1218, 1.1196, 1.1175, 1.1155, 1.1136,
1042 1.1119, 1.1101, 1.1085, 1.1070, 1.1055, 1.1040},
1043 {4.0000, 3.0900, 2.6750, 2.4290, 2.2628, 2.1415, 2.0481, 1.9736,
1044 1.9124, 1.8610, 1.8171, 1.7791, 1.7457, 1.7161, 1.6897, 1.6658,
1045 1.6442, 1.6246, 1.6065, 1.5899, 1.5745, 1.5603, 1.5470, 1.5346,
1046 1.5230, 1.5120, 1.5017, 1.4920, 1.4829, 1.4742},
1047 {9.0000, 5.9146, 4.7184, 4.0628, 3.6410, 3.3436, 3.1209, 2.9468,
1048 2.8063, 2.6902, 2.5922, 2.5082, 2.4352, 2.3711, 2.3143, 2.2635,
1049 2.2178, 2.1764, 2.1386, 2.1040, 2.0722, 2.0428, 2.0155, 1.9901,
1050 1.9665, 1.9443, 1.9235, 1.9040, 1.8855, 1.8681}};
1056 lNSig = TMath::Max(1,TMath::Min(nSig,3));
1059 return table[lNSig-1][nDoF-1];
1061 else { // approximation
1062 return ((sn[lNSig-1]+TMath::Sqrt(float(2*nDoF-3)))*
1063 (sn[lNSig-1]+TMath::Sqrt(float(2*nDoF-3))))/float(2*nDoF-2);
1068 //_____________________________________________________________________________________________
1069 Int_t AliMillePede2::SetIterations(double lChi2CutFac)
1071 // Number of iterations is calculated from lChi2CutFac
1072 fChi2CutFactor = TMath::Max(1.0, lChi2CutFac);
1074 AliInfo(Form("Initial cut factor is %f",fChi2CutFactor));
1075 fIter = 1; // Initializes the iteration process
1079 //_____________________________________________________________________________________________
1080 Double_t AliMillePede2::GetParError(int iPar) const
1082 // return error for parameter iPar
1083 if (fGloSolveStatus==gkInvert) {
1084 double res = fMatCGlo->QueryDiag(iPar);
1085 if (res>=0) return TMath::Sqrt(res);
1091 //_____________________________________________________________________________________________
1092 Int_t AliMillePede2::PrintGlobalParameters() const
1094 /// Print the final results into the logfile
1096 double lGlobalCor =0.;
1099 AliInfo(" Result of fit for global parameters");
1100 AliInfo(" ===================================");
1101 AliInfo(" I initial final differ lastcor error gcor Npnt");
1102 AliInfo("----------------------------------------------------------------------------------------------");
1104 for (int i=0; i<fNGloPar; i++) {
1105 lError = GetParError(i);
1109 if (fGloSolveStatus==gkInvert && TMath::Abs( (dg=fMatCGlo->QueryDiag(i)) *fDiagCGlo[i]) > 0) {
1110 lGlobalCor = TMath::Sqrt(TMath::Abs(1.0-1.0/(dg*fDiagCGlo[i])));
1111 AliInfo(Form("%d\t %.6f\t %.6f\t %.6f\t %.6f\t %.6f\t %.6f\t %6d",
1112 i,fInitPar[i],fInitPar[i]+fDeltaPar[i],fDeltaPar[i],fVecBGlo[i],lError,lGlobalCor,fProcPnt[i]));
1115 AliInfo(Form("%d\t %.6f\t %.6f\t %.6f\t %.6f\t OFF\t OFF\t %6d",i,fInitPar[i],fInitPar[i]+fDeltaPar[i],
1116 fDeltaPar[i],fVecBGlo[i],fProcPnt[i]));
1122 //_____________________________________________________________________________________________
1123 Bool_t AliMillePede2::IsRecordAcceptable() const
1125 // validate record according run lists set by the user
1126 static Long_t prevRunID = kMaxInt;
1127 static Bool_t prevAns = kTRUE;
1128 Long_t runID = fRecord->GetRunID();
1129 if (runID!=prevRunID) {
1132 // is run to be rejected?
1133 if (fRejRunList && (n=fRejRunList->GetSize())) {
1135 for (int i=n;i--;) if (runID == (*fRejRunList)[i]) {
1137 printf("New Run to reject: %ld -> %d\n",runID,prevAns);
1141 else if (fAccRunList && (n=fAccRunList->GetSize())) { // is run specifically selected
1143 for (int i=n;i--;) if (runID == (*fAccRunList)[i]) {prevAns = kTRUE; break;}
1151 //_____________________________________________________________________________________________
1152 void AliMillePede2::SetRejRunList(const UInt_t *runs, Int_t nruns)
1154 // set the list of runs to be rejected
1155 if (fRejRunList) delete fRejRunList;
1157 if (nruns<1 || !runs) return;
1158 fRejRunList = new TArrayL(nruns);
1159 for (int i=0;i<nruns;i++) (*fRejRunList)[i] = runs[i];
1162 //_____________________________________________________________________________________________
1163 void AliMillePede2::SetAccRunList(const UInt_t *runs, Int_t nruns)
1165 // set the list of runs to be selected
1166 if (fAccRunList) delete fAccRunList;
1168 if (nruns<1 || !runs) return;
1169 fAccRunList = new TArrayL(nruns);
1170 for (int i=0;i<nruns;i++) (*fAccRunList)[i] = runs[i];