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>
16 #include "AliMatrixSq.h"
17 #include "AliSymMatrix.h"
18 #include "AliRectMatrix.h"
19 #include "AliMatrixSparse.h"
25 ClassImp(AliMillePede2)
27 Bool_t AliMillePede2::fgInvChol = kTRUE; // Invert global matrix with Cholesky solver
28 Bool_t AliMillePede2::fgWeightSigma = kTRUE; // weight local constraint by module statistics
29 Bool_t AliMillePede2::fgIsMatGloSparse = kFALSE; // use faster dense matrix by default
30 Int_t AliMillePede2::fgMinResCondType = 1; // Jacoby preconditioner by default
31 Double_t AliMillePede2::fgMinResTol = 1.e-11; // default tolerance
32 Int_t AliMillePede2::fgMinResMaxIter = 10000; // default max number of iterations
33 Int_t AliMillePede2::fgIterSol = AliMinResSolve::kSolMinRes; // default iterative solver
34 Int_t AliMillePede2::fgNKrylovV = 240; // default number of Krylov vectors to keep
36 //_____________________________________________________________________________________________
37 AliMillePede2::AliMillePede2()
47 fNLagrangeConstraints(0),
51 fGloSolveStatus(gkFailed),
81 fDataRecFName("/tmp/mp2_data_records.root"),
87 fConstrRecFName("/tmp/mp2_constraints_records.root"),
95 //_____________________________________________________________________________________________
96 AliMillePede2::AliMillePede2(const AliMillePede2& src) :
97 TObject(src),fNLocPar(0),fNGloPar(0),fNGloSize(0),fNLocEquations(0),fIter(0),
98 fMaxIter(10),fNStdDev(3),fNGloConstraints(0),fNLagrangeConstraints(0),
99 fNLocFits(0),fNLocFitsRejected(0),
100 fNGloFix(0),fGloSolveStatus(0),fChi2CutFactor(0),fChi2CutRef(0),fResCutInit(0),
101 fResCut(0),fMinPntValid(1),fNGroupsSet(0),fParamGrID(0),fProcPnt(0),fVecBLoc(0),fDiagCGlo(0),fVecBGlo(0),
102 fInitPar(0),fDeltaPar(0),fSigmaPar(0),fIsLinear(0),fConstrUsed(0),fGlo2CGlo(0),fCGlo2Glo(0),
103 fMatCLoc(0),fMatCGlo(0),fMatCGloLoc(0),fFillIndex(0),fFillValue(0),
104 fDataRecFName(0),fRecord(0),fDataRecFile(0),
105 fTreeData(0),fRecFileStatus(0),fConstrRecFName(0),fTreeConstr(0),fConsRecFile(0),fCurrRecDataID(0),
106 fCurrRecConstrID(0),fLocFitAdd(kTRUE)
109 //_____________________________________________________________________________________________
110 AliMillePede2::~AliMillePede2()
112 CloseDataRecStorage();
113 CloseConsRecStorage();
126 delete[] fConstrUsed;
136 //_____________________________________________________________________________________________
137 Int_t AliMillePede2::InitMille(int nGlo, int nLoc, int lNStdDev,double lResCut, double lResCutInit)
140 if (nLoc>0) fNLocPar = nLoc;
141 if (nGlo>0) fNGloPar = nGlo;
142 if (lResCutInit>0) fResCutInit = lResCutInit;
143 if (lResCut>0) fResCut = lResCut;
144 if (lNStdDev>0) fNStdDev = lNStdDev;
146 fNGloSize = fNGloPar;
150 if (fgIsMatGloSparse) {fMatCGlo = new AliMatrixSparse(fNGloPar); fMatCGlo->SetSymmetric(kTRUE);}
151 else fMatCGlo = new AliSymMatrix(fNGloPar);
153 fFillIndex = new Int_t[fNGloPar];
154 fFillValue = new Double_t[fNGloPar];
156 fMatCLoc = new AliSymMatrix(fNLocPar);
157 fMatCGloLoc = new AliRectMatrix(fNGloPar,fNLocPar);
159 fParamGrID = new Int_t[fNGloPar];
160 fProcPnt = new Int_t[fNGloPar];
161 fVecBLoc = new Double_t[fNLocPar];
162 fDiagCGlo = new Double_t[fNGloPar];
164 fInitPar = new Double_t[fNGloPar];
165 fDeltaPar = new Double_t[fNGloPar];
166 fSigmaPar = new Double_t[fNGloPar];
167 fIsLinear = new Bool_t[fNGloPar];
169 fGlo2CGlo = new Int_t[fNGloPar];
170 fCGlo2Glo = new Int_t[fNGloPar];
173 AliInfo(Form("Failed to allocate the memory for %d global and %d local parameters",fNGloPar,fNLocPar));
177 memset(fVecBLoc ,0,fNLocPar*sizeof(Double_t));
178 memset(fDiagCGlo ,0,fNGloPar*sizeof(Double_t));
179 memset(fInitPar ,0,fNGloPar*sizeof(Double_t));
180 memset(fDeltaPar ,0,fNGloPar*sizeof(Double_t));
181 memset(fSigmaPar ,0,fNGloPar*sizeof(Double_t));
182 memset(fProcPnt ,0,fNGloPar*sizeof(Int_t));
184 for (int i=fNGloPar;i--;) {
185 fGlo2CGlo[i] = fCGlo2Glo[i] = -1;
186 fIsLinear[i] = kTRUE;
193 //_____________________________________________________________________________________________
194 Bool_t AliMillePede2::ImposeDataRecFile(const char* fname)
196 CloseDataRecStorage();
197 SetDataRecFName(fname);
198 return InitDataRecStorage(kTRUE); // open in read mode
201 //_____________________________________________________________________________________________
202 Bool_t AliMillePede2::ImposeConsRecFile(const char* fname)
204 CloseConsRecStorage();
205 SetConsRecFName(fname);
206 return InitConsRecStorage(kTRUE); // open in read mode
209 //_____________________________________________________________________________________________
210 Bool_t AliMillePede2::InitDataRecStorage(Bool_t read)
212 // initialize the buffer for processed measurements records
214 if (fDataRecFile) {AliInfo("Data Records File is already initialized"); return kFALSE;}
216 if (!fRecord) fRecord = new AliMillePedeRecord();
218 fDataRecFile = TFile::Open(GetDataRecFName(),read ? "":"recreate");
219 if (!fDataRecFile) {AliFatal(Form("Failed to initialize data records file %s",GetDataRecFName())); return kFALSE;}
221 AliInfo(Form("File %s used for derivatives records",GetDataRecFName()));
223 fTreeData = (TTree*)fDataRecFile->Get("AliMillePedeRecords_Data");
224 if (!fTreeData) {AliFatal(Form("Did not find data records tree in %s",GetDataRecFName())); return kFALSE;}
225 fTreeData->SetBranchAddress("Record_Data",&fRecord);
226 AliInfo(Form("Found %d derivatives records",fTreeData->GetEntries()));
229 fTreeData = new TTree("AliMillePedeRecords_Data","Data Records for AliMillePede2");
230 fTreeData->Branch("Record_Data","AliMillePedeRecord",&fRecord,32000,99);
233 fRecFileStatus = read ? 1:2;
238 //_____________________________________________________________________________________________
239 Bool_t AliMillePede2::InitConsRecStorage(Bool_t read)
241 // initialize the buffer for processed measurements records
243 if (fConsRecFile) {AliInfo("Constraints Records File is already initialized"); return kFALSE;}
245 if (!fRecord) fRecord = new AliMillePedeRecord();
247 fConsRecFile = TFile::Open(GetConsRecFName(),read ? "":"recreate");
248 if (!fConsRecFile) {AliInfo(Form("Failed to initialize constraints records file %s",GetConsRecFName())); return kFALSE;}
250 AliInfo(Form("File %s used for constraints records",GetConsRecFName()));
252 fTreeConstr = (TTree*)fConsRecFile->Get("AliMillePedeRecords_Constraints");
253 if (!fTreeConstr) {AliInfo(Form("Did not find constraints records tree in %s",GetConsRecFName())); return kFALSE;}
254 fTreeConstr->SetBranchAddress("Record_Constraints",&fRecord);
255 AliInfo(Form("Found %d constraints records",fTreeConstr->GetEntries()));
260 fTreeConstr = new TTree("AliMillePedeRecords_Constraints","Constraints Records for AliMillePede2");
261 fTreeConstr->Branch("Record_Constraints","AliMillePedeRecord",&fRecord,32000,99);
263 fCurrRecConstrID = -1;
268 //_____________________________________________________________________________________________
269 void AliMillePede2::CloseDataRecStorage()
272 if (fDataRecFile->IsWritable()) {
278 fDataRecFile->Close();
286 //_____________________________________________________________________________________________
287 void AliMillePede2::CloseConsRecStorage()
290 if (fConsRecFile->IsWritable()) {
292 fTreeConstr->Write();
296 fConsRecFile->Close();
303 //_____________________________________________________________________________________________
304 Bool_t AliMillePede2::ReadNextRecordData()
306 // read next data record (if any)
307 if (!fTreeData || ++fCurrRecDataID >= fTreeData->GetEntries()) { fCurrRecDataID--; return kFALSE;}
308 fTreeData->GetEntry(fCurrRecDataID);
312 //_____________________________________________________________________________________________
313 Bool_t AliMillePede2::ReadNextRecordConstraint()
315 // read next constraint record (if any)
316 if (!fTreeConstr || ++fCurrRecConstrID >= fTreeConstr->GetEntries()) { fCurrRecConstrID--; return kFALSE;}
317 fTreeConstr->GetEntry(fCurrRecConstrID);
321 //_____________________________________________________________________________________________
322 void AliMillePede2::SetRecordWeight(double wgh)
324 if (fRecFileStatus<2) InitDataRecStorage(); // create a buffer to store the data
325 fRecord->SetWeight(wgh);
328 //_____________________________________________________________________________________________
329 void AliMillePede2::SetLocalEquation(double *dergb, double *derlc, double lMeas, double lSigma)
331 if (fRecFileStatus<2) InitDataRecStorage(); // create a buffer to store the data
333 // write data of single measurement
334 if (lSigma<=0.0) { // If parameter is fixed, then no equation
335 for (int i=fNLocPar; i--;) derlc[i] = 0.0;
336 for (int i=fNGloPar; i--;) dergb[i] = 0.0;
340 fRecord->AddResidual(lMeas);
342 // Retrieve local param interesting indices
343 for (int i=0;i<fNLocPar;i++) if (!IsZero(derlc[i])) {fRecord->AddIndexValue(i,derlc[i]); derlc[i] = 0.0;}
345 fRecord->AddWeight( 1.0/lSigma/lSigma );
347 // Idem for global parameters
348 for (int i=0;i<fNGloPar;i++) if (!IsZero(dergb[i])) {
349 fRecord->AddIndexValue(i,dergb[i]); dergb[i] = 0.0;
350 fRecord->MarkGroup(fParamGrID[i]);
355 //_____________________________________________________________________________________________
356 void AliMillePede2::SetLocalEquation(int *indgb, double *dergb, int ngb, int *indlc,
357 double *derlc,int nlc,double lMeas,double lSigma)
359 // write data of single measurement
360 if (lSigma<=0.0) { // If parameter is fixed, then no equation
361 for (int i=nlc;i--;) derlc[i] = 0.0;
362 for (int i=ngb;i--;) dergb[i] = 0.0;
366 if (fRecFileStatus<2) InitDataRecStorage(); // create a buffer to store the data
368 fRecord->AddResidual(lMeas);
370 // Retrieve local param interesting indices
371 for (int i=0;i<nlc;i++) if (!IsZero(derlc[i])) {fRecord->AddIndexValue(indlc[i],derlc[i]); derlc[i]=0.; indlc[i]=0;}
373 fRecord->AddWeight( 1./lSigma/lSigma );
375 // Idem for global parameters
376 for (int i=0;i<ngb;i++) if (!IsZero(dergb[i])) {fRecord->AddIndexValue(indgb[i],dergb[i]); dergb[i]=0.; indgb[i]=0;}
381 //_____________________________________________________________________________________________
382 void AliMillePede2::SetGlobalConstraint(double *dergb, double val, double sigma)
384 // Define a constraint equation.
385 if (!fConsRecFile || !fConsRecFile->IsWritable()) InitConsRecStorage(); // create a buffer to store the data
388 fRecord->AddResidual(val);
389 fRecord->AddWeight(sigma);
390 for (int i=0; i<fNGloPar; i++) if (!IsZero(dergb[i])) fRecord->AddIndexValue(i,dergb[i]);
392 if (IsZero(sigma)) fNLagrangeConstraints++;
393 // printf("NewConstraint:\n"); fRecord->Print(); //RRR
394 SaveRecordConstraint();
398 //_____________________________________________________________________________________________
399 void AliMillePede2::SetGlobalConstraint(const int *indgb, double *dergb, int ngb, double val,double sigma)
401 // Define a constraint equation.
402 if (!fConsRecFile || !fConsRecFile->IsWritable()) InitConsRecStorage(); // create a buffer to store the data
404 fRecord->AddResidual(val);
405 fRecord->AddWeight(sigma); // dummy
406 for (int i=0; i<ngb; i++) if (!IsZero(dergb[i])) fRecord->AddIndexValue(indgb[i],dergb[i]);
408 if (IsZero(sigma)) fNLagrangeConstraints++;
409 SaveRecordConstraint();
413 //_____________________________________________________________________________________________
414 Int_t AliMillePede2::LocalFit(double *localParams)
417 Perform local parameters fit once all the local equations have been set
418 -----------------------------------------------------------
419 localParams = (if !=0) will contain the fitted track parameters and
422 static int nrefSize = 0;
423 // static TArrayI refLoc,refGlo,nrefLoc,nrefGlo;
424 static Int_t *refLoc=0,*refGlo=0,*nrefLoc=0,*nrefGlo=0;
427 AliSymMatrix &matCLoc = *fMatCLoc;
428 AliMatrixSq &matCGlo = *fMatCGlo;
429 AliRectMatrix &matCGloLoc = *fMatCGloLoc;
431 memset(fVecBLoc,0,fNLocPar*sizeof(double));
435 int recSz = fRecord->GetSize();
437 while(cnt<recSz) { // Transfer the measurement records to matrices
439 // extract addresses of residual, weight and pointers on local and global derivatives for each point
440 if (nrefSize<=nPoints) {
442 nrefSize = 2*(nPoints+1);
443 tmpA = refLoc; refLoc = new Int_t[nrefSize]; if (tmpA) memcpy(refLoc,tmpA,nPoints*sizeof(int));
444 tmpA = refGlo; refGlo = new Int_t[nrefSize]; if (tmpA) memcpy(refGlo,tmpA,nPoints*sizeof(int));
445 tmpA = nrefLoc; nrefLoc = new Int_t[nrefSize]; if (tmpA) memcpy(nrefLoc,tmpA,nPoints*sizeof(int));
446 tmpA = nrefGlo; nrefGlo = new Int_t[nrefSize]; if (tmpA) memcpy(nrefGlo,tmpA,nPoints*sizeof(int));
449 refLoc[nPoints] = ++cnt;
451 while(!fRecord->IsWeight(cnt)) {nLoc++; cnt++;}
452 nrefLoc[nPoints] = nLoc;
454 refGlo[nPoints] = ++cnt;
456 while(!fRecord->IsResidual(cnt) && cnt<recSz) {nGlo++; cnt++;}
457 nrefGlo[nPoints] = nGlo;
463 double gloWgh = fRecord->GetWeight(); // global weight for this set
464 Int_t maxLocUsed = 0;
466 for (int ip=nPoints;ip--;) { // Transfer the measurement records to matrices
467 double resid = fRecord->GetValue( refLoc[ip]-1 );
468 double weight = fRecord->GetValue( refGlo[ip]-1 )*gloWgh;
469 double *derLoc = fRecord->GetValue()+refLoc[ip];
470 double *derGlo = fRecord->GetValue()+refGlo[ip];
471 int *indLoc = fRecord->GetIndex()+refLoc[ip];
472 int *indGlo = fRecord->GetIndex()+refGlo[ip];
474 for (int i=nrefGlo[ip];i--;) { // suppress the global part (only relevant with iterations)
475 int iID = indGlo[i]; // Global param indice
476 if (fSigmaPar[iID]<=0.) continue; // fixed parameter RRRCheck
477 if (fIsLinear[iID]) resid -= derGlo[i]*(fInitPar[iID]+fDeltaPar[iID]); // linear parameter
478 else resid -= derGlo[i]*fDeltaPar[iID]; // nonlinear parameter
481 // Symmetric matrix, don't bother j>i coeffs
482 for (int i=nrefLoc[ip];i--;) { // Fill local matrix and vector
483 fVecBLoc[ indLoc[i] ] += weight*resid*derLoc[i];
484 if (indLoc[i]>maxLocUsed) maxLocUsed = indLoc[i];
485 for (int j=i+1;j--;) matCLoc(indLoc[i] ,indLoc[j]) += weight*derLoc[i]*derLoc[j];
488 } // end of the transfer of the measurement record to matrices
490 matCLoc.SetSizeUsed(++maxLocUsed); // data with B=0 may use less than declared nLocals
492 // first try to solve by faster Cholesky decomposition, then by Gaussian elimination
493 if (!matCLoc.SolveChol(fVecBLoc,kTRUE)) {
494 AliInfo("Failed to solve locals by Cholesky, trying Gaussian Elimination");
495 if (!matCLoc.SolveSpmInv(fVecBLoc,kTRUE)) {
496 AliInfo("Failed to solve locals by Gaussian Elimination, skip...");
498 return 0; // failed to solve
502 // If requested, store the track params and errors
503 // printf("locfit: "); for (int i=0;i<fNLocPar;i++) printf("%+e |",fVecBLoc[i]); printf("\n");
504 if (localParams) for (int i=maxLocUsed; i--;) {
505 localParams[2*i] = fVecBLoc[i];
506 localParams[2*i+1] = TMath::Sqrt(TMath::Abs(matCLoc.QueryDiag(i)));
512 for (int ip=nPoints;ip--;) { // Calculate residuals
513 double resid = fRecord->GetValue( refLoc[ip]-1 );
514 double weight = fRecord->GetValue( refGlo[ip]-1 )*gloWgh;
515 double *derLoc = fRecord->GetValue()+refLoc[ip];
516 double *derGlo = fRecord->GetValue()+refGlo[ip];
517 int *indLoc = fRecord->GetIndex()+refLoc[ip];
518 int *indGlo = fRecord->GetIndex()+refGlo[ip];
520 // Suppress local and global contribution in residuals;
521 for (int i=nrefLoc[ip];i--;) resid -= derLoc[i]*fVecBLoc[ indLoc[i] ]; // local part
523 for (int i=nrefGlo[ip];i--;) { // global part
525 if ( fSigmaPar[iID] <= 0.) continue; // fixed parameter RRRCheck
526 if (fIsLinear[iID]) resid -= derGlo[i]*(fInitPar[iID]+fDeltaPar[iID]); // linear parameter
527 else resid -= derGlo[i]*fDeltaPar[iID]; // nonlinear parameter
530 // reject the track if the residual is too large (outlier)
531 double absres = TMath::Abs(resid);
532 if ( (absres >= fResCutInit && fIter ==1 ) ||
533 (absres >= fResCut && fIter > 1)) {
534 if (fLocFitAdd) fNLocFitsRejected++;
535 // printf("reject res %5ld %+e\n",fCurrRecDataID,resid);
539 lChi2 += weight*resid*resid ; // total chi^2
540 nEq++; // number of equations
541 } // end of Calculate residuals
544 int nDoF = nEq-maxLocUsed;
545 lChi2 = (nDoF>0) ? lChi2/nDoF : 0; // Chi^2/dof
547 if (fNStdDev != 0 && nDoF>0 && lChi2 > Chi2DoFLim(fNStdDev,nDoF)*fChi2CutFactor) { // check final chi2
548 if (fLocFitAdd) fNLocFitsRejected++;
549 // printf("reject chi2 %5ld: %+e\n",fCurrRecDataID, lChi2);
555 fNLocEquations += nEq;
559 fNLocEquations -= nEq;
562 // local operations are finished, track is accepted
563 // We now update the global parameters (other matrices)
567 for (int ip=nPoints;ip--;) { // Update matrices
568 double resid = fRecord->GetValue( refLoc[ip]-1 );
569 double weight = fRecord->GetValue( refGlo[ip]-1 )*gloWgh;
570 double *derLoc = fRecord->GetValue()+refLoc[ip];
571 double *derGlo = fRecord->GetValue()+refGlo[ip];
572 int *indLoc = fRecord->GetIndex()+refLoc[ip];
573 int *indGlo = fRecord->GetIndex()+refGlo[ip];
575 for (int i=nrefGlo[ip];i--;) { // suppress the global part
576 int iID = indGlo[i]; // Global param indice
577 if ( fSigmaPar[iID] <= 0.) continue; // fixed parameter RRRCheck
578 if (fIsLinear[iID]) resid -= derGlo[i]*(fInitPar[iID]+fDeltaPar[iID]); // linear parameter
579 else resid -= derGlo[i]*fDeltaPar[iID]; // nonlinear parameter
582 for (int ig=nrefGlo[ip];ig--;) {
583 int iIDg = indGlo[ig]; // Global param indice (the matrix line)
584 if ( fSigmaPar[iIDg] <= 0.) continue; // fixed parameter RRRCheck
585 if (fLocFitAdd) fVecBGlo[ iIDg ] += weight*resid*derGlo[ig]; //!!!
586 else fVecBGlo[ iIDg ] -= weight*resid*derGlo[ig]; //!!!
588 // First of all, the global/global terms (exactly like local matrix)
590 for (int jg=ig+1;jg--;) { // matCGlo is symmetric by construction
591 int jIDg = indGlo[jg];
592 if ( fSigmaPar[jIDg] <= 0.) continue; // fixed parameter RRRCheck
593 if ( !IsZero(vl = weight*derGlo[ig]*derGlo[jg]) ) {
594 fFillIndex[nfill] = jIDg;
595 fFillValue[nfill++] = fLocFitAdd ? vl:-vl;
598 if (nfill) matCGlo.AddToRow(iIDg,fFillValue,fFillIndex,nfill);
600 // Now we have also rectangular matrices containing global/local terms.
601 int iCIDg = fGlo2CGlo[iIDg]; // compressed Index of index
603 Double_t *rowGL = matCGloLoc(nGloInFit);
604 for (int k=maxLocUsed;k--;) rowGL[k] = 0.0; // reset the row
605 iCIDg = fGlo2CGlo[iIDg] = nGloInFit;
606 fCGlo2Glo[nGloInFit++] = iIDg;
609 Double_t *rowGLIDg = matCGloLoc(iCIDg);
610 for (int il=nrefLoc[ip];il--;) rowGLIDg[ indLoc[il] ] += weight*derGlo[ig]*derLoc[il];
611 fProcPnt[iIDg] += fLocFitAdd ? 1:-1; // update counter
614 } // end of Update matrices
616 // calculate fMatCGlo -= fMatCGloLoc * fMatCLoc * fMatCGloLoc^T
617 // and fVecBGlo -= fMatCGloLoc * fVecBLoc
619 //-------------------------------------------------------------- >>>
621 for (int iCIDg=0; iCIDg<nGloInFit; iCIDg++) {
622 int iIDg = fCGlo2Glo[iCIDg];
625 Double_t *rowGLIDg = matCGloLoc(iCIDg);
626 for (int kl=0;kl<maxLocUsed;kl++) if (rowGLIDg[kl]) vl += rowGLIDg[kl]*fVecBLoc[kl];
627 if (!IsZero(vl)) fVecBGlo[iIDg] -= fLocFitAdd ? vl : -vl;
630 for (int jCIDg=0;jCIDg<=iCIDg; jCIDg++) {
631 int jIDg = fCGlo2Glo[jCIDg];
634 Double_t *rowGLJDg = matCGloLoc(jCIDg);
635 for (int kl=0;kl<maxLocUsed;kl++) {
637 if ( (!IsZero(vll=rowGLIDg[kl]*rowGLJDg[kl])) ) vl += matCLoc.QueryDiag(kl)*vll;
640 for (int ll=0;ll<kl;ll++) {
641 if ( !IsZero(vll=rowGLIDg[kl]*rowGLJDg[ll]) ) vl += matCLoc(kl,ll)*vll;
642 if ( !IsZero(vll=rowGLIDg[ll]*rowGLJDg[kl]) ) vl += matCLoc(kl,ll)*vll;
646 fFillIndex[nfill] = jIDg;
647 fFillValue[nfill++] = fLocFitAdd ? -vl : vl;
650 if (nfill) matCGlo.AddToRow(iIDg,fFillValue,fFillIndex,nfill);
653 // reset compressed index array
655 for (int i=nGloInFit;i--;) {
656 fGlo2CGlo[ fCGlo2Glo[i] ] = -1;
660 //---------------------------------------------------- <<<
664 //_____________________________________________________________________________________________
665 Int_t AliMillePede2::GlobalFit(Double_t *par, Double_t *error, Double_t *pull)
667 // performs a requested number of global iterations
670 TStopwatch sw; sw.Start();
673 AliInfo("Starting Global fit.");
674 while (fIter<=fMaxIter) {
676 res = GlobalFitIteration();
679 if (!IsZero(fChi2CutFactor-fChi2CutRef)) {
680 fChi2CutFactor = TMath::Sqrt(fChi2CutFactor);
681 if (fChi2CutFactor < 1.2*fChi2CutRef) {
682 fChi2CutFactor = fChi2CutRef;
683 fIter = fMaxIter - 1; // Last iteration
690 AliInfo(Form("Global fit %s, CPU time: %.1f",res ? "Converged":"Failed",sw.CpuTime()));
693 if (par) for (int i=fNGloPar;i--;) par[i] = fInitPar[i]+fDeltaPar[i];
695 if (fGloSolveStatus==gkInvert) { // errors on params are available
696 if (error) for (int i=fNGloPar;i--;) error[i] = fProcPnt[i]>0 ? TMath::Sqrt(TMath::Abs(fMatCGlo->QueryDiag(i))) : 0.;
697 if (pull) for (int i=fNGloPar;i--;) pull[i] = fProcPnt[i]>0 && (fSigmaPar[i]*fSigmaPar[i]-fMatCGlo->QueryDiag(i))>0. && fSigmaPar[i]>0
698 ? fDeltaPar[i]/TMath::Sqrt(fSigmaPar[i]*fSigmaPar[i]-fMatCGlo->QueryDiag(i)) : 0;
704 //_____________________________________________________________________________________________
705 Int_t AliMillePede2::GlobalFitIteration()
707 // perform global parameters fit once all the local equations have been fitted
709 AliInfo(Form("Global Fit Iteration#%2d (Local Fit Chi^2 cut factor: %.2f)",fIter,fChi2CutFactor));
711 if (!fNGloPar || !fTreeData) {
712 AliInfo("No data was stored, aborting iteration");
720 fConstrUsed = new Bool_t[fNGloConstraints];
721 memset(fConstrUsed,0,fNGloConstraints*sizeof(Bool_t));
723 // Reset all info specific for this step
724 AliMatrixSq& matCGlo = *fMatCGlo;
726 memset(fProcPnt,0,fNGloPar*sizeof(Int_t));
728 fNGloConstraints = fTreeConstr ? fTreeConstr->GetEntries() : 0;
730 // count number of Lagrange constraints: they need new row/cols to be added
731 fNLagrangeConstraints = 0;
732 for (int i=0; i<fNGloConstraints; i++) {
733 ReadRecordConstraint(i);
734 if ( IsZero(fRecord->GetValue(1)) ) fNLagrangeConstraints++; // exact constraint (no error) -> Lagrange multiplier
737 // if needed, readjust the size of the global vector (for matrices this is done automatically)
738 if (!fVecBGlo || fNGloSize!=fNGloPar+fNLagrangeConstraints) {
739 delete[] fVecBGlo; // in case some constraint was added between the two manual iterations
740 fNGloSize = fNGloPar+fNLagrangeConstraints;
741 fVecBGlo = new Double_t[fNGloSize];
743 memset(fVecBGlo,0,fNGloSize*sizeof(double));
746 fNLocFitsRejected = 0;
749 // Process data records and build the matrices
750 Long_t ndr = fTreeData->GetEntries();
751 AliInfo(Form("Building the Global matrix from %d data records",ndr));
754 TStopwatch swt; swt.Start();
755 fLocFitAdd = kTRUE; // add contributions of matching tracks
756 for (Long_t i=0;i<ndr;i++) {
759 if ( (i%int(0.2*ndr)) == 0) printf("%.1f%% of local fits done\n", double(100.*i)/ndr);
762 printf("%ld local fits done: ", ndr);
767 // ---------------------- Reject parameters with low statistics ------------>>
769 if (fMinPntValid>1 && fNGroupsSet) {
771 printf("Checking parameters with statistics < %d\n",fMinPntValid);
774 // 1) build the list of parameters to fix
775 Int_t fixArrSize = 10;
776 Int_t nFixedGroups = 0;
777 TArrayI fixGroups(fixArrSize);
781 double oldMin = 1.e20;
782 double oldMax =-1.e20;
784 for (int i=fNGloPar;i--;) { // // Reset row and column of fixed params and add 1/sig^2 to free ones
785 int grID = fParamGrID[i];
786 if (grID<0) continue; // not in the group
788 if (grID!=grIDold) { // starting new group
789 if (grIDold>=0) { // decide if the group has enough statistics
790 if (oldMin<fMinPntValid && oldMax<2*fMinPntValid) { // suppress group
791 for (int iold=oldStart;iold>i;iold--) fProcPnt[iold] = 0;
792 Bool_t fnd = kFALSE; // check if the group is already accounted
793 for (int j=nFixedGroups;j--;) if (fixGroups[j]==grIDold) {fnd=kTRUE; break;}
795 if (nFixedGroups>=fixArrSize) {fixArrSize*=2; fixGroups.Set(fixArrSize);}
796 fixGroups[nFixedGroups++] = grIDold; // add group to fix
800 grIDold = grID; // mark the start of the new group
805 if (oldMin>fProcPnt[i]) oldMin = fProcPnt[i];
806 if (oldMax<fProcPnt[i]) oldMax = fProcPnt[i];
809 // extra check for the last group
810 if (grIDold>=0 && oldMin<fMinPntValid && oldMax<2*fMinPntValid) { // suppress group
811 for (int iold=oldStart;iold--;) fProcPnt[iold] = 0;
812 Bool_t fnd = kFALSE; // check if the group is already accounted
813 for (int j=nFixedGroups;j--;) if (fixGroups[j]==grIDold) {fnd=kTRUE; break;}
815 if (nFixedGroups>=fixArrSize) {fixArrSize*=2; fixGroups.Set(fixArrSize);}
816 fixGroups[nFixedGroups++] = grIDold; // add group to fix
820 // 2) loop over records and add contributions of fixed groups with negative sign
823 for (Long_t i=0;i<ndr;i++) {
825 Bool_t suppr = kFALSE;
826 for (int ifx=nFixedGroups;ifx--;)if (fRecord->IsGroupPresent(fixGroups[ifx])) suppr = kTRUE;
827 if (suppr) LocalFit();
832 printf("Suppressed contributions of groups with NPoints<%d :\n",fMinPntValid);
833 for (int i=0;i<nFixedGroups;i++) printf("%d ",fixGroups[i]); printf("\n");
838 // ---------------------- Reject parameters with low statistics ------------<<
840 // add large number to diagonal of fixed params
842 for (int i=fNGloPar;i--;) { // // Reset row and column of fixed params and add 1/sig^2 to free ones
843 // printf("#%3d : Nproc : %5d grp: %d\n",i,fProcPnt[i],fParamGrID[i]);
847 matCGlo.DiagElem(i) = 1.;//float(fNLocEquations*fNLocEquations);
849 else matCGlo.DiagElem(i) += (fgWeightSigma ? fProcPnt[i] : 1.)/(fSigmaPar[i]*fSigmaPar[i]);
852 for (int i=fNGloPar;i--;) fDiagCGlo[i] = matCGlo.QueryDiag(i); // save the diagonal elements
854 // add constraint equations
855 int nVar = fNGloPar; // Current size of global matrix
856 for (int i=0; i<fNGloConstraints; i++) {
857 ReadRecordConstraint(i);
858 double val = fRecord->GetValue(0);
859 double sig = fRecord->GetValue(1);
860 int *indV = fRecord->GetIndex()+2;
861 double *der = fRecord->GetValue()+2;
862 int csize = fRecord->GetSize()-2;
864 // check if after suppression of fixed variables there are non-0 derivatives
865 // and determine the max statistics of involved params
868 for (int j=csize;j--;) {
869 if (fProcPnt[indV[j]]<1) nSuppressed++;
871 maxStat = TMath::Max(maxStat,fProcPnt[indV[j]]);
875 if (nSuppressed==csize) {
876 // AliInfo(Form("Neglecting constraint %d of %d derivatives since no free parameters left",i,csize));
878 // was this constraint ever created ?
879 if ( sig==0 && fConstrUsed[i] ) { // this is needed only for constraints with Lagrange multiplier
880 // to avoid empty row impose dummy constraint on "Lagrange multiplier"
881 matCGlo.DiagElem(nVar) = 1.;
882 fVecBGlo[nVar++] = 0;
887 // account for already accumulated corrections
888 for (int j=csize; j--;) val -= der[j]*(fInitPar[ indV[j] ]+fDeltaPar[ indV[j] ]);
890 if (sig > 0) { // this is a gaussian constriant: no Lagrange multipliers are added
892 double sig2i = (fgWeightSigma ? TMath::Sqrt(maxStat) : 1.)/sig/sig;
893 for (int ir=0;ir<csize;ir++) {
895 for (int ic=0;ic<=ir;ic++) { // matrix is symmetric
897 double vl = der[ir]*der[ic]*sig2i;
898 if (!IsZero(vl)) matCGlo(iID,jID) += vl;
900 fVecBGlo[iID] += val*der[ir]*sig2i;
903 else { // this is exact constriant: Lagrange multipliers must be added
904 for (int j=csize; j--;) {
906 if (fProcPnt[jID]<1) continue; // this parameter was fixed, don't put it into constraint
907 matCGlo(nVar,jID) = float(fNLocEquations)*der[j]; // fMatCGlo is symmetric, only lower triangle is filled
910 if (matCGlo.QueryDiag(nVar)) matCGlo.DiagElem(nVar) = 0.0;
911 fVecBGlo[nVar++] = float(fNLocEquations)*val; //RS ? should we use here fNLocFits ?
912 fConstrUsed[i] = kTRUE;
916 AliInfo(Form("Obtained %-7ld equations from %-7ld records (%-7ld rejected). Fixed %-4d globals",
917 fNLocEquations,fNLocFits,fNLocFitsRejected,fNGloFix));
921 fGloSolveStatus = SolveGlobalMatEq(); // obtain solution for this step
923 printf("Solve %d |",fIter); sws.Print();
926 AliInfo(Form("Iteration#%2d %s. CPU time: %.1f",fIter,fGloSolveStatus==gkFailed ? "Failed":"Converged",sw.CpuTime()));
927 if (fGloSolveStatus==gkFailed) return 0;
929 for (int i=fNGloPar;i--;) fDeltaPar[i] += fVecBGlo[i]; // Update global parameters values (for iterations)
931 // PrintGlobalParameters();
935 //_____________________________________________________________________________________________
936 Int_t AliMillePede2::SolveGlobalMatEq()
939 // solve global matrix equation MatCGlob*X=VecBGlo and store the result in the VecBGlo
942 printf("GlobalMatrix\n");
945 for (int i=0;i<fNGloPar;i++) printf("%d %+e\n",i,fVecBGlo[i]);
948 if (!fgIsMatGloSparse) {
950 if (fNLagrangeConstraints==0) { // pos-def systems are faster to solve by Cholesky
951 if ( ((AliSymMatrix*)fMatCGlo)->SolveChol(fVecBGlo, fgInvChol) ) return fgInvChol ? gkInvert:gkNoInversion;
952 else AliInfo("Solution of Global Dense System by Cholesky failed, trying Gaussian Elimiation");
955 if (((AliSymMatrix*)fMatCGlo)->SolveSpmInv(fVecBGlo, kTRUE)) return gkInvert;
956 else AliInfo("Solution of Global Dense System by Gaussian Elimiation failed, trying iterative methods");
958 // try to solve by minres
959 TVectorD sol(fNGloSize);
961 AliMinResSolve *slv = new AliMinResSolve(fMatCGlo,fVecBGlo);
962 if (!slv) return gkFailed;
965 if (fgIterSol == AliMinResSolve::kSolMinRes)
966 res = slv->SolveMinRes(sol,fgMinResCondType,fgMinResMaxIter,fgMinResTol);
967 else if (fgIterSol == AliMinResSolve::kSolFGMRes)
968 res = slv->SolveFGMRES(sol,fgMinResCondType,fgMinResMaxIter,fgMinResTol,fgNKrylovV);
970 AliInfo(Form("Undefined Iteritive Solver ID=%d, only %d are defined",fgIterSol,AliMinResSolve::kNSolvers));
972 if (!res) return gkFailed;
973 for (int i=fNGloSize;i--;) fVecBGlo[i] = sol[i];
974 return gkNoInversion;
978 //_____________________________________________________________________________________________
979 Float_t AliMillePede2::Chi2DoFLim(int nSig, int nDoF) const
981 /// return the limit in chi^2/nd for n sigmas stdev authorized
982 // Only n=1, 2, and 3 are expected in input
984 float sn[3] = {0.47523, 1.690140, 2.782170};
985 float table[3][30] = {{1.0000, 1.1479, 1.1753, 1.1798, 1.1775, 1.1730, 1.1680, 1.1630,
986 1.1581, 1.1536, 1.1493, 1.1454, 1.1417, 1.1383, 1.1351, 1.1321,
987 1.1293, 1.1266, 1.1242, 1.1218, 1.1196, 1.1175, 1.1155, 1.1136,
988 1.1119, 1.1101, 1.1085, 1.1070, 1.1055, 1.1040},
989 {4.0000, 3.0900, 2.6750, 2.4290, 2.2628, 2.1415, 2.0481, 1.9736,
990 1.9124, 1.8610, 1.8171, 1.7791, 1.7457, 1.7161, 1.6897, 1.6658,
991 1.6442, 1.6246, 1.6065, 1.5899, 1.5745, 1.5603, 1.5470, 1.5346,
992 1.5230, 1.5120, 1.5017, 1.4920, 1.4829, 1.4742},
993 {9.0000, 5.9146, 4.7184, 4.0628, 3.6410, 3.3436, 3.1209, 2.9468,
994 2.8063, 2.6902, 2.5922, 2.5082, 2.4352, 2.3711, 2.3143, 2.2635,
995 2.2178, 2.1764, 2.1386, 2.1040, 2.0722, 2.0428, 2.0155, 1.9901,
996 1.9665, 1.9443, 1.9235, 1.9040, 1.8855, 1.8681}};
1002 lNSig = TMath::Max(1,TMath::Min(nSig,3));
1005 return table[lNSig-1][nDoF-1];
1007 else { // approximation
1008 return ((sn[lNSig-1]+TMath::Sqrt(float(2*nDoF-3)))*
1009 (sn[lNSig-1]+TMath::Sqrt(float(2*nDoF-3))))/float(2*nDoF-2);
1014 //_____________________________________________________________________________________________
1015 Int_t AliMillePede2::SetIterations(double lChi2CutFac)
1017 // Number of iterations is calculated from lChi2CutFac
1018 fChi2CutFactor = TMath::Max(1.0, lChi2CutFac);
1020 AliInfo(Form("Initial cut factor is %f",fChi2CutFactor));
1021 fIter = 1; // Initializes the iteration process
1025 //_____________________________________________________________________________________________
1026 Double_t AliMillePede2::GetParError(int iPar) const
1028 // return error for parameter iPar
1029 if (fGloSolveStatus==gkInvert) {
1030 double res = fMatCGlo->QueryDiag(iPar);
1031 if (res>=0) return TMath::Sqrt(res);
1037 //_____________________________________________________________________________________________
1038 Int_t AliMillePede2::PrintGlobalParameters() const
1040 /// Print the final results into the logfile
1042 double lGlobalCor =0.;
1045 AliInfo(" Result of fit for global parameters");
1046 AliInfo(" ===================================");
1047 AliInfo(" I initial final differ lastcor error gcor Npnt");
1048 AliInfo("----------------------------------------------------------------------------------------------");
1050 for (int i=0; i<fNGloPar; i++) {
1051 lError = GetParError(i);
1055 if (fGloSolveStatus==gkInvert && TMath::Abs( (dg=fMatCGlo->QueryDiag(i)) *fDiagCGlo[i]) > 0) {
1056 lGlobalCor = TMath::Sqrt(TMath::Abs(1.0-1.0/(dg*fDiagCGlo[i])));
1057 AliInfo(Form("%d\t %.6f\t %.6f\t %.6f\t %.6f\t %.6f\t %.6f\t %6d",
1058 i,fInitPar[i],fInitPar[i]+fDeltaPar[i],fDeltaPar[i],fVecBGlo[i],lError,lGlobalCor,fProcPnt[i]));
1061 AliInfo(Form("%d\t %.6f\t %.6f\t %.6f\t %.6f\t OFF\t OFF\t %6d",i,fInitPar[i],fInitPar[i]+fDeltaPar[i],
1062 fDeltaPar[i],fVecBGlo[i],fProcPnt[i]));