1 /**********************************************************************************************/
2 /* General class for alignment with large number of degrees of freedom */
3 /* Based on the original milliped2 by Volker Blobel */
4 /* and AliMillepede class by Javier */
5 /* Allows operations with large sparse matrices */
6 /* http://www.desy.de/~blobel/mptalks.html */
8 /* Author: ruben.shahoyan@cern.ch */
10 /**********************************************************************************************/
12 #include "AliMillePede2.h"
14 #include <TStopwatch.h>
22 #include "AliMatrixSq.h"
23 #include "AliSymMatrix.h"
24 #include "AliRectMatrix.h"
25 #include "AliMatrixSparse.h"
29 #include <sys/types.h>
34 //#define _DUMP_EQ_BEFORE_
35 //#define _DUMP_EQ_AFTER_
37 //#define _DUMPEQ_BEFORE_
38 //#define _DUMPEQ_AFTER_
41 ClassImp(AliMillePede2)
43 Bool_t AliMillePede2::fgInvChol = kTRUE; // Invert global matrix with Cholesky solver
44 Bool_t AliMillePede2::fgWeightSigma = kTRUE; // weight local constraint by module statistics
45 Bool_t AliMillePede2::fgIsMatGloSparse = kFALSE; // use faster dense matrix by default
46 Int_t AliMillePede2::fgMinResCondType = 1; // Jacoby preconditioner by default
47 Double_t AliMillePede2::fgMinResTol = 1.e-11; // default tolerance
48 Int_t AliMillePede2::fgMinResMaxIter = 10000; // default max number of iterations
49 Int_t AliMillePede2::fgIterSol = AliMinResSolve::kSolMinRes; // default iterative solver
50 Int_t AliMillePede2::fgNKrylovV = 240; // default number of Krylov vectors to keep
52 //_____________________________________________________________________________________________
53 AliMillePede2::AliMillePede2()
64 fNLagrangeConstraints(0),
68 fGloSolveStatus(kFailed),
98 fRecDataTreeName("AliMillePedeRecords_Data"),
99 fRecConsTreeName("AliMillePedeRecords_Consaints"),
100 fRecDataBranchName("Record_Data"),
101 fRecConsBranchName("Record_Consaints"),
103 fDataRecFName("/tmp/mp2_data_records.root"),
109 fConstrRecFName("/tmp/mp2_constraints_records.root"),
115 fUseRecordWeight(kTRUE),
125 fWghScl[0] = fWghScl[1] = -1;
128 //_____________________________________________________________________________________________
129 AliMillePede2::AliMillePede2(const AliMillePede2& src) :
130 TObject(src),fNLocPar(0),fNGloPar(0),fNGloParIni(0),fNGloSize(0),fNLocEquations(0),fIter(0),
131 fMaxIter(10),fNStdDev(3),fNGloConstraints(0),fNLagrangeConstraints(0),
132 fNLocFits(0),fNLocFitsRejected(0),
133 fNGloFix(0),fGloSolveStatus(0),fChi2CutFactor(0),fChi2CutRef(0),fResCutInit(0),
134 fResCut(0),fMinPntValid(1),fNGroupsSet(0),fParamGrID(0),fProcPnt(0),fVecBLoc(0),fDiagCGlo(0),fVecBGlo(0),
135 fInitPar(0),fDeltaPar(0),fSigmaPar(0),fIsLinear(0),fConstrUsed(0),fGlo2CGlo(0),fCGlo2Glo(0),
136 fMatCLoc(0),fMatCGlo(0),fMatCGloLoc(0),fFillIndex(0),fFillValue(0),
137 fRecDataTreeName(0),fRecConsTreeName(0),fRecDataBranchName(0),fRecConsBranchName(0),
138 fDataRecFName(0),fRecord(0),fDataRecFile(0),
139 fTreeData(0),fRecFileStatus(0),fConstrRecFName(0),fTreeConstr(0),fConsRecFile(0),fCurrRecDataID(0),
140 fCurrRecConstrID(0),fLocFitAdd(kTRUE),
141 fUseRecordWeight(kTRUE),
151 fWghScl[0] = src.fWghScl[0];
152 fWghScl[1] = src.fWghScl[1];
156 //_____________________________________________________________________________________________
157 AliMillePede2::~AliMillePede2()
160 CloseDataRecStorage();
161 CloseConsRecStorage();
174 delete[] fConstrUsed;
184 delete fAccRunListWgh;
187 //_____________________________________________________________________________________________
188 Int_t AliMillePede2::InitMille(int nGlo, int nLoc, int lNStdDev,double lResCut, double lResCutInit, const Int_t* regroup)
193 if (regroup) { // regrouping is requested
195 int ng = 0; // recalculate N globals
197 for (int i=0;i<nGlo;i++) if (regroup[i]>=0) {ng++; if (regroup[i]>maxPID) maxPID = regroup[i];}
199 AliInfo(Form("Regrouping is requested: from %d raw to %d formal globals grouped to %d real globals",nGlo,ng,maxPID));
202 if (nLoc>0) fNLocPar = nLoc;
203 if (nGlo>0) fNGloPar = nGlo;
204 if (lResCutInit>0) fResCutInit = lResCutInit;
205 if (lResCut>0) fResCut = lResCut;
206 if (lNStdDev>0) fNStdDev = lNStdDev;
208 AliInfo(Form("NLoc: %d NGlo: %d",fNLocPar,fNGloPar));
210 fNGloSize = fNGloPar;
212 if (fgIsMatGloSparse) {fMatCGlo = new AliMatrixSparse(fNGloPar); fMatCGlo->SetSymmetric(kTRUE);}
213 else fMatCGlo = new AliSymMatrix(fNGloPar);
215 fFillIndex = new Int_t[fNGloPar];
216 fFillValue = new Double_t[fNGloPar];
218 fMatCLoc = new AliSymMatrix(fNLocPar);
219 fMatCGloLoc = new AliRectMatrix(fNGloPar,fNLocPar);
221 fParamGrID = new Int_t[fNGloPar];
222 fProcPnt = new Int_t[fNGloPar];
223 fVecBLoc = new Double_t[fNLocPar];
224 fDiagCGlo = new Double_t[fNGloPar];
226 fInitPar = new Double_t[fNGloPar];
227 fDeltaPar = new Double_t[fNGloPar];
228 fSigmaPar = new Double_t[fNGloPar];
229 fIsLinear = new Bool_t[fNGloPar];
231 fGlo2CGlo = new Int_t[fNGloPar];
232 fCGlo2Glo = new Int_t[fNGloPar];
234 memset(fVecBLoc ,0,fNLocPar*sizeof(Double_t));
235 memset(fDiagCGlo ,0,fNGloPar*sizeof(Double_t));
236 memset(fInitPar ,0,fNGloPar*sizeof(Double_t));
237 memset(fDeltaPar ,0,fNGloPar*sizeof(Double_t));
238 memset(fSigmaPar ,0,fNGloPar*sizeof(Double_t));
239 memset(fProcPnt ,0,fNGloPar*sizeof(Int_t));
241 for (int i=fNGloPar;i--;) {
242 fGlo2CGlo[i] = fCGlo2Glo[i] = -1;
243 fIsLinear[i] = kTRUE;
252 //_____________________________________________________________________________________________
253 Bool_t AliMillePede2::ImposeDataRecFile(const char* fname)
255 // set filename for records
256 CloseDataRecStorage();
257 SetDataRecFName(fname);
258 return InitDataRecStorage(kTRUE); // open in read mode
261 //_____________________________________________________________________________________________
262 Bool_t AliMillePede2::ImposeConsRecFile(const char* fname)
264 // set filename for constraints
265 CloseConsRecStorage();
266 SetConsRecFName(fname);
267 return InitConsRecStorage(kTRUE); // open in read mode
270 //_____________________________________________________________________________________________
271 Bool_t AliMillePede2::InitDataRecStorage(Bool_t read)
273 // initialize the buffer for processed measurements records
275 if (fTreeData) {AliInfo("Data Records File is already initialized"); return kFALSE;}
277 if (!fRecord) fRecord = new AliMillePedeRecord();
279 if (!read) { // write mode: cannot use chain
280 fDataRecFile = TFile::Open(GetDataRecFName(),"recreate");
281 if (!fDataRecFile) {AliFatal(Form("Failed to initialize data records file %s",GetDataRecFName())); return kFALSE;}
282 AliInfo(Form("File %s used for derivatives records",GetDataRecFName()));
283 fTreeData = new TTree(GetRecDataTreeName(),"Data Records for AliMillePede2");
284 fTreeData->Branch(GetRecDataBranchName(),"AliMillePedeRecord",&fRecord,32000,99);
287 TChain* ch = new TChain(GetRecDataTreeName());
289 if (fDataRecFName.EndsWith(".root")) ch->AddFile(fDataRecFName);
290 else { // assume text file with list of filenames
292 ifstream inpf(fDataRecFName.Data());
293 if (!inpf.good()) {AliInfo(Form("Failed on input records list %s\n",fDataRecFName.Data())); return kFALSE;}
296 while ( !(recfName.ReadLine(inpf)).eof() ) {
297 recfName = recfName.Strip(TString::kBoth,' ');
298 if (recfName.BeginsWith("//") || recfName.BeginsWith("#") || !recfName.EndsWith(".root")) { // comment
299 AliInfo(Form("Skip %s\n",recfName.Data()));
303 recfName = recfName.Strip(TString::kBoth,',');
304 recfName = recfName.Strip(TString::kBoth,'"');
305 gSystem->ExpandPathName(recfName);
306 printf("Adding %s\n",recfName.Data());
307 ch->AddFile(recfName.Data());
311 Long64_t nent = ch->GetEntries();
312 if (nent<1) { AliInfo("Obtained chain is empty"); return kFALSE;}
314 fTreeData->SetBranchAddress(GetRecDataBranchName(),&fRecord);
315 AliInfo(Form("Found %lld derivatives records",nent));
318 fRecFileStatus = read ? 1:2;
323 //_____________________________________________________________________________________________
324 Bool_t AliMillePede2::InitConsRecStorage(Bool_t read)
326 // initialize the buffer for processed measurements records
328 if (fConsRecFile) {AliInfo("Constraints Records File is already initialized"); return kFALSE;}
330 if (!fRecord) fRecord = new AliMillePedeRecord();
332 fConsRecFile = TFile::Open(GetConsRecFName(),read ? "":"recreate");
333 if (!fConsRecFile) {AliInfo(Form("Failed to initialize constraints records file %s",GetConsRecFName())); return kFALSE;}
335 AliInfo(Form("File %s used for constraints records",GetConsRecFName()));
337 fTreeConstr = (TTree*)fConsRecFile->Get(GetRecConsTreeName());
338 if (!fTreeConstr) {AliInfo(Form("Did not find constraints records tree in %s",GetConsRecFName())); return kFALSE;}
339 fTreeConstr->SetBranchAddress(GetRecConsBranchName(),&fRecord);
340 AliInfo(Form("Found %lld constraints records",fTreeConstr->GetEntries()));
345 fTreeConstr = new TTree(GetRecConsTreeName(),"Constraints Records for AliMillePede2");
346 fTreeConstr->Branch(GetRecConsBranchName(),"AliMillePedeRecord",&fRecord,32000,99);
348 fCurrRecConstrID = -1;
353 //_____________________________________________________________________________________________
354 void AliMillePede2::CloseDataRecStorage()
356 // close records file
358 if (fDataRecFile && fDataRecFile->IsWritable()) {
365 fDataRecFile->Close();
374 //_____________________________________________________________________________________________
375 void AliMillePede2::CloseConsRecStorage()
377 // close constraints file
379 if (fConsRecFile->IsWritable()) {
381 fTreeConstr->Write();
385 fConsRecFile->Close();
392 //_____________________________________________________________________________________________
393 Bool_t AliMillePede2::ReadNextRecordData()
395 // read next data record (if any)
396 if (!fTreeData || ++fCurrRecDataID >= fTreeData->GetEntries()) { fCurrRecDataID--; return kFALSE;}
397 fTreeData->GetEntry(fCurrRecDataID);
401 //_____________________________________________________________________________________________
402 Bool_t AliMillePede2::ReadNextRecordConstraint()
404 // read next constraint record (if any)
405 if (!fTreeConstr || ++fCurrRecConstrID >= fTreeConstr->GetEntries()) { fCurrRecConstrID--; return kFALSE;}
406 fTreeConstr->GetEntry(fCurrRecConstrID);
410 //_____________________________________________________________________________________________
411 void AliMillePede2::SetRecordWeight(double wgh)
414 if (fRecFileStatus<2) InitDataRecStorage(); // create a buffer to store the data
415 fRecord->SetWeight(wgh);
418 //_____________________________________________________________________________________________
419 void AliMillePede2::SetRecordRun(Int_t run)
422 if (fRecFileStatus<2) InitDataRecStorage(); // create a buffer to store the data
423 fRecord->SetRunID(run);
426 //_____________________________________________________________________________________________
427 void AliMillePede2::SetLocalEquation(double *dergb, double *derlc, double lMeas, double lSigma)
429 // assing derivs of loc.eq.
430 if (fRecFileStatus<2) InitDataRecStorage(); // create a buffer to store the data
432 // write data of single measurement
433 if (lSigma<=0.0) { // If parameter is fixed, then no equation
434 for (int i=fNLocPar; i--;) derlc[i] = 0.0;
435 for (int i=fNGloParIni; i--;) dergb[i] = 0.0;
439 fRecord->AddResidual(lMeas);
441 // Retrieve local param interesting indices
442 for (int i=0;i<fNLocPar;i++) if (!IsZero(derlc[i])) {fRecord->AddIndexValue(i,derlc[i]); derlc[i] = 0.0;}
444 fRecord->AddWeight( 1.0/lSigma/lSigma );
446 // Idem for global parameters
447 for (int i=0;i<fNGloParIni;i++) if (!IsZero(dergb[i])) {
448 fRecord->AddIndexValue(i,dergb[i]); dergb[i] = 0.0;
449 int idrg = GetRGId(i);
450 fRecord->MarkGroup(idrg<0 ? -1 : fParamGrID[i]);
456 //_____________________________________________________________________________________________
457 void AliMillePede2::SetLocalEquation(int *indgb, double *dergb, int ngb, int *indlc,
458 double *derlc,int nlc,double lMeas,double lSigma)
460 // write data of single measurement. Note: the records ignore regrouping, store direct parameters
461 if (lSigma<=0.0) { // If parameter is fixed, then no equation
462 for (int i=nlc;i--;) derlc[i] = 0.0;
463 for (int i=ngb;i--;) dergb[i] = 0.0;
467 if (fRecFileStatus<2) InitDataRecStorage(); // create a buffer to store the data
469 fRecord->AddResidual(lMeas);
471 // Retrieve local param interesting indices
472 for (int i=0;i<nlc;i++) if (!IsZero(derlc[i])) {fRecord->AddIndexValue(indlc[i],derlc[i]); derlc[i]=0.; indlc[i]=0;}
474 fRecord->AddWeight( 1./lSigma/lSigma );
476 // Idem for global parameters
477 for (int i=0;i<ngb;i++) if (!IsZero(dergb[i])) {fRecord->AddIndexValue(indgb[i],dergb[i]); dergb[i]=0.; indgb[i]=0;}
482 //_____________________________________________________________________________________________
483 void AliMillePede2::SetGlobalConstraint(const double *dergb, double val, double sigma)
485 // Define a constraint equation.
486 if (!fConsRecFile || !fConsRecFile->IsWritable()) InitConsRecStorage(); // create a buffer to store the data
489 fRecord->AddResidual(val);
490 fRecord->AddWeight(sigma);
491 for (int i=0; i<fNGloParIni; i++) if (!IsZero(dergb[i])) fRecord->AddIndexValue(i,dergb[i]);
493 if (IsZero(sigma)) fNLagrangeConstraints++;
494 // printf("NewConstraint:\n"); fRecord->Print(); //RRR
495 SaveRecordConstraint();
499 //_____________________________________________________________________________________________
500 void AliMillePede2::SetGlobalConstraint(const int *indgb, const double *dergb, int ngb, double val,double sigma)
502 // Define a constraint equation.
503 if (!fConsRecFile || !fConsRecFile->IsWritable()) InitConsRecStorage(); // create a buffer to store the data
505 fRecord->AddResidual(val);
506 fRecord->AddWeight(sigma); // dummy
507 for (int i=0; i<ngb; i++) if (!IsZero(dergb[i])) fRecord->AddIndexValue(indgb[i],dergb[i]);
509 if (IsZero(sigma)) fNLagrangeConstraints++;
510 SaveRecordConstraint();
514 //_____________________________________________________________________________________________
515 Int_t AliMillePede2::LocalFit(double *localParams)
518 Perform local parameters fit once all the local equations have been set
519 -----------------------------------------------------------
520 localParams = (if !=0) will contain the fitted track parameters and
523 static int nrefSize = 0;
524 // static TArrayI refLoc,refGlo,nrefLoc,nrefGlo;
525 static Int_t *refLoc=0,*refGlo=0,*nrefLoc=0,*nrefGlo=0;
528 AliSymMatrix &matCLoc = *fMatCLoc;
529 AliMatrixSq &matCGlo = *fMatCGlo;
530 AliRectMatrix &matCGloLoc = *fMatCGloLoc;
532 memset(fVecBLoc,0,fNLocPar*sizeof(double));
536 int recSz = fRecord->GetSize();
538 while(cnt<recSz) { // Transfer the measurement records to matrices
540 // extract addresses of residual, weight and pointers on local and global derivatives for each point
541 if (nrefSize<=nPoints) {
543 nrefSize = 2*(nPoints+1);
544 tmpA = refLoc; refLoc = new Int_t[nrefSize]; if (tmpA) memcpy(refLoc,tmpA,nPoints*sizeof(int));
545 tmpA = refGlo; refGlo = new Int_t[nrefSize]; if (tmpA) memcpy(refGlo,tmpA,nPoints*sizeof(int));
546 tmpA = nrefLoc; nrefLoc = new Int_t[nrefSize]; if (tmpA) memcpy(nrefLoc,tmpA,nPoints*sizeof(int));
547 tmpA = nrefGlo; nrefGlo = new Int_t[nrefSize]; if (tmpA) memcpy(nrefGlo,tmpA,nPoints*sizeof(int));
550 refLoc[nPoints] = ++cnt;
552 while(!fRecord->IsWeight(cnt)) {nLoc++; cnt++;}
553 nrefLoc[nPoints] = nLoc;
555 refGlo[nPoints] = ++cnt;
557 while(!fRecord->IsResidual(cnt) && cnt<recSz) {nGlo++; cnt++;}
558 nrefGlo[nPoints] = nGlo;
562 if (fMinRecordLength>0 && nPoints < fMinRecordLength) return 0; // ignore
566 double gloWgh = fRunWgh;
567 if (fUseRecordWeight) gloWgh *= fRecord->GetWeight(); // global weight for this set
568 Int_t maxLocUsed = 0;
570 for (int ip=nPoints;ip--;) { // Transfer the measurement records to matrices
571 double resid = fRecord->GetValue( refLoc[ip]-1 );
572 double weight = fRecord->GetValue( refGlo[ip]-1 )*gloWgh;
574 if (fWghScl[odd]>0) weight *= fWghScl[odd];
575 double *derLoc = fRecord->GetValue()+refLoc[ip];
576 double *derGlo = fRecord->GetValue()+refGlo[ip];
577 int *indLoc = fRecord->GetIndex()+refLoc[ip];
578 int *indGlo = fRecord->GetIndex()+refGlo[ip];
580 for (int i=nrefGlo[ip];i--;) { // suppress the global part (only relevant with iterations)
582 // if regrouping was requested, do it here
584 int idtmp = fkReGroup[ indGlo[i] ];
585 if (idtmp == kFixParID) indGlo[i] = kFixParID; // fixed param in regrouping
586 else indGlo[i] = idtmp;
589 int iID = indGlo[i]; // Global param indice
590 if (iID<0 || fSigmaPar[iID]<=0.) continue; // fixed parameter RRRCheck
591 if (fIsLinear[iID]) resid -= derGlo[i]*(fInitPar[iID]+fDeltaPar[iID]); // linear parameter
592 else resid -= derGlo[i]*fDeltaPar[iID]; // nonlinear parameter
595 // Symmetric matrix, don't bother j>i coeffs
596 for (int i=nrefLoc[ip];i--;) { // Fill local matrix and vector
597 fVecBLoc[ indLoc[i] ] += weight*resid*derLoc[i];
598 if (indLoc[i]>maxLocUsed) maxLocUsed = indLoc[i];
599 for (int j=i+1;j--;) matCLoc(indLoc[i] ,indLoc[j]) += weight*derLoc[i]*derLoc[j];
602 } // end of the transfer of the measurement record to matrices
604 matCLoc.SetSizeUsed(++maxLocUsed); // data with B=0 may use less than declared nLocals
608 printf("\nBefore\nLocalMatrix: "); matCLoc.Print("l");
609 printf("RHSLoc: "); for (int i=0;i<fNLocPar;i++) printf("%+e |",fVecBLoc[i]); printf("\n");
611 // first try to solve by faster Cholesky decomposition, then by Gaussian elimination
612 if (!matCLoc.SolveChol(fVecBLoc,kTRUE)) {
613 AliInfo("Failed to solve locals by Cholesky, trying Gaussian Elimination");
614 if (!matCLoc.SolveSpmInv(fVecBLoc,kTRUE)) {
615 AliInfo("Failed to solve locals by Gaussian Elimination, skip...");
617 return 0; // failed to solve
621 // If requested, store the track params and errors
622 //RRR printf("locfit: "); for (int i=0;i<fNLocPar;i++) printf("%+e |",fVecBLoc[i]); printf("\n");
624 if (localParams) for (int i=maxLocUsed; i--;) {
625 localParams[2*i] = fVecBLoc[i];
626 localParams[2*i+1] = TMath::Sqrt(TMath::Abs(matCLoc.QueryDiag(i)));
632 for (int ip=nPoints;ip--;) { // Calculate residuals
633 double resid = fRecord->GetValue( refLoc[ip]-1 );
634 double weight = fRecord->GetValue( refGlo[ip]-1 )*gloWgh;
636 if (fWghScl[odd]>0) weight *= fWghScl[odd];
637 double *derLoc = fRecord->GetValue()+refLoc[ip];
638 double *derGlo = fRecord->GetValue()+refGlo[ip];
639 int *indLoc = fRecord->GetIndex()+refLoc[ip];
640 int *indGlo = fRecord->GetIndex()+refGlo[ip];
642 // Suppress local and global contribution in residuals;
643 for (int i=nrefLoc[ip];i--;) resid -= derLoc[i]*fVecBLoc[ indLoc[i] ]; // local part
645 for (int i=nrefGlo[ip];i--;) { // global part
647 if ( iID<0 || fSigmaPar[iID] <= 0.) continue; // fixed parameter RRRCheck
648 if (fIsLinear[iID]) resid -= derGlo[i]*(fInitPar[iID]+fDeltaPar[iID]); // linear parameter
649 else resid -= derGlo[i]*fDeltaPar[iID]; // nonlinear parameter
652 // reject the track if the residual is too large (outlier)
653 double absres = TMath::Abs(resid);
654 if ( (absres >= fResCutInit && fIter ==1 ) ||
655 (absres >= fResCut && fIter > 1)) {
656 if (fLocFitAdd) fNLocFitsRejected++;
657 // printf("reject res %5ld %+e\n",fCurrRecDataID,resid);
661 lChi2 += weight*resid*resid ; // total chi^2
662 nEq++; // number of equations
663 } // end of Calculate residuals
666 int nDoF = nEq-maxLocUsed;
667 lChi2 = (nDoF>0) ? lChi2/nDoF : 0; // Chi^2/dof
669 if (fNStdDev != 0 && nDoF>0 && lChi2 > Chi2DoFLim(fNStdDev,nDoF)*fChi2CutFactor) { // check final chi2
670 if (fLocFitAdd) fNLocFitsRejected++;
671 // printf("reject chi2 %5ld: %+e\n",fCurrRecDataID, lChi2);
677 fNLocEquations += nEq;
681 fNLocEquations -= nEq;
684 // local operations are finished, track is accepted
685 // We now update the global parameters (other matrices)
689 for (int ip=nPoints;ip--;) { // Update matrices
690 double resid = fRecord->GetValue( refLoc[ip]-1 );
691 double weight = fRecord->GetValue( refGlo[ip]-1 )*gloWgh;
693 if (fWghScl[odd]>0) weight *= fWghScl[odd];
694 double *derLoc = fRecord->GetValue()+refLoc[ip];
695 double *derGlo = fRecord->GetValue()+refGlo[ip];
696 int *indLoc = fRecord->GetIndex()+refLoc[ip];
697 int *indGlo = fRecord->GetIndex()+refGlo[ip];
699 for (int i=nrefGlo[ip];i--;) { // suppress the global part
700 int iID = indGlo[i]; // Global param indice
701 if ( iID<0 || fSigmaPar[iID] <= 0.) continue; // fixed parameter RRRCheck
702 if (fIsLinear[iID]) resid -= derGlo[i]*(fInitPar[iID]+fDeltaPar[iID]); // linear parameter
703 else resid -= derGlo[i]*fDeltaPar[iID]; // nonlinear parameter
706 for (int ig=nrefGlo[ip];ig--;) {
707 int iIDg = indGlo[ig]; // Global param indice (the matrix line)
708 if ( iIDg<0 || fSigmaPar[iIDg] <= 0.) continue; // fixed parameter RRRCheck
709 if (fLocFitAdd) fVecBGlo[ iIDg ] += weight*resid*derGlo[ig]; //!!!
710 else fVecBGlo[ iIDg ] -= weight*resid*derGlo[ig]; //!!!
712 // First of all, the global/global terms (exactly like local matrix)
714 for (int jg=ig+1;jg--;) { // matCGlo is symmetric by construction
715 int jIDg = indGlo[jg];
716 if ( jIDg<0 || fSigmaPar[jIDg] <= 0.) continue; // fixed parameter RRRCheck
717 if ( !IsZero(vl = weight*derGlo[ig]*derGlo[jg]) ) {
718 fFillIndex[nfill] = jIDg;
719 fFillValue[nfill++] = fLocFitAdd ? vl:-vl;
722 if (nfill) matCGlo.AddToRow(iIDg,fFillValue,fFillIndex,nfill);
724 // Now we have also rectangular matrices containing global/local terms.
725 int iCIDg = fGlo2CGlo[iIDg]; // compressed Index of index
727 Double_t *rowGL = matCGloLoc(nGloInFit);
728 for (int k=maxLocUsed;k--;) rowGL[k] = 0.0; // reset the row
729 iCIDg = fGlo2CGlo[iIDg] = nGloInFit;
730 fCGlo2Glo[nGloInFit++] = iIDg;
733 Double_t *rowGLIDg = matCGloLoc(iCIDg);
734 for (int il=nrefLoc[ip];il--;) rowGLIDg[ indLoc[il] ] += weight*derGlo[ig]*derLoc[il];
735 fProcPnt[iIDg] += fLocFitAdd ? 1:-1; // update counter
738 } // end of Update matrices
741 printf("After GLO\n");
742 printf("MatCLoc: "); fMatCLoc->Print("l");
743 printf("MatCGlo: "); fMatCGlo->Print("l");
744 printf("MatCGlLc:"); fMatCGloLoc->Print("l");
745 printf("BGlo: "); for (int i=0; i<fNGloPar; i++) printf("%+e |",fVecBGlo[i]); printf("\n");
747 // calculate fMatCGlo -= fMatCGloLoc * fMatCLoc * fMatCGloLoc^T
748 // and fVecBGlo -= fMatCGloLoc * fVecBLoc
750 //-------------------------------------------------------------- >>>
752 for (int iCIDg=0; iCIDg<nGloInFit; iCIDg++) {
753 int iIDg = fCGlo2Glo[iCIDg];
756 Double_t *rowGLIDg = matCGloLoc(iCIDg);
757 for (int kl=0;kl<maxLocUsed;kl++) if (rowGLIDg[kl]) vl += rowGLIDg[kl]*fVecBLoc[kl];
758 if (!IsZero(vl)) fVecBGlo[iIDg] -= fLocFitAdd ? vl : -vl;
761 for (int jCIDg=0;jCIDg<=iCIDg; jCIDg++) {
762 int jIDg = fCGlo2Glo[jCIDg];
765 Double_t *rowGLJDg = matCGloLoc(jCIDg);
766 for (int kl=0;kl<maxLocUsed;kl++) {
768 if ( (!IsZero(vll=rowGLIDg[kl]*rowGLJDg[kl])) ) vl += matCLoc.QueryDiag(kl)*vll;
771 for (int ll=0;ll<kl;ll++) {
772 if ( !IsZero(vll=rowGLIDg[kl]*rowGLJDg[ll]) ) vl += matCLoc(kl,ll)*vll;
773 if ( !IsZero(vll=rowGLIDg[ll]*rowGLJDg[kl]) ) vl += matCLoc(kl,ll)*vll;
777 fFillIndex[nfill] = jIDg;
778 fFillValue[nfill++] = fLocFitAdd ? -vl : vl;
781 if (nfill) matCGlo.AddToRow(iIDg,fFillValue,fFillIndex,nfill);
784 // reset compressed index array
787 printf("After GLOLoc\n");
788 printf("MatCGlo: "); fMatCGlo->Print("");
789 printf("BGlo: "); for (int i=0; i<fNGloPar; i++) printf("%+e |",fVecBGlo[i]); printf("\n");
791 for (int i=nGloInFit;i--;) {
792 fGlo2CGlo[ fCGlo2Glo[i] ] = -1;
796 //---------------------------------------------------- <<<
800 //_____________________________________________________________________________________________
801 Int_t AliMillePede2::GlobalFit(Double_t *par, Double_t *error, Double_t *pull)
803 // performs a requested number of global iterations
806 TStopwatch sw; sw.Start();
809 AliInfo("Starting Global fit.");
810 while (fIter<=fMaxIter) {
812 res = GlobalFitIteration();
815 if (!IsZero(fChi2CutFactor-fChi2CutRef)) {
816 fChi2CutFactor = TMath::Sqrt(fChi2CutFactor);
817 if (fChi2CutFactor < 1.2*fChi2CutRef) {
818 fChi2CutFactor = fChi2CutRef;
819 //RRR fIter = fMaxIter - 1; // Last iteration
826 AliInfo(Form("Global fit %s, CPU time: %.1f",res ? "Converged":"Failed",sw.CpuTime()));
829 if (par) for (int i=fNGloParIni;i--;) par[i] = GetFinalParam(i);
831 if (fGloSolveStatus==kInvert) { // errors on params are available
832 if (error) for (int i=fNGloParIni;i--;) error[i] = GetFinalError(i);
833 if (pull) for (int i=fNGloParIni;i--;) pull[i] = GetPull(i);
839 //_____________________________________________________________________________________________
840 Int_t AliMillePede2::GlobalFitIteration()
842 // perform global parameters fit once all the local equations have been fitted
844 AliInfo(Form("Global Fit Iteration#%2d (Local Fit Chi^2 cut factor: %.2f)",fIter,fChi2CutFactor));
846 if (!fNGloPar || !fTreeData) {
847 AliInfo("No data was stored, stopping iteration");
855 fConstrUsed = new Bool_t[fNGloConstraints];
856 memset(fConstrUsed,0,fNGloConstraints*sizeof(Bool_t));
858 // Reset all info specific for this step
859 AliMatrixSq& matCGlo = *fMatCGlo;
861 memset(fProcPnt,0,fNGloPar*sizeof(Int_t));
863 fNGloConstraints = fTreeConstr ? fTreeConstr->GetEntries() : 0;
865 // count number of Lagrange constraints: they need new row/cols to be added
866 fNLagrangeConstraints = 0;
867 for (int i=0; i<fNGloConstraints; i++) {
868 ReadRecordConstraint(i);
869 if ( IsZero(fRecord->GetValue(1)) ) fNLagrangeConstraints++; // exact constraint (no error) -> Lagrange multiplier
872 // if needed, readjust the size of the global vector (for matrices this is done automatically)
873 if (!fVecBGlo || fNGloSize!=fNGloPar+fNLagrangeConstraints) {
874 delete[] fVecBGlo; // in case some constraint was added between the two manual iterations
875 fNGloSize = fNGloPar+fNLagrangeConstraints;
876 fVecBGlo = new Double_t[fNGloSize];
878 memset(fVecBGlo,0,fNGloSize*sizeof(double));
881 fNLocFitsRejected = 0;
884 // Process data records and build the matrices
885 Long_t ndr = fTreeData->GetEntries();
886 Long_t first = fSelFirst>0 ? fSelFirst : 0;
887 Long_t last = fSelLast<1 ? ndr : (fSelLast>=ndr ? ndr : fSelLast+Long_t(1));
890 AliInfo(Form("Building the Global matrix from data records %ld : %ld",first,last));
893 TStopwatch swt; swt.Start();
894 fLocFitAdd = kTRUE; // add contributions of matching tracks
895 for (Long_t i=0;i<ndr;i++) {
896 Long_t iev = i+first;
898 if (!IsRecordAcceptable()) continue;
900 if ( (i%int(0.2*ndr)) == 0) printf("%.1f%% of local fits done\n", double(100.*i)/ndr);
903 printf("%ld local fits done: ", ndr);
905 printf("MatCGlo: "); fMatCGlo->Print("l");
906 printf("BGlo: "); for (int i=0; i<fNGloPar; i++) printf("%+e |",fVecBGlo[i]); printf("\n");
912 // ---------------------- Reject parameters with low statistics ------------>>
914 if (fMinPntValid>1 && fNGroupsSet) {
916 printf("Checking parameters with statistics < %d\n",fMinPntValid);
919 // 1) build the list of parameters to fix
920 Int_t fixArrSize = 10;
921 Int_t nFixedGroups = 0;
922 TArrayI fixGroups(fixArrSize);
926 double oldMin = 1.e20;
927 double oldMax =-1.e20;
929 for (int i=fNGloPar;i--;) { // // Reset row and column of fixed params and add 1/sig^2 to free ones
930 int grID = fParamGrID[i];
931 if (grID<0) continue; // not in the group
933 if (grID!=grIDold) { // starting new group
934 if (grIDold>=0) { // decide if the group has enough statistics
935 if (oldMin<fMinPntValid && oldMax<2*fMinPntValid) { // suppress group
936 for (int iold=oldStart;iold>i;iold--) fProcPnt[iold] = 0;
937 Bool_t fnd = kFALSE; // check if the group is already accounted
938 for (int j=nFixedGroups;j--;) if (fixGroups[j]==grIDold) {fnd=kTRUE; break;}
940 if (nFixedGroups>=fixArrSize) {fixArrSize*=2; fixGroups.Set(fixArrSize);}
941 fixGroups[nFixedGroups++] = grIDold; // add group to fix
945 grIDold = grID; // mark the start of the new group
950 if (oldMin>fProcPnt[i]) oldMin = fProcPnt[i];
951 if (oldMax<fProcPnt[i]) oldMax = fProcPnt[i];
954 // extra check for the last group
955 if (grIDold>=0 && oldMin<fMinPntValid && oldMax<2*fMinPntValid) { // suppress group
956 for (int iold=oldStart;iold--;) fProcPnt[iold] = 0;
957 Bool_t fnd = kFALSE; // check if the group is already accounted
958 for (int j=nFixedGroups;j--;) if (fixGroups[j]==grIDold) {fnd=kTRUE; break;}
960 if (nFixedGroups>=fixArrSize) {fixArrSize*=2; fixGroups.Set(fixArrSize);}
961 fixGroups[nFixedGroups++] = grIDold; // add group to fix
965 // 2) loop over records and add contributions of fixed groups with negative sign
968 for (Long_t i=0;i<ndr;i++) {
969 Long_t iev = i+first;
971 if (!IsRecordAcceptable()) continue;
972 Bool_t suppr = kFALSE;
973 for (int ifx=nFixedGroups;ifx--;)if (fRecord->IsGroupPresent(fixGroups[ifx])) suppr = kTRUE;
974 if (suppr) LocalFit();
979 printf("Suppressed contributions of groups with NPoints<%d :\n",fMinPntValid);
980 for (int i=0;i<nFixedGroups;i++) printf("%d ",fixGroups[i]); printf("\n");
985 // ---------------------- Reject parameters with low statistics ------------<<
987 // add large number to diagonal of fixed params
989 for (int i=fNGloPar;i--;) { // // Reset row and column of fixed params and add 1/sig^2 to free ones
990 // printf("#%3d : Nproc : %5d grp: %d\n",i,fProcPnt[i],fParamGrID[i]);
994 matCGlo.DiagElem(i) = 1.;//float(fNLocEquations*fNLocEquations);
995 // matCGlo.DiagElem(i) = float(fNLocEquations*fNLocEquations);
997 else matCGlo.DiagElem(i) += (fgWeightSigma ? fProcPnt[i] : 1.)/(fSigmaPar[i]*fSigmaPar[i]);
1000 for (int i=fNGloPar;i--;) fDiagCGlo[i] = matCGlo.QueryDiag(i); // save the diagonal elements
1002 // add constraint equations
1003 int nVar = fNGloPar; // Current size of global matrix
1004 for (int i=0; i<fNGloConstraints; i++) {
1005 ReadRecordConstraint(i);
1006 double val = fRecord->GetValue(0);
1007 double sig = fRecord->GetValue(1);
1008 int *indV = fRecord->GetIndex()+2;
1009 double *der = fRecord->GetValue()+2;
1010 int csize = fRecord->GetSize()-2;
1013 for (int jp=csize;jp--;) {
1015 if (fkReGroup[idp]<0) AliFatal(Form("Constain is requested for suppressed parameter #%d",indV[jp]));
1019 // check if after suppression of fixed variables there are non-0 derivatives
1020 // and determine the max statistics of involved params
1021 int nSuppressed = 0;
1023 for (int j=csize;j--;) {
1024 if (fProcPnt[indV[j]]<1) nSuppressed++;
1026 maxStat = TMath::Max(maxStat,fProcPnt[indV[j]]);
1030 if (nSuppressed==csize) {
1031 // AliInfo(Form("Neglecting constraint %d of %d derivatives since no free parameters left",i,csize));
1033 // was this constraint ever created ?
1034 if ( sig==0 && fConstrUsed[i] ) { // this is needed only for constraints with Lagrange multiplier
1035 // to avoid empty row impose dummy constraint on "Lagrange multiplier"
1036 matCGlo.DiagElem(nVar) = 1.;
1037 fVecBGlo[nVar++] = 0;
1042 // account for already accumulated corrections
1043 for (int j=csize; j--;) val -= der[j]*(fInitPar[ indV[j] ]+fDeltaPar[ indV[j] ]);
1045 if (sig > 0) { // this is a gaussian constriant: no Lagrange multipliers are added
1047 double sig2i = (fgWeightSigma ? TMath::Sqrt(maxStat) : 1.)/sig/sig;
1048 for (int ir=0;ir<csize;ir++) {
1050 for (int ic=0;ic<=ir;ic++) { // matrix is symmetric
1052 double vl = der[ir]*der[ic]*sig2i;
1053 if (!IsZero(vl)) matCGlo(iID,jID) += vl;
1055 fVecBGlo[iID] += val*der[ir]*sig2i;
1058 else { // this is exact constriant: Lagrange multipliers must be added
1059 for (int j=csize; j--;) {
1061 if (fProcPnt[jID]<1) continue; // this parameter was fixed, don't put it into constraint
1062 matCGlo(nVar,jID) = float(fNLocEquations)*der[j]; // fMatCGlo is symmetric, only lower triangle is filled
1065 if (matCGlo.QueryDiag(nVar)) matCGlo.DiagElem(nVar) = 0.0;
1066 fVecBGlo[nVar++] = float(fNLocEquations)*val; //RS ? should we use here fNLocFits ?
1067 fConstrUsed[i] = kTRUE;
1071 AliInfo(Form("Obtained %-7ld equations from %-7ld records (%-7ld rejected). Fixed %-4d globals",
1072 fNLocEquations,fNLocFits,fNLocFitsRejected,fNGloFix));
1077 #ifdef _DUMP_EQ_BEFORE_
1078 const char* faildumpB = Form("mp2eq_before%d.dat",fIter);
1079 int defoutB = dup(1);
1080 if (defoutB<0) {AliFatal("Failed on dup"); exit(1);}
1081 int slvDumpB = open(faildumpB, O_RDWR|O_CREAT, 0666);
1084 printf("Solving%d for %d params\n",fIter,fNGloSize);
1085 matCGlo.Print("10");
1086 for (int i=0;i<fNGloSize;i++) printf("b%2d : %+.10f\n",i,fVecBGlo[i]);
1094 printf("Solving:\n");
1096 for (int i=0;i<fNGloSize;i++) printf("b%2d : %+e\n",i,fVecBGlo[i]);
1098 #ifdef _DUMPEQ_BEFORE_
1099 const char* faildumpB = Form("mp2eq_before%d.dat",fIter);
1100 int defoutB = dup(1);
1101 int slvDumpB = open(faildumpB, O_RDWR|O_CREAT, 0666);
1104 printf("#Equation before step %d\n",fIter);
1105 fMatCGlo->Print("10");
1106 printf("#RHS/STAT : NGlo:%d NGloSize:%d\n",fNGloPar,fNGloSize);
1107 for (int i=0;i<fNGloSize;i++) printf("%d %+.10f %d\n",i,fVecBGlo[i],fProcPnt[i]);
1114 fGloSolveStatus = SolveGlobalMatEq(); // obtain solution for this step
1115 #ifdef _DUMPEQ_AFTER_
1116 const char* faildumpA = Form("mp2eq_after%d.dat",fIter);
1117 int defoutA = dup(1);
1118 int slvDumpA = open(faildumpA, O_RDWR|O_CREAT, 0666);
1121 printf("#Matrix after step %d\n",fIter);
1122 fMatCGlo->Print("10");
1123 printf("#RHS/STAT : NGlo:%d NGloSize:%d\n",fNGloPar,fNGloSize);
1124 for (int i=0;i<fNGloSize;i++) printf("%d %+.10f %d\n",i,fVecBGlo[i],fProcPnt[i]);
1132 printf("Solve %d |",fIter); sws.Print();
1135 AliInfo(Form("Iteration#%2d %s. CPU time: %.1f",fIter,fGloSolveStatus==kFailed ? "Failed":"Converged",sw.CpuTime()));
1136 if (fGloSolveStatus==kFailed) return 0;
1138 for (int i=fNGloPar;i--;) fDeltaPar[i] += fVecBGlo[i]; // Update global parameters values (for iterations)
1140 #ifdef _DUMP_EQ_AFTER_
1141 const char* faildumpA = Form("mp2eq_after%d.dat",fIter);
1142 int defoutA = dup(1);
1143 if (defoutA<0) {AliFatal("Failed on dup"); exit(1);}
1144 int slvDumpA = open(faildumpA, O_RDWR|O_CREAT, 0666);
1147 printf("Solving%d for %d params\n",fIter,fNGloSize);
1148 matCGlo.Print("10");
1149 for (int i=0;i<fNGloSize;i++) printf("b%2d : %+.10f\n",i,fVecBGlo[i]);
1157 printf("Solved:\n");
1159 for (int i=0;i<fNGloSize;i++) printf("b%2d : %+e (->%+e)\n",i,fVecBGlo[i], fDeltaPar[i]);
1162 PrintGlobalParameters();
1166 //_____________________________________________________________________________________________
1167 Int_t AliMillePede2::SolveGlobalMatEq()
1170 // solve global matrix equation MatCGlob*X=VecBGlo and store the result in the VecBGlo
1173 printf("GlobalMatrix\n");
1174 fMatCGlo->Print("l");
1176 for (int i=0;i<fNGloPar;i++) printf("%d %+e\n",i,fVecBGlo[i]);
1179 if (!fgIsMatGloSparse) {
1181 if (fNLagrangeConstraints==0) { // pos-def systems are faster to solve by Cholesky
1182 if ( ((AliSymMatrix*)fMatCGlo)->SolveChol(fVecBGlo, fgInvChol) ) return fgInvChol ? kInvert:kNoInversion;
1183 else AliInfo("Solution of Global Dense System by Cholesky failed, trying Gaussian Elimiation");
1186 if (((AliSymMatrix*)fMatCGlo)->SolveSpmInv(fVecBGlo, kTRUE)) return kInvert;
1187 else AliInfo("Solution of Global Dense System by Gaussian Elimination failed, trying iterative methods");
1189 // try to solve by minres
1190 TVectorD sol(fNGloSize);
1192 AliMinResSolve *slv = new AliMinResSolve(fMatCGlo,fVecBGlo);
1193 if (!slv) return kFailed;
1195 Bool_t res = kFALSE;
1196 if (fgIterSol == AliMinResSolve::kSolMinRes)
1197 res = slv->SolveMinRes(sol,fgMinResCondType,fgMinResMaxIter,fgMinResTol);
1198 else if (fgIterSol == AliMinResSolve::kSolFGMRes)
1199 res = slv->SolveFGMRES(sol,fgMinResCondType,fgMinResMaxIter,fgMinResTol,fgNKrylovV);
1201 AliInfo(Form("Undefined Iteritive Solver ID=%d, only %d are defined",fgIterSol,AliMinResSolve::kNSolvers));
1204 const char* faildump = "fgmr_failed.dat";
1205 int defout = dup(1);
1207 AliInfo("Failed on dup");
1210 int slvDump = open(faildump, O_RDWR|O_CREAT, 0666);
1214 printf("#Failed to solve using solver %d with PreCond: %d MaxIter: %d Tol: %e NKrylov: %d\n",
1215 fgIterSol,fgMinResCondType,fgMinResMaxIter,fgMinResTol,fgNKrylovV);
1216 printf("#Dump of matrix:\n");
1217 fMatCGlo->Print("10");
1218 printf("#Dump of RHS:\n");
1219 for (int i=0;i<fNGloSize;i++) printf("%d %+.10f\n",i,fVecBGlo[i]);
1224 printf("#Dumped failed matrix and RHS to %s\n",faildump);
1226 else AliInfo("Failed on file open for matrix dumping");
1230 for (int i=fNGloSize;i--;) fVecBGlo[i] = sol[i];
1232 return kNoInversion;
1236 //_____________________________________________________________________________________________
1237 Float_t AliMillePede2::Chi2DoFLim(int nSig, int nDoF) const
1239 /// return the limit in chi^2/nd for n sigmas stdev authorized
1240 // Only n=1, 2, and 3 are expected in input
1242 float sn[3] = {0.47523, 1.690140, 2.782170};
1243 float table[3][30] = {{1.0000, 1.1479, 1.1753, 1.1798, 1.1775, 1.1730, 1.1680, 1.1630,
1244 1.1581, 1.1536, 1.1493, 1.1454, 1.1417, 1.1383, 1.1351, 1.1321,
1245 1.1293, 1.1266, 1.1242, 1.1218, 1.1196, 1.1175, 1.1155, 1.1136,
1246 1.1119, 1.1101, 1.1085, 1.1070, 1.1055, 1.1040},
1247 {4.0000, 3.0900, 2.6750, 2.4290, 2.2628, 2.1415, 2.0481, 1.9736,
1248 1.9124, 1.8610, 1.8171, 1.7791, 1.7457, 1.7161, 1.6897, 1.6658,
1249 1.6442, 1.6246, 1.6065, 1.5899, 1.5745, 1.5603, 1.5470, 1.5346,
1250 1.5230, 1.5120, 1.5017, 1.4920, 1.4829, 1.4742},
1251 {9.0000, 5.9146, 4.7184, 4.0628, 3.6410, 3.3436, 3.1209, 2.9468,
1252 2.8063, 2.6902, 2.5922, 2.5082, 2.4352, 2.3711, 2.3143, 2.2635,
1253 2.2178, 2.1764, 2.1386, 2.1040, 2.0722, 2.0428, 2.0155, 1.9901,
1254 1.9665, 1.9443, 1.9235, 1.9040, 1.8855, 1.8681}};
1260 lNSig = TMath::Max(1,TMath::Min(nSig,3));
1263 return table[lNSig-1][nDoF-1];
1265 else { // approximation
1266 return ((sn[lNSig-1]+TMath::Sqrt(float(2*nDoF-3)))*
1267 (sn[lNSig-1]+TMath::Sqrt(float(2*nDoF-3))))/float(2*nDoF-2);
1272 //_____________________________________________________________________________________________
1273 Int_t AliMillePede2::SetIterations(double lChi2CutFac)
1275 // Number of iterations is calculated from lChi2CutFac
1276 fChi2CutFactor = TMath::Max(1.0, lChi2CutFac);
1278 AliInfo(Form("Initial cut factor is %f",fChi2CutFactor));
1279 fIter = 1; // Initializes the iteration process
1283 //_____________________________________________________________________________________________
1284 Double_t AliMillePede2::GetParError(int iPar) const
1286 // return error for parameter iPar
1287 if (fGloSolveStatus==kInvert) {
1288 if (fkReGroup) iPar = fkReGroup[iPar];
1290 // AliDebug(2,Form("Parameter %d was suppressed in the regrouping",iPar));
1293 double res = fMatCGlo->QueryDiag(iPar);
1294 if (res>=0) return TMath::Sqrt(res);
1299 //_____________________________________________________________________________________________
1300 Double_t AliMillePede2::GetPull(int iPar) const
1302 // return pull for parameter iPar
1303 if (fGloSolveStatus==kInvert) {
1304 if (fkReGroup) iPar = fkReGroup[iPar];
1306 // AliDebug(2,Form("Parameter %d was suppressed in the regrouping",iPar));
1310 return fProcPnt[iPar]>0 && (fSigmaPar[iPar]*fSigmaPar[iPar]-fMatCGlo->QueryDiag(iPar))>0. && fSigmaPar[iPar]>0
1311 ? fDeltaPar[iPar]/TMath::Sqrt(fSigmaPar[iPar]*fSigmaPar[iPar]-fMatCGlo->QueryDiag(iPar)) : 0;
1317 //_____________________________________________________________________________________________
1318 Int_t AliMillePede2::PrintGlobalParameters() const
1320 /// Print the final results into the logfile
1322 double lGlobalCor =0.;
1324 printf("\nMillePede2 output\n");
1325 printf(" Result of fit for global parameters\n");
1326 printf(" ===================================\n");
1327 printf(" I initial final differ lastcor error gcor Npnt\n");
1328 printf("----------------------------------------------------------------------------------------------\n");
1330 int lastPrintedId = -1;
1331 for (int i0=0; i0<fNGloParIni; i0++) {
1332 int i = GetRGId(i0); if (i<0) continue;
1333 if (i!=i0 && lastPrintedId>=0 && i<=lastPrintedId) continue; // grouped param
1335 lError = GetParError(i0);
1339 if (fGloSolveStatus==kInvert && TMath::Abs( (dg=fMatCGlo->QueryDiag(i)) *fDiagCGlo[i]) > 0) {
1340 lGlobalCor = TMath::Sqrt(TMath::Abs(1.0-1.0/(dg*fDiagCGlo[i])));
1341 printf("%4d(%4d)\t %+.6f\t %+.6f\t %+.6f\t %.6f\t %.6f\t %.6f\t %6d\n",
1342 i,i0,fInitPar[i],fInitPar[i]+fDeltaPar[i],fDeltaPar[i],fVecBGlo[i],lError,lGlobalCor,fProcPnt[i]);
1345 printf("%4d (%4d)\t %+.6f\t %+.6f\t %+.6f\t %.6f\t OFF\t OFF\t %6d\n",i,i0,fInitPar[i],fInitPar[i]+fDeltaPar[i],
1346 fDeltaPar[i],fVecBGlo[i],fProcPnt[i]);
1352 //_____________________________________________________________________________________________
1353 Bool_t AliMillePede2::IsRecordAcceptable()
1355 // validate record according run lists set by the user
1356 static Long_t prevRunID = kMaxInt;
1357 static Bool_t prevAns = kTRUE;
1358 Long_t runID = fRecord->GetRunID();
1359 if (runID!=prevRunID) {
1363 // is run to be rejected?
1364 if (fRejRunList && (n=fRejRunList->GetSize())) {
1366 for (int i=n;i--;) if (runID == (*fRejRunList)[i]) {
1368 AliInfo(Form("New Run to reject: %ld",runID));
1372 else if (fAccRunList && (n=fAccRunList->GetSize())) { // is run specifically selected
1374 for (int i=n;i--;) if (runID == (*fAccRunList)[i]) {
1376 if (fAccRunListWgh) fRunWgh = (*fAccRunListWgh)[i];
1377 AliInfo(Form("New Run to accept explicitly: %ld, weight=%f",runID,fRunWgh));
1380 if (!prevAns) AliInfo(Form("New Run is not in the list to accept: %ld",runID));
1388 //_____________________________________________________________________________________________
1389 void AliMillePede2::SetRejRunList(const UInt_t *runs, Int_t nruns)
1391 // set the list of runs to be rejected
1392 if (fRejRunList) delete fRejRunList;
1394 if (nruns<1 || !runs) return;
1395 fRejRunList = new TArrayL(nruns);
1396 for (int i=0;i<nruns;i++) (*fRejRunList)[i] = runs[i];
1399 //_____________________________________________________________________________________________
1400 void AliMillePede2::SetAccRunList(const UInt_t *runs, Int_t nruns, const Float_t* wghList)
1402 // set the list of runs to be selected
1403 if (fAccRunList) delete fAccRunList;
1404 if (fAccRunListWgh) delete fAccRunListWgh;
1406 if (nruns<1 || !runs) return;
1407 fAccRunList = new TArrayL(nruns);
1408 fAccRunListWgh = new TArrayF(nruns);
1409 for (int i=0;i<nruns;i++) {
1410 (*fAccRunList)[i] = runs[i];
1411 (*fAccRunListWgh)[i] =wghList ? wghList[i] : 1.0;
1415 //_____________________________________________________________________________________________
1416 void AliMillePede2::SetInitPars(const Double_t* par)
1418 // initialize parameters, account for eventual grouping
1419 for (int i=0;i<fNGloParIni;i++) {
1420 int id = GetRGId(i); if (id<0) continue;
1421 fInitPar[id] = par[i];
1425 //_____________________________________________________________________________________________
1426 void AliMillePede2::SetSigmaPars(const Double_t* par)
1428 // initialize sigmas, account for eventual grouping
1429 for (int i=0;i<fNGloParIni;i++) {
1430 int id = GetRGId(i); if (id<0) continue;
1431 fSigmaPar[id] = par[i];
1435 //_____________________________________________________________________________________________
1436 void AliMillePede2::SetInitPar(Int_t i,Double_t par)
1438 // initialize param, account for eventual grouping
1439 int id = GetRGId(i); if (id<0) return;
1443 //_____________________________________________________________________________________________
1444 void AliMillePede2::SetSigmaPar(Int_t i,Double_t par)
1446 // initialize sigma, account for eventual grouping
1447 int id = GetRGId(i); if (id<0) return;
1448 fSigmaPar[id] = par;