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>
19 #include "AliMatrixSq.h"
20 #include "AliSymMatrix.h"
21 #include "AliRectMatrix.h"
22 #include "AliMatrixSparse.h"
26 #include <sys/types.h>
34 ClassImp(AliMillePede2)
36 Bool_t AliMillePede2::fgInvChol = kTRUE; // Invert global matrix with Cholesky solver
37 Bool_t AliMillePede2::fgWeightSigma = kTRUE; // weight local constraint by module statistics
38 Bool_t AliMillePede2::fgIsMatGloSparse = kFALSE; // use faster dense matrix by default
39 Int_t AliMillePede2::fgMinResCondType = 1; // Jacoby preconditioner by default
40 Double_t AliMillePede2::fgMinResTol = 1.e-11; // default tolerance
41 Int_t AliMillePede2::fgMinResMaxIter = 10000; // default max number of iterations
42 Int_t AliMillePede2::fgIterSol = AliMinResSolve::kSolMinRes; // default iterative solver
43 Int_t AliMillePede2::fgNKrylovV = 240; // default number of Krylov vectors to keep
45 //_____________________________________________________________________________________________
46 AliMillePede2::AliMillePede2()
56 fNLagrangeConstraints(0),
60 fGloSolveStatus(gkFailed),
90 fDataRecFName("/tmp/mp2_data_records.root"),
96 fConstrRecFName("/tmp/mp2_constraints_records.root"),
108 //_____________________________________________________________________________________________
109 AliMillePede2::AliMillePede2(const AliMillePede2& src) :
110 TObject(src),fNLocPar(0),fNGloPar(0),fNGloSize(0),fNLocEquations(0),fIter(0),
111 fMaxIter(10),fNStdDev(3),fNGloConstraints(0),fNLagrangeConstraints(0),
112 fNLocFits(0),fNLocFitsRejected(0),
113 fNGloFix(0),fGloSolveStatus(0),fChi2CutFactor(0),fChi2CutRef(0),fResCutInit(0),
114 fResCut(0),fMinPntValid(1),fNGroupsSet(0),fParamGrID(0),fProcPnt(0),fVecBLoc(0),fDiagCGlo(0),fVecBGlo(0),
115 fInitPar(0),fDeltaPar(0),fSigmaPar(0),fIsLinear(0),fConstrUsed(0),fGlo2CGlo(0),fCGlo2Glo(0),
116 fMatCLoc(0),fMatCGlo(0),fMatCGloLoc(0),fFillIndex(0),fFillValue(0),
117 fDataRecFName(0),fRecord(0),fDataRecFile(0),
118 fTreeData(0),fRecFileStatus(0),fConstrRecFName(0),fTreeConstr(0),fConsRecFile(0),fCurrRecDataID(0),
119 fCurrRecConstrID(0),fLocFitAdd(kTRUE),
126 //_____________________________________________________________________________________________
127 AliMillePede2::~AliMillePede2()
129 CloseDataRecStorage();
130 CloseConsRecStorage();
143 delete[] fConstrUsed;
155 //_____________________________________________________________________________________________
156 Int_t AliMillePede2::InitMille(int nGlo, int nLoc, int lNStdDev,double lResCut, double lResCutInit)
159 if (nLoc>0) fNLocPar = nLoc;
160 if (nGlo>0) fNGloPar = nGlo;
161 if (lResCutInit>0) fResCutInit = lResCutInit;
162 if (lResCut>0) fResCut = lResCut;
163 if (lNStdDev>0) fNStdDev = lNStdDev;
165 fNGloSize = fNGloPar;
169 if (fgIsMatGloSparse) {fMatCGlo = new AliMatrixSparse(fNGloPar); fMatCGlo->SetSymmetric(kTRUE);}
170 else fMatCGlo = new AliSymMatrix(fNGloPar);
172 fFillIndex = new Int_t[fNGloPar];
173 fFillValue = new Double_t[fNGloPar];
175 fMatCLoc = new AliSymMatrix(fNLocPar);
176 fMatCGloLoc = new AliRectMatrix(fNGloPar,fNLocPar);
178 fParamGrID = new Int_t[fNGloPar];
179 fProcPnt = new Int_t[fNGloPar];
180 fVecBLoc = new Double_t[fNLocPar];
181 fDiagCGlo = new Double_t[fNGloPar];
183 fInitPar = new Double_t[fNGloPar];
184 fDeltaPar = new Double_t[fNGloPar];
185 fSigmaPar = new Double_t[fNGloPar];
186 fIsLinear = new Bool_t[fNGloPar];
188 fGlo2CGlo = new Int_t[fNGloPar];
189 fCGlo2Glo = new Int_t[fNGloPar];
192 AliInfo(Form("Failed to allocate the memory for %d global and %d local parameters",fNGloPar,fNLocPar));
196 memset(fVecBLoc ,0,fNLocPar*sizeof(Double_t));
197 memset(fDiagCGlo ,0,fNGloPar*sizeof(Double_t));
198 memset(fInitPar ,0,fNGloPar*sizeof(Double_t));
199 memset(fDeltaPar ,0,fNGloPar*sizeof(Double_t));
200 memset(fSigmaPar ,0,fNGloPar*sizeof(Double_t));
201 memset(fProcPnt ,0,fNGloPar*sizeof(Int_t));
203 for (int i=fNGloPar;i--;) {
204 fGlo2CGlo[i] = fCGlo2Glo[i] = -1;
205 fIsLinear[i] = kTRUE;
212 //_____________________________________________________________________________________________
213 Bool_t AliMillePede2::ImposeDataRecFile(const char* fname)
215 CloseDataRecStorage();
216 SetDataRecFName(fname);
217 return InitDataRecStorage(kTRUE); // open in read mode
220 //_____________________________________________________________________________________________
221 Bool_t AliMillePede2::ImposeConsRecFile(const char* fname)
223 CloseConsRecStorage();
224 SetConsRecFName(fname);
225 return InitConsRecStorage(kTRUE); // open in read mode
228 //_____________________________________________________________________________________________
229 Bool_t AliMillePede2::InitDataRecStorage(Bool_t read)
231 // initialize the buffer for processed measurements records
233 if (fTreeData) {AliInfo("Data Records File is already initialized"); return kFALSE;}
235 if (!fRecord) fRecord = new AliMillePedeRecord();
237 if (!read) { // write mode: cannot use chain
238 fDataRecFile = TFile::Open(GetDataRecFName(),"recreate");
239 if (!fDataRecFile) {AliFatal(Form("Failed to initialize data records file %s",GetDataRecFName())); return kFALSE;}
240 AliInfo(Form("File %s used for derivatives records",GetDataRecFName()));
241 fTreeData = new TTree("AliMillePedeRecords_Data","Data Records for AliMillePede2");
242 fTreeData->Branch("Record_Data","AliMillePedeRecord",&fRecord,32000,99);
245 TChain* ch = new TChain("AliMillePedeRecords_Data");
247 if (fDataRecFName.EndsWith(".root")) ch->AddFile(fDataRecFName);
248 else { // assume text file with list of filenames
250 ifstream inpf(fDataRecFName.Data());
251 if (!inpf.good()) {AliInfo(Form("Failed on input records list %s\n",fDataRecFName.Data())); return kFALSE;}
254 recfName.ReadLine(inpf);
255 while ( !(recfName.ReadLine(inpf)).eof() ) {
256 recfName = recfName.Strip(TString::kBoth,' ');
257 if (recfName.BeginsWith("//") || recfName.BeginsWith("#") || !recfName.EndsWith(".root")) { // comment
258 AliInfo(Form("Skip %s\n",recfName.Data()));
262 recfName = recfName.Strip(TString::kBoth,',');
263 recfName = recfName.Strip(TString::kBoth,'"');
264 gSystem->ExpandPathName(recfName);
265 printf("Adding %s\n",recfName.Data());
266 ch->AddFile(recfName.Data());
270 Long64_t nent = ch->GetEntries();
271 if (nent<1) { AliInfo("Obtained chain is empty"); return kFALSE;}
273 fTreeData->SetBranchAddress("Record_Data",&fRecord);
274 AliInfo(Form("Found %lld derivatives records",nent));
277 fRecFileStatus = read ? 1:2;
282 //_____________________________________________________________________________________________
283 Bool_t AliMillePede2::InitConsRecStorage(Bool_t read)
285 // initialize the buffer for processed measurements records
287 if (fConsRecFile) {AliInfo("Constraints Records File is already initialized"); return kFALSE;}
289 if (!fRecord) fRecord = new AliMillePedeRecord();
291 fConsRecFile = TFile::Open(GetConsRecFName(),read ? "":"recreate");
292 if (!fConsRecFile) {AliInfo(Form("Failed to initialize constraints records file %s",GetConsRecFName())); return kFALSE;}
294 AliInfo(Form("File %s used for constraints records",GetConsRecFName()));
296 fTreeConstr = (TTree*)fConsRecFile->Get("AliMillePedeRecords_Constraints");
297 if (!fTreeConstr) {AliInfo(Form("Did not find constraints records tree in %s",GetConsRecFName())); return kFALSE;}
298 fTreeConstr->SetBranchAddress("Record_Constraints",&fRecord);
299 AliInfo(Form("Found %lld constraints records",fTreeConstr->GetEntries()));
304 fTreeConstr = new TTree("AliMillePedeRecords_Constraints","Constraints Records for AliMillePede2");
305 fTreeConstr->Branch("Record_Constraints","AliMillePedeRecord",&fRecord,32000,99);
307 fCurrRecConstrID = -1;
312 //_____________________________________________________________________________________________
313 void AliMillePede2::CloseDataRecStorage()
316 if (fDataRecFile && fDataRecFile->IsWritable()) {
323 fDataRecFile->Close();
332 //_____________________________________________________________________________________________
333 void AliMillePede2::CloseConsRecStorage()
336 if (fConsRecFile->IsWritable()) {
338 fTreeConstr->Write();
342 fConsRecFile->Close();
349 //_____________________________________________________________________________________________
350 Bool_t AliMillePede2::ReadNextRecordData()
352 // read next data record (if any)
353 if (!fTreeData || ++fCurrRecDataID >= fTreeData->GetEntries()) { fCurrRecDataID--; return kFALSE;}
354 fTreeData->GetEntry(fCurrRecDataID);
358 //_____________________________________________________________________________________________
359 Bool_t AliMillePede2::ReadNextRecordConstraint()
361 // read next constraint record (if any)
362 if (!fTreeConstr || ++fCurrRecConstrID >= fTreeConstr->GetEntries()) { fCurrRecConstrID--; return kFALSE;}
363 fTreeConstr->GetEntry(fCurrRecConstrID);
367 //_____________________________________________________________________________________________
368 void AliMillePede2::SetRecordWeight(double wgh)
370 if (fRecFileStatus<2) InitDataRecStorage(); // create a buffer to store the data
371 fRecord->SetWeight(wgh);
374 //_____________________________________________________________________________________________
375 void AliMillePede2::SetRecordRun(Int_t run)
377 if (fRecFileStatus<2) InitDataRecStorage(); // create a buffer to store the data
378 fRecord->SetRunID(run);
381 //_____________________________________________________________________________________________
382 void AliMillePede2::SetLocalEquation(double *dergb, double *derlc, double lMeas, double lSigma)
384 if (fRecFileStatus<2) InitDataRecStorage(); // create a buffer to store the data
386 // write data of single measurement
387 if (lSigma<=0.0) { // If parameter is fixed, then no equation
388 for (int i=fNLocPar; i--;) derlc[i] = 0.0;
389 for (int i=fNGloPar; i--;) dergb[i] = 0.0;
393 fRecord->AddResidual(lMeas);
395 // Retrieve local param interesting indices
396 for (int i=0;i<fNLocPar;i++) if (!IsZero(derlc[i])) {fRecord->AddIndexValue(i,derlc[i]); derlc[i] = 0.0;}
398 fRecord->AddWeight( 1.0/lSigma/lSigma );
400 // Idem for global parameters
401 for (int i=0;i<fNGloPar;i++) if (!IsZero(dergb[i])) {
402 fRecord->AddIndexValue(i,dergb[i]); dergb[i] = 0.0;
403 fRecord->MarkGroup(fParamGrID[i]);
408 //_____________________________________________________________________________________________
409 void AliMillePede2::SetLocalEquation(int *indgb, double *dergb, int ngb, int *indlc,
410 double *derlc,int nlc,double lMeas,double lSigma)
412 // write data of single measurement
413 if (lSigma<=0.0) { // If parameter is fixed, then no equation
414 for (int i=nlc;i--;) derlc[i] = 0.0;
415 for (int i=ngb;i--;) dergb[i] = 0.0;
419 if (fRecFileStatus<2) InitDataRecStorage(); // create a buffer to store the data
421 fRecord->AddResidual(lMeas);
423 // Retrieve local param interesting indices
424 for (int i=0;i<nlc;i++) if (!IsZero(derlc[i])) {fRecord->AddIndexValue(indlc[i],derlc[i]); derlc[i]=0.; indlc[i]=0;}
426 fRecord->AddWeight( 1./lSigma/lSigma );
428 // Idem for global parameters
429 for (int i=0;i<ngb;i++) if (!IsZero(dergb[i])) {fRecord->AddIndexValue(indgb[i],dergb[i]); dergb[i]=0.; indgb[i]=0;}
434 //_____________________________________________________________________________________________
435 void AliMillePede2::SetGlobalConstraint(double *dergb, double val, double sigma)
437 // Define a constraint equation.
438 if (!fConsRecFile || !fConsRecFile->IsWritable()) InitConsRecStorage(); // create a buffer to store the data
441 fRecord->AddResidual(val);
442 fRecord->AddWeight(sigma);
443 for (int i=0; i<fNGloPar; i++) if (!IsZero(dergb[i])) fRecord->AddIndexValue(i,dergb[i]);
445 if (IsZero(sigma)) fNLagrangeConstraints++;
446 // printf("NewConstraint:\n"); fRecord->Print(); //RRR
447 SaveRecordConstraint();
451 //_____________________________________________________________________________________________
452 void AliMillePede2::SetGlobalConstraint(const int *indgb, double *dergb, int ngb, double val,double sigma)
454 // Define a constraint equation.
455 if (!fConsRecFile || !fConsRecFile->IsWritable()) InitConsRecStorage(); // create a buffer to store the data
457 fRecord->AddResidual(val);
458 fRecord->AddWeight(sigma); // dummy
459 for (int i=0; i<ngb; i++) if (!IsZero(dergb[i])) fRecord->AddIndexValue(indgb[i],dergb[i]);
461 if (IsZero(sigma)) fNLagrangeConstraints++;
462 SaveRecordConstraint();
466 //_____________________________________________________________________________________________
467 Int_t AliMillePede2::LocalFit(double *localParams)
470 Perform local parameters fit once all the local equations have been set
471 -----------------------------------------------------------
472 localParams = (if !=0) will contain the fitted track parameters and
475 static int nrefSize = 0;
476 // static TArrayI refLoc,refGlo,nrefLoc,nrefGlo;
477 static Int_t *refLoc=0,*refGlo=0,*nrefLoc=0,*nrefGlo=0;
480 AliSymMatrix &matCLoc = *fMatCLoc;
481 AliMatrixSq &matCGlo = *fMatCGlo;
482 AliRectMatrix &matCGloLoc = *fMatCGloLoc;
484 memset(fVecBLoc,0,fNLocPar*sizeof(double));
488 int recSz = fRecord->GetSize();
490 while(cnt<recSz) { // Transfer the measurement records to matrices
492 // extract addresses of residual, weight and pointers on local and global derivatives for each point
493 if (nrefSize<=nPoints) {
495 nrefSize = 2*(nPoints+1);
496 tmpA = refLoc; refLoc = new Int_t[nrefSize]; if (tmpA) memcpy(refLoc,tmpA,nPoints*sizeof(int));
497 tmpA = refGlo; refGlo = new Int_t[nrefSize]; if (tmpA) memcpy(refGlo,tmpA,nPoints*sizeof(int));
498 tmpA = nrefLoc; nrefLoc = new Int_t[nrefSize]; if (tmpA) memcpy(nrefLoc,tmpA,nPoints*sizeof(int));
499 tmpA = nrefGlo; nrefGlo = new Int_t[nrefSize]; if (tmpA) memcpy(nrefGlo,tmpA,nPoints*sizeof(int));
502 refLoc[nPoints] = ++cnt;
504 while(!fRecord->IsWeight(cnt)) {nLoc++; cnt++;}
505 nrefLoc[nPoints] = nLoc;
507 refGlo[nPoints] = ++cnt;
509 while(!fRecord->IsResidual(cnt) && cnt<recSz) {nGlo++; cnt++;}
510 nrefGlo[nPoints] = nGlo;
516 double gloWgh = fRecord->GetWeight(); // global weight for this set
517 Int_t maxLocUsed = 0;
519 for (int ip=nPoints;ip--;) { // Transfer the measurement records to matrices
520 double resid = fRecord->GetValue( refLoc[ip]-1 );
521 double weight = fRecord->GetValue( refGlo[ip]-1 )*gloWgh;
522 double *derLoc = fRecord->GetValue()+refLoc[ip];
523 double *derGlo = fRecord->GetValue()+refGlo[ip];
524 int *indLoc = fRecord->GetIndex()+refLoc[ip];
525 int *indGlo = fRecord->GetIndex()+refGlo[ip];
527 for (int i=nrefGlo[ip];i--;) { // suppress the global part (only relevant with iterations)
528 int iID = indGlo[i]; // Global param indice
529 if (fSigmaPar[iID]<=0.) continue; // fixed parameter RRRCheck
530 if (fIsLinear[iID]) resid -= derGlo[i]*(fInitPar[iID]+fDeltaPar[iID]); // linear parameter
531 else resid -= derGlo[i]*fDeltaPar[iID]; // nonlinear parameter
534 // Symmetric matrix, don't bother j>i coeffs
535 for (int i=nrefLoc[ip];i--;) { // Fill local matrix and vector
536 fVecBLoc[ indLoc[i] ] += weight*resid*derLoc[i];
537 if (indLoc[i]>maxLocUsed) maxLocUsed = indLoc[i];
538 for (int j=i+1;j--;) matCLoc(indLoc[i] ,indLoc[j]) += weight*derLoc[i]*derLoc[j];
541 } // end of the transfer of the measurement record to matrices
543 matCLoc.SetSizeUsed(++maxLocUsed); // data with B=0 may use less than declared nLocals
547 printf("\nBefore\nLocalMatrix: "); matCLoc.Print("l");
548 printf("RHSLoc: "); for (int i=0;i<fNLocPar;i++) printf("%+e |",fVecBLoc[i]); printf("\n");
550 // first try to solve by faster Cholesky decomposition, then by Gaussian elimination
551 if (!matCLoc.SolveChol(fVecBLoc,kTRUE)) {
552 AliInfo("Failed to solve locals by Cholesky, trying Gaussian Elimination");
553 if (!matCLoc.SolveSpmInv(fVecBLoc,kTRUE)) {
554 AliInfo("Failed to solve locals by Gaussian Elimination, skip...");
556 return 0; // failed to solve
560 // If requested, store the track params and errors
561 //RRR printf("locfit: "); for (int i=0;i<fNLocPar;i++) printf("%+e |",fVecBLoc[i]); printf("\n");
563 if (localParams) for (int i=maxLocUsed; i--;) {
564 localParams[2*i] = fVecBLoc[i];
565 localParams[2*i+1] = TMath::Sqrt(TMath::Abs(matCLoc.QueryDiag(i)));
571 for (int ip=nPoints;ip--;) { // Calculate residuals
572 double resid = fRecord->GetValue( refLoc[ip]-1 );
573 double weight = fRecord->GetValue( refGlo[ip]-1 )*gloWgh;
574 double *derLoc = fRecord->GetValue()+refLoc[ip];
575 double *derGlo = fRecord->GetValue()+refGlo[ip];
576 int *indLoc = fRecord->GetIndex()+refLoc[ip];
577 int *indGlo = fRecord->GetIndex()+refGlo[ip];
579 // Suppress local and global contribution in residuals;
580 for (int i=nrefLoc[ip];i--;) resid -= derLoc[i]*fVecBLoc[ indLoc[i] ]; // local part
582 for (int i=nrefGlo[ip];i--;) { // global part
584 if ( fSigmaPar[iID] <= 0.) continue; // fixed parameter RRRCheck
585 if (fIsLinear[iID]) resid -= derGlo[i]*(fInitPar[iID]+fDeltaPar[iID]); // linear parameter
586 else resid -= derGlo[i]*fDeltaPar[iID]; // nonlinear parameter
589 // reject the track if the residual is too large (outlier)
590 double absres = TMath::Abs(resid);
591 if ( (absres >= fResCutInit && fIter ==1 ) ||
592 (absres >= fResCut && fIter > 1)) {
593 if (fLocFitAdd) fNLocFitsRejected++;
594 // printf("reject res %5ld %+e\n",fCurrRecDataID,resid);
598 lChi2 += weight*resid*resid ; // total chi^2
599 nEq++; // number of equations
600 } // end of Calculate residuals
603 int nDoF = nEq-maxLocUsed;
604 lChi2 = (nDoF>0) ? lChi2/nDoF : 0; // Chi^2/dof
606 if (fNStdDev != 0 && nDoF>0 && lChi2 > Chi2DoFLim(fNStdDev,nDoF)*fChi2CutFactor) { // check final chi2
607 if (fLocFitAdd) fNLocFitsRejected++;
608 // printf("reject chi2 %5ld: %+e\n",fCurrRecDataID, lChi2);
614 fNLocEquations += nEq;
618 fNLocEquations -= nEq;
621 // local operations are finished, track is accepted
622 // We now update the global parameters (other matrices)
626 for (int ip=nPoints;ip--;) { // Update matrices
627 double resid = fRecord->GetValue( refLoc[ip]-1 );
628 double weight = fRecord->GetValue( refGlo[ip]-1 )*gloWgh;
629 double *derLoc = fRecord->GetValue()+refLoc[ip];
630 double *derGlo = fRecord->GetValue()+refGlo[ip];
631 int *indLoc = fRecord->GetIndex()+refLoc[ip];
632 int *indGlo = fRecord->GetIndex()+refGlo[ip];
634 for (int i=nrefGlo[ip];i--;) { // suppress the global part
635 int iID = indGlo[i]; // Global param indice
636 if ( fSigmaPar[iID] <= 0.) continue; // fixed parameter RRRCheck
637 if (fIsLinear[iID]) resid -= derGlo[i]*(fInitPar[iID]+fDeltaPar[iID]); // linear parameter
638 else resid -= derGlo[i]*fDeltaPar[iID]; // nonlinear parameter
641 for (int ig=nrefGlo[ip];ig--;) {
642 int iIDg = indGlo[ig]; // Global param indice (the matrix line)
643 if ( fSigmaPar[iIDg] <= 0.) continue; // fixed parameter RRRCheck
644 if (fLocFitAdd) fVecBGlo[ iIDg ] += weight*resid*derGlo[ig]; //!!!
645 else fVecBGlo[ iIDg ] -= weight*resid*derGlo[ig]; //!!!
647 // First of all, the global/global terms (exactly like local matrix)
649 for (int jg=ig+1;jg--;) { // matCGlo is symmetric by construction
650 int jIDg = indGlo[jg];
651 if ( fSigmaPar[jIDg] <= 0.) continue; // fixed parameter RRRCheck
652 if ( !IsZero(vl = weight*derGlo[ig]*derGlo[jg]) ) {
653 fFillIndex[nfill] = jIDg;
654 fFillValue[nfill++] = fLocFitAdd ? vl:-vl;
657 if (nfill) matCGlo.AddToRow(iIDg,fFillValue,fFillIndex,nfill);
659 // Now we have also rectangular matrices containing global/local terms.
660 int iCIDg = fGlo2CGlo[iIDg]; // compressed Index of index
662 Double_t *rowGL = matCGloLoc(nGloInFit);
663 for (int k=maxLocUsed;k--;) rowGL[k] = 0.0; // reset the row
664 iCIDg = fGlo2CGlo[iIDg] = nGloInFit;
665 fCGlo2Glo[nGloInFit++] = iIDg;
668 Double_t *rowGLIDg = matCGloLoc(iCIDg);
669 for (int il=nrefLoc[ip];il--;) rowGLIDg[ indLoc[il] ] += weight*derGlo[ig]*derLoc[il];
670 fProcPnt[iIDg] += fLocFitAdd ? 1:-1; // update counter
673 } // end of Update matrices
676 printf("After GLO\n");
677 printf("MatCLoc: "); fMatCLoc->Print("l");
678 printf("MatCGlo: "); fMatCGlo->Print("l");
679 printf("MatCGlLc:"); fMatCGloLoc->Print("l");
680 printf("BGlo: "); for (int i=0; i<fNGloPar; i++) printf("%+e |",fVecBGlo[i]); printf("\n");
682 // calculate fMatCGlo -= fMatCGloLoc * fMatCLoc * fMatCGloLoc^T
683 // and fVecBGlo -= fMatCGloLoc * fVecBLoc
685 //-------------------------------------------------------------- >>>
687 for (int iCIDg=0; iCIDg<nGloInFit; iCIDg++) {
688 int iIDg = fCGlo2Glo[iCIDg];
691 Double_t *rowGLIDg = matCGloLoc(iCIDg);
692 for (int kl=0;kl<maxLocUsed;kl++) if (rowGLIDg[kl]) vl += rowGLIDg[kl]*fVecBLoc[kl];
693 if (!IsZero(vl)) fVecBGlo[iIDg] -= fLocFitAdd ? vl : -vl;
696 for (int jCIDg=0;jCIDg<=iCIDg; jCIDg++) {
697 int jIDg = fCGlo2Glo[jCIDg];
700 Double_t *rowGLJDg = matCGloLoc(jCIDg);
701 for (int kl=0;kl<maxLocUsed;kl++) {
703 if ( (!IsZero(vll=rowGLIDg[kl]*rowGLJDg[kl])) ) vl += matCLoc.QueryDiag(kl)*vll;
706 for (int ll=0;ll<kl;ll++) {
707 if ( !IsZero(vll=rowGLIDg[kl]*rowGLJDg[ll]) ) vl += matCLoc(kl,ll)*vll;
708 if ( !IsZero(vll=rowGLIDg[ll]*rowGLJDg[kl]) ) vl += matCLoc(kl,ll)*vll;
712 fFillIndex[nfill] = jIDg;
713 fFillValue[nfill++] = fLocFitAdd ? -vl : vl;
716 if (nfill) matCGlo.AddToRow(iIDg,fFillValue,fFillIndex,nfill);
719 // reset compressed index array
722 printf("After GLOLoc\n");
723 printf("MatCGlo: "); fMatCGlo->Print("");
724 printf("BGlo: "); for (int i=0; i<fNGloPar; i++) printf("%+e |",fVecBGlo[i]); printf("\n");
726 for (int i=nGloInFit;i--;) {
727 fGlo2CGlo[ fCGlo2Glo[i] ] = -1;
731 //---------------------------------------------------- <<<
735 //_____________________________________________________________________________________________
736 Int_t AliMillePede2::GlobalFit(Double_t *par, Double_t *error, Double_t *pull)
738 // performs a requested number of global iterations
741 TStopwatch sw; sw.Start();
744 AliInfo("Starting Global fit.");
745 while (fIter<=fMaxIter) {
747 res = GlobalFitIteration();
750 if (!IsZero(fChi2CutFactor-fChi2CutRef)) {
751 fChi2CutFactor = TMath::Sqrt(fChi2CutFactor);
752 if (fChi2CutFactor < 1.2*fChi2CutRef) {
753 fChi2CutFactor = fChi2CutRef;
754 //RRR fIter = fMaxIter - 1; // Last iteration
761 AliInfo(Form("Global fit %s, CPU time: %.1f",res ? "Converged":"Failed",sw.CpuTime()));
764 if (par) for (int i=fNGloPar;i--;) par[i] = fInitPar[i]+fDeltaPar[i];
766 if (fGloSolveStatus==gkInvert) { // errors on params are available
767 if (error) for (int i=fNGloPar;i--;) error[i] = fProcPnt[i]>0 ? TMath::Sqrt(TMath::Abs(fMatCGlo->QueryDiag(i))) : 0.;
768 if (pull) for (int i=fNGloPar;i--;) pull[i] = fProcPnt[i]>0 && (fSigmaPar[i]*fSigmaPar[i]-fMatCGlo->QueryDiag(i))>0. && fSigmaPar[i]>0
769 ? fDeltaPar[i]/TMath::Sqrt(fSigmaPar[i]*fSigmaPar[i]-fMatCGlo->QueryDiag(i)) : 0;
775 //_____________________________________________________________________________________________
776 Int_t AliMillePede2::GlobalFitIteration()
778 // perform global parameters fit once all the local equations have been fitted
780 AliInfo(Form("Global Fit Iteration#%2d (Local Fit Chi^2 cut factor: %.2f)",fIter,fChi2CutFactor));
782 if (!fNGloPar || !fTreeData) {
783 AliInfo("No data was stored, aborting iteration");
791 fConstrUsed = new Bool_t[fNGloConstraints];
792 memset(fConstrUsed,0,fNGloConstraints*sizeof(Bool_t));
794 // Reset all info specific for this step
795 AliMatrixSq& matCGlo = *fMatCGlo;
797 memset(fProcPnt,0,fNGloPar*sizeof(Int_t));
799 fNGloConstraints = fTreeConstr ? fTreeConstr->GetEntries() : 0;
801 // count number of Lagrange constraints: they need new row/cols to be added
802 fNLagrangeConstraints = 0;
803 for (int i=0; i<fNGloConstraints; i++) {
804 ReadRecordConstraint(i);
805 if ( IsZero(fRecord->GetValue(1)) ) fNLagrangeConstraints++; // exact constraint (no error) -> Lagrange multiplier
808 // if needed, readjust the size of the global vector (for matrices this is done automatically)
809 if (!fVecBGlo || fNGloSize!=fNGloPar+fNLagrangeConstraints) {
810 delete[] fVecBGlo; // in case some constraint was added between the two manual iterations
811 fNGloSize = fNGloPar+fNLagrangeConstraints;
812 fVecBGlo = new Double_t[fNGloSize];
814 memset(fVecBGlo,0,fNGloSize*sizeof(double));
817 fNLocFitsRejected = 0;
820 // Process data records and build the matrices
821 Long_t ndr = fTreeData->GetEntries();
822 Long_t first = fSelFirst>0 ? fSelFirst : 0;
823 Long_t last = fSelLast<1 ? ndr : (fSelLast>=ndr ? ndr : fSelLast+Long_t(1));
826 AliInfo(Form("Building the Global matrix from data records %ld : %ld",first,last));
829 TStopwatch swt; swt.Start();
830 fLocFitAdd = kTRUE; // add contributions of matching tracks
831 for (Long_t i=0;i<ndr;i++) {
832 Long_t iev = i+first;
834 if (!IsRecordAcceptable()) continue;
836 if ( (i%int(0.2*ndr)) == 0) printf("%.1f%% of local fits done\n", double(100.*i)/ndr);
839 printf("%ld local fits done: ", ndr);
841 printf("MatCGlo: "); fMatCGlo->Print("l");
842 printf("BGlo: "); for (int i=0; i<fNGloPar; i++) printf("%+e |",fVecBGlo[i]); printf("\n");
848 // ---------------------- Reject parameters with low statistics ------------>>
850 if (fMinPntValid>1 && fNGroupsSet) {
852 printf("Checking parameters with statistics < %d\n",fMinPntValid);
855 // 1) build the list of parameters to fix
856 Int_t fixArrSize = 10;
857 Int_t nFixedGroups = 0;
858 TArrayI fixGroups(fixArrSize);
862 double oldMin = 1.e20;
863 double oldMax =-1.e20;
865 for (int i=fNGloPar;i--;) { // // Reset row and column of fixed params and add 1/sig^2 to free ones
866 int grID = fParamGrID[i];
867 if (grID<0) continue; // not in the group
869 if (grID!=grIDold) { // starting new group
870 if (grIDold>=0) { // decide if the group has enough statistics
871 if (oldMin<fMinPntValid && oldMax<2*fMinPntValid) { // suppress group
872 for (int iold=oldStart;iold>i;iold--) fProcPnt[iold] = 0;
873 Bool_t fnd = kFALSE; // check if the group is already accounted
874 for (int j=nFixedGroups;j--;) if (fixGroups[j]==grIDold) {fnd=kTRUE; break;}
876 if (nFixedGroups>=fixArrSize) {fixArrSize*=2; fixGroups.Set(fixArrSize);}
877 fixGroups[nFixedGroups++] = grIDold; // add group to fix
881 grIDold = grID; // mark the start of the new group
886 if (oldMin>fProcPnt[i]) oldMin = fProcPnt[i];
887 if (oldMax<fProcPnt[i]) oldMax = fProcPnt[i];
890 // extra check for the last group
891 if (grIDold>=0 && oldMin<fMinPntValid && oldMax<2*fMinPntValid) { // suppress group
892 for (int iold=oldStart;iold--;) fProcPnt[iold] = 0;
893 Bool_t fnd = kFALSE; // check if the group is already accounted
894 for (int j=nFixedGroups;j--;) if (fixGroups[j]==grIDold) {fnd=kTRUE; break;}
896 if (nFixedGroups>=fixArrSize) {fixArrSize*=2; fixGroups.Set(fixArrSize);}
897 fixGroups[nFixedGroups++] = grIDold; // add group to fix
901 // 2) loop over records and add contributions of fixed groups with negative sign
904 for (Long_t i=0;i<ndr;i++) {
905 Long_t iev = i+first;
907 if (!IsRecordAcceptable()) continue;
908 Bool_t suppr = kFALSE;
909 for (int ifx=nFixedGroups;ifx--;)if (fRecord->IsGroupPresent(fixGroups[ifx])) suppr = kTRUE;
910 if (suppr) LocalFit();
915 printf("Suppressed contributions of groups with NPoints<%d :\n",fMinPntValid);
916 for (int i=0;i<nFixedGroups;i++) printf("%d ",fixGroups[i]); printf("\n");
921 // ---------------------- Reject parameters with low statistics ------------<<
923 // add large number to diagonal of fixed params
925 for (int i=fNGloPar;i--;) { // // Reset row and column of fixed params and add 1/sig^2 to free ones
926 // printf("#%3d : Nproc : %5d grp: %d\n",i,fProcPnt[i],fParamGrID[i]);
930 matCGlo.DiagElem(i) = 1.;//float(fNLocEquations*fNLocEquations);
931 // matCGlo.DiagElem(i) = float(fNLocEquations*fNLocEquations);
933 else matCGlo.DiagElem(i) += (fgWeightSigma ? fProcPnt[i] : 1.)/(fSigmaPar[i]*fSigmaPar[i]);
936 for (int i=fNGloPar;i--;) fDiagCGlo[i] = matCGlo.QueryDiag(i); // save the diagonal elements
938 // add constraint equations
939 int nVar = fNGloPar; // Current size of global matrix
940 for (int i=0; i<fNGloConstraints; i++) {
941 ReadRecordConstraint(i);
942 double val = fRecord->GetValue(0);
943 double sig = fRecord->GetValue(1);
944 int *indV = fRecord->GetIndex()+2;
945 double *der = fRecord->GetValue()+2;
946 int csize = fRecord->GetSize()-2;
948 // check if after suppression of fixed variables there are non-0 derivatives
949 // and determine the max statistics of involved params
952 for (int j=csize;j--;) {
953 if (fProcPnt[indV[j]]<1) nSuppressed++;
955 maxStat = TMath::Max(maxStat,fProcPnt[indV[j]]);
959 if (nSuppressed==csize) {
960 // AliInfo(Form("Neglecting constraint %d of %d derivatives since no free parameters left",i,csize));
962 // was this constraint ever created ?
963 if ( sig==0 && fConstrUsed[i] ) { // this is needed only for constraints with Lagrange multiplier
964 // to avoid empty row impose dummy constraint on "Lagrange multiplier"
965 matCGlo.DiagElem(nVar) = 1.;
966 fVecBGlo[nVar++] = 0;
971 // account for already accumulated corrections
972 for (int j=csize; j--;) val -= der[j]*(fInitPar[ indV[j] ]+fDeltaPar[ indV[j] ]);
974 if (sig > 0) { // this is a gaussian constriant: no Lagrange multipliers are added
976 double sig2i = (fgWeightSigma ? TMath::Sqrt(maxStat) : 1.)/sig/sig;
977 for (int ir=0;ir<csize;ir++) {
979 for (int ic=0;ic<=ir;ic++) { // matrix is symmetric
981 double vl = der[ir]*der[ic]*sig2i;
982 if (!IsZero(vl)) matCGlo(iID,jID) += vl;
984 fVecBGlo[iID] += val*der[ir]*sig2i;
987 else { // this is exact constriant: Lagrange multipliers must be added
988 for (int j=csize; j--;) {
990 if (fProcPnt[jID]<1) continue; // this parameter was fixed, don't put it into constraint
991 matCGlo(nVar,jID) = float(fNLocEquations)*der[j]; // fMatCGlo is symmetric, only lower triangle is filled
994 if (matCGlo.QueryDiag(nVar)) matCGlo.DiagElem(nVar) = 0.0;
995 fVecBGlo[nVar++] = float(fNLocEquations)*val; //RS ? should we use here fNLocFits ?
996 fConstrUsed[i] = kTRUE;
1000 AliInfo(Form("Obtained %-7ld equations from %-7ld records (%-7ld rejected). Fixed %-4d globals",
1001 fNLocEquations,fNLocFits,fNLocFitsRejected,fNGloFix));
1007 printf("Solving:\n");
1009 for (int i=0;i<fNGloSize;i++) printf("b%2d : %+e\n",i,fVecBGlo[i]);
1011 fGloSolveStatus = SolveGlobalMatEq(); // obtain solution for this step
1013 printf("Solve %d |",fIter); sws.Print();
1016 AliInfo(Form("Iteration#%2d %s. CPU time: %.1f",fIter,fGloSolveStatus==gkFailed ? "Failed":"Converged",sw.CpuTime()));
1017 if (fGloSolveStatus==gkFailed) return 0;
1019 for (int i=fNGloPar;i--;) fDeltaPar[i] += fVecBGlo[i]; // Update global parameters values (for iterations)
1021 // PrintGlobalParameters();
1025 //_____________________________________________________________________________________________
1026 Int_t AliMillePede2::SolveGlobalMatEq()
1029 // solve global matrix equation MatCGlob*X=VecBGlo and store the result in the VecBGlo
1032 printf("GlobalMatrix\n");
1035 for (int i=0;i<fNGloPar;i++) printf("%d %+e\n",i,fVecBGlo[i]);
1038 if (!fgIsMatGloSparse) {
1040 if (fNLagrangeConstraints==0) { // pos-def systems are faster to solve by Cholesky
1041 if ( ((AliSymMatrix*)fMatCGlo)->SolveChol(fVecBGlo, fgInvChol) ) return fgInvChol ? gkInvert:gkNoInversion;
1042 else AliInfo("Solution of Global Dense System by Cholesky failed, trying Gaussian Elimiation");
1045 if (((AliSymMatrix*)fMatCGlo)->SolveSpmInv(fVecBGlo, kTRUE)) return gkInvert;
1046 else AliInfo("Solution of Global Dense System by Gaussian Elimination failed, trying iterative methods");
1048 // try to solve by minres
1049 TVectorD sol(fNGloSize);
1051 AliMinResSolve *slv = new AliMinResSolve(fMatCGlo,fVecBGlo);
1052 if (!slv) return gkFailed;
1054 Bool_t res = kFALSE;
1055 if (fgIterSol == AliMinResSolve::kSolMinRes)
1056 res = slv->SolveMinRes(sol,fgMinResCondType,fgMinResMaxIter,fgMinResTol);
1057 else if (fgIterSol == AliMinResSolve::kSolFGMRes)
1058 res = slv->SolveFGMRES(sol,fgMinResCondType,fgMinResMaxIter,fgMinResTol,fgNKrylovV);
1060 AliInfo(Form("Undefined Iteritive Solver ID=%d, only %d are defined",fgIterSol,AliMinResSolve::kNSolvers));
1063 const char* faildump = "fgmr_failed.dat";
1064 int defout = dup(1);
1066 AliInfo("Failed on dup");
1069 int slvDump = open(faildump, O_RDWR|O_CREAT, 0666);
1073 printf("#Failed to solve using solver %d with PreCond: %d MaxIter: %d Tol: %e NKrylov: %d\n",
1074 fgIterSol,fgMinResCondType,fgMinResMaxIter,fgMinResTol,fgNKrylovV);
1075 printf("#Dump of matrix:\n");
1076 fMatCGlo->Print("10");
1077 printf("#Dump of RHS:\n");
1078 for (int i=0;i<fNGloSize;i++) printf("%d %+.10f\n",i,fVecBGlo[i]);
1083 printf("#Dumped failed matrix and RHS to %s\n",faildump);
1085 else AliInfo("Failed on file open for matrix dumping");
1088 for (int i=fNGloSize;i--;) fVecBGlo[i] = sol[i];
1089 return gkNoInversion;
1093 //_____________________________________________________________________________________________
1094 Float_t AliMillePede2::Chi2DoFLim(int nSig, int nDoF) const
1096 /// return the limit in chi^2/nd for n sigmas stdev authorized
1097 // Only n=1, 2, and 3 are expected in input
1099 float sn[3] = {0.47523, 1.690140, 2.782170};
1100 float table[3][30] = {{1.0000, 1.1479, 1.1753, 1.1798, 1.1775, 1.1730, 1.1680, 1.1630,
1101 1.1581, 1.1536, 1.1493, 1.1454, 1.1417, 1.1383, 1.1351, 1.1321,
1102 1.1293, 1.1266, 1.1242, 1.1218, 1.1196, 1.1175, 1.1155, 1.1136,
1103 1.1119, 1.1101, 1.1085, 1.1070, 1.1055, 1.1040},
1104 {4.0000, 3.0900, 2.6750, 2.4290, 2.2628, 2.1415, 2.0481, 1.9736,
1105 1.9124, 1.8610, 1.8171, 1.7791, 1.7457, 1.7161, 1.6897, 1.6658,
1106 1.6442, 1.6246, 1.6065, 1.5899, 1.5745, 1.5603, 1.5470, 1.5346,
1107 1.5230, 1.5120, 1.5017, 1.4920, 1.4829, 1.4742},
1108 {9.0000, 5.9146, 4.7184, 4.0628, 3.6410, 3.3436, 3.1209, 2.9468,
1109 2.8063, 2.6902, 2.5922, 2.5082, 2.4352, 2.3711, 2.3143, 2.2635,
1110 2.2178, 2.1764, 2.1386, 2.1040, 2.0722, 2.0428, 2.0155, 1.9901,
1111 1.9665, 1.9443, 1.9235, 1.9040, 1.8855, 1.8681}};
1117 lNSig = TMath::Max(1,TMath::Min(nSig,3));
1120 return table[lNSig-1][nDoF-1];
1122 else { // approximation
1123 return ((sn[lNSig-1]+TMath::Sqrt(float(2*nDoF-3)))*
1124 (sn[lNSig-1]+TMath::Sqrt(float(2*nDoF-3))))/float(2*nDoF-2);
1129 //_____________________________________________________________________________________________
1130 Int_t AliMillePede2::SetIterations(double lChi2CutFac)
1132 // Number of iterations is calculated from lChi2CutFac
1133 fChi2CutFactor = TMath::Max(1.0, lChi2CutFac);
1135 AliInfo(Form("Initial cut factor is %f",fChi2CutFactor));
1136 fIter = 1; // Initializes the iteration process
1140 //_____________________________________________________________________________________________
1141 Double_t AliMillePede2::GetParError(int iPar) const
1143 // return error for parameter iPar
1144 if (fGloSolveStatus==gkInvert) {
1145 double res = fMatCGlo->QueryDiag(iPar);
1146 if (res>=0) return TMath::Sqrt(res);
1152 //_____________________________________________________________________________________________
1153 Int_t AliMillePede2::PrintGlobalParameters() const
1155 /// Print the final results into the logfile
1157 double lGlobalCor =0.;
1160 AliInfo(" Result of fit for global parameters");
1161 AliInfo(" ===================================");
1162 AliInfo(" I initial final differ lastcor error gcor Npnt");
1163 AliInfo("----------------------------------------------------------------------------------------------");
1165 for (int i=0; i<fNGloPar; i++) {
1166 lError = GetParError(i);
1170 if (fGloSolveStatus==gkInvert && TMath::Abs( (dg=fMatCGlo->QueryDiag(i)) *fDiagCGlo[i]) > 0) {
1171 lGlobalCor = TMath::Sqrt(TMath::Abs(1.0-1.0/(dg*fDiagCGlo[i])));
1172 AliInfo(Form("%d\t %.6f\t %.6f\t %.6f\t %.6f\t %.6f\t %.6f\t %6d",
1173 i,fInitPar[i],fInitPar[i]+fDeltaPar[i],fDeltaPar[i],fVecBGlo[i],lError,lGlobalCor,fProcPnt[i]));
1176 AliInfo(Form("%d\t %.6f\t %.6f\t %.6f\t %.6f\t OFF\t OFF\t %6d",i,fInitPar[i],fInitPar[i]+fDeltaPar[i],
1177 fDeltaPar[i],fVecBGlo[i],fProcPnt[i]));
1183 //_____________________________________________________________________________________________
1184 Bool_t AliMillePede2::IsRecordAcceptable() const
1186 // validate record according run lists set by the user
1187 static Long_t prevRunID = kMaxInt;
1188 static Bool_t prevAns = kTRUE;
1189 Long_t runID = fRecord->GetRunID();
1190 if (runID!=prevRunID) {
1193 // is run to be rejected?
1194 if (fRejRunList && (n=fRejRunList->GetSize())) {
1196 for (int i=n;i--;) if (runID == (*fRejRunList)[i]) {
1198 printf("New Run to reject: %ld -> %d\n",runID,prevAns);
1202 else if (fAccRunList && (n=fAccRunList->GetSize())) { // is run specifically selected
1204 for (int i=n;i--;) if (runID == (*fAccRunList)[i]) {prevAns = kTRUE; break;}
1212 //_____________________________________________________________________________________________
1213 void AliMillePede2::SetRejRunList(const UInt_t *runs, Int_t nruns)
1215 // set the list of runs to be rejected
1216 if (fRejRunList) delete fRejRunList;
1218 if (nruns<1 || !runs) return;
1219 fRejRunList = new TArrayL(nruns);
1220 for (int i=0;i<nruns;i++) (*fRejRunList)[i] = runs[i];
1223 //_____________________________________________________________________________________________
1224 void AliMillePede2::SetAccRunList(const UInt_t *runs, Int_t nruns)
1226 // set the list of runs to be selected
1227 if (fAccRunList) delete fAccRunList;
1229 if (nruns<1 || !runs) return;
1230 fAccRunList = new TArrayL(nruns);
1231 for (int i=0;i<nruns;i++) (*fAccRunList)[i] = runs[i];