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_
40 //#define _DUMPEQ_BEFORE_
41 //#define _DUMPEQ_AFTER_
43 ClassImp(AliMillePede2)
45 Bool_t AliMillePede2::fgInvChol = kTRUE; // Invert global matrix with Cholesky solver
46 Bool_t AliMillePede2::fgWeightSigma = kTRUE; // weight local constraint by module statistics
47 Bool_t AliMillePede2::fgIsMatGloSparse = kFALSE; // use faster dense matrix by default
48 Int_t AliMillePede2::fgMinResCondType = 1; // Jacoby preconditioner by default
49 Double_t AliMillePede2::fgMinResTol = 1.e-11; // default tolerance
50 Int_t AliMillePede2::fgMinResMaxIter = 10000; // default max number of iterations
51 Int_t AliMillePede2::fgIterSol = AliMinResSolve::kSolMinRes; // default iterative solver
52 Int_t AliMillePede2::fgNKrylovV = 240; // default number of Krylov vectors to keep
54 //_____________________________________________________________________________________________
55 AliMillePede2::AliMillePede2()
66 fNLagrangeConstraints(0),
70 fGloSolveStatus(kFailed),
100 fRecDataTreeName("AliMillePedeRecords_Data"),
101 fRecConsTreeName("AliMillePedeRecords_Consaints"),
102 fRecDataBranchName("Record_Data"),
103 fRecConsBranchName("Record_Consaints"),
105 fDataRecFName("/tmp/mp2_data_records.root"),
111 fConstrRecFName("/tmp/mp2_constraints_records.root"),
117 fUseRecordWeight(kTRUE),
127 fWghScl[0] = fWghScl[1] = -1;
130 //_____________________________________________________________________________________________
131 AliMillePede2::AliMillePede2(const AliMillePede2& src) :
132 TObject(src),fNLocPar(0),fNGloPar(0),fNGloParIni(0),fNGloSize(0),fNLocEquations(0),fIter(0),
133 fMaxIter(10),fNStdDev(3),fNGloConstraints(0),fNLagrangeConstraints(0),
134 fNLocFits(0),fNLocFitsRejected(0),
135 fNGloFix(0),fGloSolveStatus(0),fChi2CutFactor(0),fChi2CutRef(0),fResCutInit(0),
136 fResCut(0),fMinPntValid(1),fNGroupsSet(0),fParamGrID(0),fProcPnt(0),fVecBLoc(0),fDiagCGlo(0),fVecBGlo(0),
137 fInitPar(0),fDeltaPar(0),fSigmaPar(0),fIsLinear(0),fConstrUsed(0),fGlo2CGlo(0),fCGlo2Glo(0),
138 fMatCLoc(0),fMatCGlo(0),fMatCGloLoc(0),fFillIndex(0),fFillValue(0),
139 fRecDataTreeName(0),fRecConsTreeName(0),fRecDataBranchName(0),fRecConsBranchName(0),
140 fDataRecFName(0),fRecord(0),fDataRecFile(0),
141 fTreeData(0),fRecFileStatus(0),fConstrRecFName(0),fTreeConstr(0),fConsRecFile(0),fCurrRecDataID(0),
142 fCurrRecConstrID(0),fLocFitAdd(kTRUE),
143 fUseRecordWeight(kTRUE),
153 fWghScl[0] = src.fWghScl[0];
154 fWghScl[1] = src.fWghScl[1];
158 //_____________________________________________________________________________________________
159 AliMillePede2::~AliMillePede2()
162 CloseDataRecStorage();
163 CloseConsRecStorage();
176 delete[] fConstrUsed;
186 delete fAccRunListWgh;
189 //_____________________________________________________________________________________________
190 Int_t AliMillePede2::InitMille(int nGlo, int nLoc, int lNStdDev,double lResCut, double lResCutInit, const Int_t* regroup)
195 if (regroup) { // regrouping is requested
197 int ng = 0; // recalculate N globals
199 for (int i=0;i<nGlo;i++) if (regroup[i]>=0) {ng++; if (regroup[i]>maxPID) maxPID = regroup[i];}
201 AliInfo(Form("Regrouping is requested: from %d raw to %d formal globals grouped to %d real globals",nGlo,ng,maxPID));
204 if (nLoc>0) fNLocPar = nLoc;
205 if (nGlo>0) fNGloPar = nGlo;
206 if (lResCutInit>0) fResCutInit = lResCutInit;
207 if (lResCut>0) fResCut = lResCut;
208 if (lNStdDev>0) fNStdDev = lNStdDev;
210 AliInfo(Form("NLoc: %d NGlo: %d",fNLocPar,fNGloPar));
212 fNGloSize = fNGloPar;
214 if (fgIsMatGloSparse) {fMatCGlo = new AliMatrixSparse(fNGloPar); fMatCGlo->SetSymmetric(kTRUE);}
215 else fMatCGlo = new AliSymMatrix(fNGloPar);
217 fFillIndex = new Int_t[fNGloPar];
218 fFillValue = new Double_t[fNGloPar];
220 fMatCLoc = new AliSymMatrix(fNLocPar);
221 fMatCGloLoc = new AliRectMatrix(fNGloPar,fNLocPar);
223 fParamGrID = new Int_t[fNGloPar];
224 fProcPnt = new Int_t[fNGloPar];
225 fVecBLoc = new Double_t[fNLocPar];
226 fDiagCGlo = new Double_t[fNGloPar];
228 fInitPar = new Double_t[fNGloPar];
229 fDeltaPar = new Double_t[fNGloPar];
230 fSigmaPar = new Double_t[fNGloPar];
231 fIsLinear = new Bool_t[fNGloPar];
233 fGlo2CGlo = new Int_t[fNGloPar];
234 fCGlo2Glo = new Int_t[fNGloPar];
236 memset(fVecBLoc ,0,fNLocPar*sizeof(Double_t));
237 memset(fDiagCGlo ,0,fNGloPar*sizeof(Double_t));
238 memset(fInitPar ,0,fNGloPar*sizeof(Double_t));
239 memset(fDeltaPar ,0,fNGloPar*sizeof(Double_t));
240 memset(fSigmaPar ,0,fNGloPar*sizeof(Double_t));
241 memset(fProcPnt ,0,fNGloPar*sizeof(Int_t));
243 for (int i=fNGloPar;i--;) {
244 fGlo2CGlo[i] = fCGlo2Glo[i] = -1;
245 fIsLinear[i] = kTRUE;
254 //_____________________________________________________________________________________________
255 Bool_t AliMillePede2::ImposeDataRecFile(const char* fname)
257 // set filename for records
258 CloseDataRecStorage();
259 SetDataRecFName(fname);
260 return InitDataRecStorage(kTRUE); // open in read mode
263 //_____________________________________________________________________________________________
264 Bool_t AliMillePede2::ImposeConsRecFile(const char* fname)
266 // set filename for constraints
267 CloseConsRecStorage();
268 SetConsRecFName(fname);
269 return InitConsRecStorage(kTRUE); // open in read mode
272 //_____________________________________________________________________________________________
273 Bool_t AliMillePede2::InitDataRecStorage(Bool_t read)
275 // initialize the buffer for processed measurements records
277 if (fTreeData) {AliInfo("Data Records File is already initialized"); return kFALSE;}
279 if (!fRecord) fRecord = new AliMillePedeRecord();
281 if (!read) { // write mode: cannot use chain
282 fDataRecFile = TFile::Open(GetDataRecFName(),"recreate");
283 if (!fDataRecFile) {AliFatal(Form("Failed to initialize data records file %s",GetDataRecFName())); return kFALSE;}
284 AliInfo(Form("File %s used for derivatives records",GetDataRecFName()));
285 fTreeData = new TTree(GetRecDataTreeName(),"Data Records for AliMillePede2");
286 fTreeData->Branch(GetRecDataBranchName(),"AliMillePedeRecord",&fRecord,32000,99);
289 TChain* ch = new TChain(GetRecDataTreeName());
291 if (fDataRecFName.EndsWith(".root")) ch->AddFile(fDataRecFName);
292 else { // assume text file with list of filenames
294 ifstream inpf(fDataRecFName.Data());
295 if (!inpf.good()) {AliInfo(Form("Failed on input records list %s\n",fDataRecFName.Data())); return kFALSE;}
298 while ( !(recfName.ReadLine(inpf)).eof() ) {
299 recfName = recfName.Strip(TString::kBoth,' ');
300 if (recfName.BeginsWith("//") || recfName.BeginsWith("#") || !recfName.EndsWith(".root")) { // comment
301 AliInfo(Form("Skip %s\n",recfName.Data()));
305 recfName = recfName.Strip(TString::kBoth,',');
306 recfName = recfName.Strip(TString::kBoth,'"');
307 gSystem->ExpandPathName(recfName);
308 printf("Adding %s\n",recfName.Data());
309 ch->AddFile(recfName.Data());
313 Long64_t nent = ch->GetEntries();
314 if (nent<1) { AliInfo("Obtained chain is empty"); return kFALSE;}
316 fTreeData->SetBranchAddress(GetRecDataBranchName(),&fRecord);
317 AliInfo(Form("Found %lld derivatives records",nent));
320 fRecFileStatus = read ? 1:2;
325 //_____________________________________________________________________________________________
326 Bool_t AliMillePede2::InitConsRecStorage(Bool_t read)
328 // initialize the buffer for processed measurements records
330 if (fConsRecFile) {AliInfo("Constraints Records File is already initialized"); return kFALSE;}
332 if (!fRecord) fRecord = new AliMillePedeRecord();
334 fConsRecFile = TFile::Open(GetConsRecFName(),read ? "":"recreate");
335 if (!fConsRecFile) {AliInfo(Form("Failed to initialize constraints records file %s",GetConsRecFName())); return kFALSE;}
337 AliInfo(Form("File %s used for constraints records",GetConsRecFName()));
339 fTreeConstr = (TTree*)fConsRecFile->Get(GetRecConsTreeName());
340 if (!fTreeConstr) {AliInfo(Form("Did not find constraints records tree in %s",GetConsRecFName())); return kFALSE;}
341 fTreeConstr->SetBranchAddress(GetRecConsBranchName(),&fRecord);
342 AliInfo(Form("Found %lld constraints records",fTreeConstr->GetEntries()));
347 fTreeConstr = new TTree(GetRecConsTreeName(),"Constraints Records for AliMillePede2");
348 fTreeConstr->Branch(GetRecConsBranchName(),"AliMillePedeRecord",&fRecord,32000,99);
350 fCurrRecConstrID = -1;
355 //_____________________________________________________________________________________________
356 void AliMillePede2::CloseDataRecStorage()
358 // close records file
360 if (fDataRecFile && fDataRecFile->IsWritable()) {
367 fDataRecFile->Close();
376 //_____________________________________________________________________________________________
377 void AliMillePede2::CloseConsRecStorage()
379 // close constraints file
381 if (fConsRecFile->IsWritable()) {
383 fTreeConstr->Write();
387 fConsRecFile->Close();
394 //_____________________________________________________________________________________________
395 Bool_t AliMillePede2::ReadNextRecordData()
397 // read next data record (if any)
398 if (!fTreeData || ++fCurrRecDataID >= fTreeData->GetEntries()) { fCurrRecDataID--; return kFALSE;}
399 fTreeData->GetEntry(fCurrRecDataID);
403 //_____________________________________________________________________________________________
404 Bool_t AliMillePede2::ReadNextRecordConstraint()
406 // read next constraint record (if any)
407 if (!fTreeConstr || ++fCurrRecConstrID >= fTreeConstr->GetEntries()) { fCurrRecConstrID--; return kFALSE;}
408 fTreeConstr->GetEntry(fCurrRecConstrID);
412 //_____________________________________________________________________________________________
413 void AliMillePede2::SetRecordWeight(double wgh)
416 if (fRecFileStatus<2) InitDataRecStorage(); // create a buffer to store the data
417 fRecord->SetWeight(wgh);
420 //_____________________________________________________________________________________________
421 void AliMillePede2::SetRecordRun(Int_t run)
424 if (fRecFileStatus<2) InitDataRecStorage(); // create a buffer to store the data
425 fRecord->SetRunID(run);
428 //_____________________________________________________________________________________________
429 void AliMillePede2::SetLocalEquation(double *dergb, double *derlc, double lMeas, double lSigma)
431 // assing derivs of loc.eq.
432 if (fRecFileStatus<2) InitDataRecStorage(); // create a buffer to store the data
434 // write data of single measurement
435 if (lSigma<=0.0) { // If parameter is fixed, then no equation
436 for (int i=fNLocPar; i--;) derlc[i] = 0.0;
437 for (int i=fNGloParIni; i--;) dergb[i] = 0.0;
441 fRecord->AddResidual(lMeas);
443 // Retrieve local param interesting indices
444 for (int i=0;i<fNLocPar;i++) if (!IsZero(derlc[i])) {fRecord->AddIndexValue(i,derlc[i]); derlc[i] = 0.0;}
446 fRecord->AddWeight( 1.0/lSigma/lSigma );
448 // Idem for global parameters
449 for (int i=0;i<fNGloParIni;i++) if (!IsZero(dergb[i])) {
450 fRecord->AddIndexValue(i,dergb[i]); dergb[i] = 0.0;
451 int idrg = GetRGId(i);
452 fRecord->MarkGroup(idrg<0 ? -1 : fParamGrID[i]);
458 //_____________________________________________________________________________________________
459 void AliMillePede2::SetLocalEquation(int *indgb, double *dergb, int ngb, int *indlc,
460 double *derlc,int nlc,double lMeas,double lSigma)
462 // write data of single measurement. Note: the records ignore regrouping, store direct parameters
463 if (lSigma<=0.0) { // If parameter is fixed, then no equation
464 for (int i=nlc;i--;) derlc[i] = 0.0;
465 for (int i=ngb;i--;) dergb[i] = 0.0;
469 if (fRecFileStatus<2) InitDataRecStorage(); // create a buffer to store the data
471 fRecord->AddResidual(lMeas);
473 // Retrieve local param interesting indices
474 for (int i=0;i<nlc;i++) if (!IsZero(derlc[i])) {fRecord->AddIndexValue(indlc[i],derlc[i]); derlc[i]=0.; indlc[i]=0;}
476 fRecord->AddWeight( 1./lSigma/lSigma );
478 // Idem for global parameters
479 for (int i=0;i<ngb;i++) if (!IsZero(dergb[i])) {fRecord->AddIndexValue(indgb[i],dergb[i]); dergb[i]=0.; indgb[i]=0;}
484 //_____________________________________________________________________________________________
485 void AliMillePede2::SetGlobalConstraint(const double *dergb, double val, double sigma)
487 // Define a constraint equation.
488 if (!fConsRecFile || !fConsRecFile->IsWritable()) InitConsRecStorage(); // create a buffer to store the data
491 fRecord->AddResidual(val);
492 fRecord->AddWeight(sigma);
493 for (int i=0; i<fNGloParIni; i++) if (!IsZero(dergb[i])) fRecord->AddIndexValue(i,dergb[i]);
495 if (IsZero(sigma)) fNLagrangeConstraints++;
496 // printf("NewConstraint:\n"); fRecord->Print(); //RRR
497 SaveRecordConstraint();
501 //_____________________________________________________________________________________________
502 void AliMillePede2::SetGlobalConstraint(const int *indgb, const double *dergb, int ngb, double val,double sigma)
504 // Define a constraint equation.
505 if (!fConsRecFile || !fConsRecFile->IsWritable()) InitConsRecStorage(); // create a buffer to store the data
507 fRecord->AddResidual(val);
508 fRecord->AddWeight(sigma); // dummy
509 for (int i=0; i<ngb; i++) if (!IsZero(dergb[i])) fRecord->AddIndexValue(indgb[i],dergb[i]);
511 if (IsZero(sigma)) fNLagrangeConstraints++;
512 SaveRecordConstraint();
516 //_____________________________________________________________________________________________
517 Int_t AliMillePede2::LocalFit(double *localParams)
520 Perform local parameters fit once all the local equations have been set
521 -----------------------------------------------------------
522 localParams = (if !=0) will contain the fitted track parameters and
525 static int nrefSize = 0;
526 // static TArrayI refLoc,refGlo,nrefLoc,nrefGlo;
527 static Int_t *refLoc=0,*refGlo=0,*nrefLoc=0,*nrefGlo=0;
530 AliSymMatrix &matCLoc = *fMatCLoc;
531 AliMatrixSq &matCGlo = *fMatCGlo;
532 AliRectMatrix &matCGloLoc = *fMatCGloLoc;
534 memset(fVecBLoc,0,fNLocPar*sizeof(double));
538 int recSz = fRecord->GetSize();
540 while(cnt<recSz) { // Transfer the measurement records to matrices
542 // extract addresses of residual, weight and pointers on local and global derivatives for each point
543 if (nrefSize<=nPoints) {
545 nrefSize = 2*(nPoints+1);
546 tmpA = refLoc; refLoc = new Int_t[nrefSize]; if (tmpA) memcpy(refLoc,tmpA,nPoints*sizeof(int));
547 tmpA = refGlo; refGlo = new Int_t[nrefSize]; if (tmpA) memcpy(refGlo,tmpA,nPoints*sizeof(int));
548 tmpA = nrefLoc; nrefLoc = new Int_t[nrefSize]; if (tmpA) memcpy(nrefLoc,tmpA,nPoints*sizeof(int));
549 tmpA = nrefGlo; nrefGlo = new Int_t[nrefSize]; if (tmpA) memcpy(nrefGlo,tmpA,nPoints*sizeof(int));
552 refLoc[nPoints] = ++cnt;
554 while(!fRecord->IsWeight(cnt)) {nLoc++; cnt++;}
555 nrefLoc[nPoints] = nLoc;
557 refGlo[nPoints] = ++cnt;
559 while(!fRecord->IsResidual(cnt) && cnt<recSz) {nGlo++; cnt++;}
560 nrefGlo[nPoints] = nGlo;
564 if (fMinRecordLength>0 && nPoints < fMinRecordLength) return 0; // ignore
568 double gloWgh = fRunWgh;
569 if (fUseRecordWeight) gloWgh *= fRecord->GetWeight(); // global weight for this set
570 Int_t maxLocUsed = 0;
572 for (int ip=nPoints;ip--;) { // Transfer the measurement records to matrices
573 double resid = fRecord->GetValue( refLoc[ip]-1 );
574 double weight = fRecord->GetValue( refGlo[ip]-1 )*gloWgh;
576 if (fWghScl[odd]>0) weight *= fWghScl[odd];
577 double *derLoc = fRecord->GetValue()+refLoc[ip];
578 double *derGlo = fRecord->GetValue()+refGlo[ip];
579 int *indLoc = fRecord->GetIndex()+refLoc[ip];
580 int *indGlo = fRecord->GetIndex()+refGlo[ip];
582 for (int i=nrefGlo[ip];i--;) { // suppress the global part (only relevant with iterations)
584 // if regrouping was requested, do it here
586 int idtmp = fkReGroup[ indGlo[i] ];
587 if (idtmp == kFixParID) indGlo[i] = kFixParID; // fixed param in regrouping
588 else indGlo[i] = idtmp;
591 int iID = indGlo[i]; // Global param indice
592 if (iID<0 || fSigmaPar[iID]<=0.) continue; // fixed parameter RRRCheck
593 if (fIsLinear[iID]) resid -= derGlo[i]*(fInitPar[iID]+fDeltaPar[iID]); // linear parameter
594 else resid -= derGlo[i]*fDeltaPar[iID]; // nonlinear parameter
597 // Symmetric matrix, don't bother j>i coeffs
598 for (int i=nrefLoc[ip];i--;) { // Fill local matrix and vector
599 fVecBLoc[ indLoc[i] ] += weight*resid*derLoc[i];
600 if (indLoc[i]>maxLocUsed) maxLocUsed = indLoc[i];
601 for (int j=i+1;j--;) matCLoc(indLoc[i] ,indLoc[j]) += weight*derLoc[i]*derLoc[j];
604 } // end of the transfer of the measurement record to matrices
606 matCLoc.SetSizeUsed(++maxLocUsed); // data with B=0 may use less than declared nLocals
610 printf("\nBefore\nLocalMatrix: "); matCLoc.Print("l");
611 printf("RHSLoc: "); for (int i=0;i<fNLocPar;i++) printf("%+e |",fVecBLoc[i]); printf("\n");
613 // first try to solve by faster Cholesky decomposition, then by Gaussian elimination
614 if (!matCLoc.SolveChol(fVecBLoc,kTRUE)) {
615 AliInfo("Failed to solve locals by Cholesky, trying Gaussian Elimination");
616 if (!matCLoc.SolveSpmInv(fVecBLoc,kTRUE)) {
617 AliInfo("Failed to solve locals by Gaussian Elimination, skip...");
619 return 0; // failed to solve
623 // If requested, store the track params and errors
624 //RRR printf("locfit: "); for (int i=0;i<fNLocPar;i++) printf("%+e |",fVecBLoc[i]); printf("\n");
626 if (localParams) for (int i=maxLocUsed; i--;) {
627 localParams[2*i] = fVecBLoc[i];
628 localParams[2*i+1] = TMath::Sqrt(TMath::Abs(matCLoc.QueryDiag(i)));
634 for (int ip=nPoints;ip--;) { // Calculate residuals
635 double resid = fRecord->GetValue( refLoc[ip]-1 );
636 double weight = fRecord->GetValue( refGlo[ip]-1 )*gloWgh;
638 if (fWghScl[odd]>0) weight *= fWghScl[odd];
639 double *derLoc = fRecord->GetValue()+refLoc[ip];
640 double *derGlo = fRecord->GetValue()+refGlo[ip];
641 int *indLoc = fRecord->GetIndex()+refLoc[ip];
642 int *indGlo = fRecord->GetIndex()+refGlo[ip];
644 // Suppress local and global contribution in residuals;
645 for (int i=nrefLoc[ip];i--;) resid -= derLoc[i]*fVecBLoc[ indLoc[i] ]; // local part
647 for (int i=nrefGlo[ip];i--;) { // global part
649 if ( iID<0 || fSigmaPar[iID] <= 0.) continue; // fixed parameter RRRCheck
650 if (fIsLinear[iID]) resid -= derGlo[i]*(fInitPar[iID]+fDeltaPar[iID]); // linear parameter
651 else resid -= derGlo[i]*fDeltaPar[iID]; // nonlinear parameter
654 // reject the track if the residual is too large (outlier)
655 double absres = TMath::Abs(resid);
656 if ( (absres >= fResCutInit && fIter ==1 ) ||
657 (absres >= fResCut && fIter > 1)) {
658 if (fLocFitAdd) fNLocFitsRejected++;
659 // printf("reject res %5ld %+e\n",fCurrRecDataID,resid);
663 lChi2 += weight*resid*resid ; // total chi^2
664 nEq++; // number of equations
665 } // end of Calculate residuals
668 int nDoF = nEq-maxLocUsed;
669 lChi2 = (nDoF>0) ? lChi2/nDoF : 0; // Chi^2/dof
671 if (fNStdDev != 0 && nDoF>0 && lChi2 > Chi2DoFLim(fNStdDev,nDoF)*fChi2CutFactor) { // check final chi2
672 if (fLocFitAdd) fNLocFitsRejected++;
673 // printf("reject chi2 %5ld: %+e\n",fCurrRecDataID, lChi2);
679 fNLocEquations += nEq;
683 fNLocEquations -= nEq;
686 // local operations are finished, track is accepted
687 // We now update the global parameters (other matrices)
691 for (int ip=nPoints;ip--;) { // Update matrices
692 double resid = fRecord->GetValue( refLoc[ip]-1 );
693 double weight = fRecord->GetValue( refGlo[ip]-1 )*gloWgh;
695 if (fWghScl[odd]>0) weight *= fWghScl[odd];
696 double *derLoc = fRecord->GetValue()+refLoc[ip];
697 double *derGlo = fRecord->GetValue()+refGlo[ip];
698 int *indLoc = fRecord->GetIndex()+refLoc[ip];
699 int *indGlo = fRecord->GetIndex()+refGlo[ip];
701 for (int i=nrefGlo[ip];i--;) { // suppress the global part
702 int iID = indGlo[i]; // Global param indice
703 if ( iID<0 || fSigmaPar[iID] <= 0.) continue; // fixed parameter RRRCheck
704 if (fIsLinear[iID]) resid -= derGlo[i]*(fInitPar[iID]+fDeltaPar[iID]); // linear parameter
705 else resid -= derGlo[i]*fDeltaPar[iID]; // nonlinear parameter
708 for (int ig=nrefGlo[ip];ig--;) {
709 int iIDg = indGlo[ig]; // Global param indice (the matrix line)
710 if ( iIDg<0 || fSigmaPar[iIDg] <= 0.) continue; // fixed parameter RRRCheck
711 if (fLocFitAdd) fVecBGlo[ iIDg ] += weight*resid*derGlo[ig]; //!!!
712 else fVecBGlo[ iIDg ] -= weight*resid*derGlo[ig]; //!!!
714 // First of all, the global/global terms (exactly like local matrix)
716 for (int jg=ig+1;jg--;) { // matCGlo is symmetric by construction
717 int jIDg = indGlo[jg];
718 if ( jIDg<0 || fSigmaPar[jIDg] <= 0.) continue; // fixed parameter RRRCheck
719 if ( !IsZero(vl = weight*derGlo[ig]*derGlo[jg]) ) {
720 fFillIndex[nfill] = jIDg;
721 fFillValue[nfill++] = fLocFitAdd ? vl:-vl;
724 if (nfill) matCGlo.AddToRow(iIDg,fFillValue,fFillIndex,nfill);
726 // Now we have also rectangular matrices containing global/local terms.
727 int iCIDg = fGlo2CGlo[iIDg]; // compressed Index of index
729 Double_t *rowGL = matCGloLoc(nGloInFit);
730 for (int k=maxLocUsed;k--;) rowGL[k] = 0.0; // reset the row
731 iCIDg = fGlo2CGlo[iIDg] = nGloInFit;
732 fCGlo2Glo[nGloInFit++] = iIDg;
735 Double_t *rowGLIDg = matCGloLoc(iCIDg);
736 for (int il=nrefLoc[ip];il--;) rowGLIDg[ indLoc[il] ] += weight*derGlo[ig]*derLoc[il];
737 fProcPnt[iIDg] += fLocFitAdd ? 1:-1; // update counter
740 } // end of Update matrices
743 printf("After GLO\n");
744 printf("MatCLoc: "); fMatCLoc->Print("l");
745 printf("MatCGlo: "); fMatCGlo->Print("l");
746 printf("MatCGlLc:"); fMatCGloLoc->Print("l");
747 printf("BGlo: "); for (int i=0; i<fNGloPar; i++) printf("%+e |",fVecBGlo[i]); printf("\n");
749 // calculate fMatCGlo -= fMatCGloLoc * fMatCLoc * fMatCGloLoc^T
750 // and fVecBGlo -= fMatCGloLoc * fVecBLoc
752 //-------------------------------------------------------------- >>>
754 for (int iCIDg=0; iCIDg<nGloInFit; iCIDg++) {
755 int iIDg = fCGlo2Glo[iCIDg];
758 Double_t *rowGLIDg = matCGloLoc(iCIDg);
759 for (int kl=0;kl<maxLocUsed;kl++) if (rowGLIDg[kl]) vl += rowGLIDg[kl]*fVecBLoc[kl];
760 if (!IsZero(vl)) fVecBGlo[iIDg] -= fLocFitAdd ? vl : -vl;
763 for (int jCIDg=0;jCIDg<=iCIDg; jCIDg++) {
764 int jIDg = fCGlo2Glo[jCIDg];
767 Double_t *rowGLJDg = matCGloLoc(jCIDg);
768 for (int kl=0;kl<maxLocUsed;kl++) {
770 if ( (!IsZero(vll=rowGLIDg[kl]*rowGLJDg[kl])) ) vl += matCLoc.QueryDiag(kl)*vll;
773 for (int ll=0;ll<kl;ll++) {
774 if ( !IsZero(vll=rowGLIDg[kl]*rowGLJDg[ll]) ) vl += matCLoc(kl,ll)*vll;
775 if ( !IsZero(vll=rowGLIDg[ll]*rowGLJDg[kl]) ) vl += matCLoc(kl,ll)*vll;
779 fFillIndex[nfill] = jIDg;
780 fFillValue[nfill++] = fLocFitAdd ? -vl : vl;
783 if (nfill) matCGlo.AddToRow(iIDg,fFillValue,fFillIndex,nfill);
786 // reset compressed index array
789 printf("After GLOLoc\n");
790 printf("MatCGlo: "); fMatCGlo->Print("");
791 printf("BGlo: "); for (int i=0; i<fNGloPar; i++) printf("%+e |",fVecBGlo[i]); printf("\n");
793 for (int i=nGloInFit;i--;) {
794 fGlo2CGlo[ fCGlo2Glo[i] ] = -1;
798 //---------------------------------------------------- <<<
802 //_____________________________________________________________________________________________
803 Int_t AliMillePede2::GlobalFit(Double_t *par, Double_t *error, Double_t *pull)
805 // performs a requested number of global iterations
808 TStopwatch sw; sw.Start();
811 AliInfo("Starting Global fit.");
812 while (fIter<=fMaxIter) {
814 res = GlobalFitIteration();
817 if (!IsZero(fChi2CutFactor-fChi2CutRef)) {
818 fChi2CutFactor = TMath::Sqrt(fChi2CutFactor);
819 if (fChi2CutFactor < 1.2*fChi2CutRef) {
820 fChi2CutFactor = fChi2CutRef;
821 //RRR fIter = fMaxIter - 1; // Last iteration
828 AliInfo(Form("Global fit %s, CPU time: %.1f",res ? "Converged":"Failed",sw.CpuTime()));
831 if (par) for (int i=fNGloParIni;i--;) par[i] = GetFinalParam(i);
833 if (fGloSolveStatus==kInvert) { // errors on params are available
834 if (error) for (int i=fNGloParIni;i--;) error[i] = GetFinalError(i);
835 if (pull) for (int i=fNGloParIni;i--;) pull[i] = GetPull(i);
841 //_____________________________________________________________________________________________
842 Int_t AliMillePede2::GlobalFitIteration()
844 // perform global parameters fit once all the local equations have been fitted
846 AliInfo(Form("Global Fit Iteration#%2d (Local Fit Chi^2 cut factor: %.2f)",fIter,fChi2CutFactor));
848 if (!fNGloPar || !fTreeData) {
849 AliInfo("No data was stored, stopping iteration");
857 fConstrUsed = new Bool_t[fNGloConstraints];
858 memset(fConstrUsed,0,fNGloConstraints*sizeof(Bool_t));
860 // Reset all info specific for this step
861 AliMatrixSq& matCGlo = *fMatCGlo;
863 memset(fProcPnt,0,fNGloPar*sizeof(Int_t));
865 fNGloConstraints = fTreeConstr ? fTreeConstr->GetEntries() : 0;
867 // count number of Lagrange constraints: they need new row/cols to be added
868 fNLagrangeConstraints = 0;
869 for (int i=0; i<fNGloConstraints; i++) {
870 ReadRecordConstraint(i);
871 if ( IsZero(fRecord->GetValue(1)) ) fNLagrangeConstraints++; // exact constraint (no error) -> Lagrange multiplier
874 // if needed, readjust the size of the global vector (for matrices this is done automatically)
875 if (!fVecBGlo || fNGloSize!=fNGloPar+fNLagrangeConstraints) {
876 delete[] fVecBGlo; // in case some constraint was added between the two manual iterations
877 fNGloSize = fNGloPar+fNLagrangeConstraints;
878 fVecBGlo = new Double_t[fNGloSize];
880 memset(fVecBGlo,0,fNGloSize*sizeof(double));
883 fNLocFitsRejected = 0;
886 // Process data records and build the matrices
887 Long_t ndr = fTreeData->GetEntries();
888 Long_t first = fSelFirst>0 ? fSelFirst : 0;
889 Long_t last = fSelLast<1 ? ndr : (fSelLast>=ndr ? ndr : fSelLast+Long_t(1));
892 AliInfo(Form("Building the Global matrix from data records %ld : %ld",first,last));
895 TStopwatch swt; swt.Start();
896 fLocFitAdd = kTRUE; // add contributions of matching tracks
897 for (Long_t i=0;i<ndr;i++) {
898 Long_t iev = i+first;
900 if (!IsRecordAcceptable()) continue;
902 if ( (i%int(0.2*ndr)) == 0) printf("%.1f%% of local fits done\n", double(100.*i)/ndr);
905 printf("%ld local fits done: ", ndr);
907 printf("MatCGlo: "); fMatCGlo->Print("l");
908 printf("BGlo: "); for (int i=0; i<fNGloPar; i++) printf("%+e |",fVecBGlo[i]); printf("\n");
914 // ---------------------- Reject parameters with low statistics ------------>>
916 if (fMinPntValid>1 && fNGroupsSet) {
918 printf("Checking parameters with statistics < %d\n",fMinPntValid);
921 // 1) build the list of parameters to fix
922 Int_t fixArrSize = 10;
923 Int_t nFixedGroups = 0;
924 TArrayI fixGroups(fixArrSize);
928 double oldMin = 1.e20;
929 double oldMax =-1.e20;
931 for (int i=fNGloPar;i--;) { // // Reset row and column of fixed params and add 1/sig^2 to free ones
932 int grID = fParamGrID[i];
933 if (grID<0) continue; // not in the group
935 if (grID!=grIDold) { // starting new group
936 if (grIDold>=0) { // decide if the group has enough statistics
937 if (oldMin<fMinPntValid && oldMax<2*fMinPntValid) { // suppress group
938 for (int iold=oldStart;iold>i;iold--) fProcPnt[iold] = 0;
939 Bool_t fnd = kFALSE; // check if the group is already accounted
940 for (int j=nFixedGroups;j--;) if (fixGroups[j]==grIDold) {fnd=kTRUE; break;}
942 if (nFixedGroups>=fixArrSize) {fixArrSize*=2; fixGroups.Set(fixArrSize);}
943 fixGroups[nFixedGroups++] = grIDold; // add group to fix
947 grIDold = grID; // mark the start of the new group
952 if (oldMin>fProcPnt[i]) oldMin = fProcPnt[i];
953 if (oldMax<fProcPnt[i]) oldMax = fProcPnt[i];
956 // extra check for the last group
957 if (grIDold>=0 && oldMin<fMinPntValid && oldMax<2*fMinPntValid) { // suppress group
958 for (int iold=oldStart;iold--;) fProcPnt[iold] = 0;
959 Bool_t fnd = kFALSE; // check if the group is already accounted
960 for (int j=nFixedGroups;j--;) if (fixGroups[j]==grIDold) {fnd=kTRUE; break;}
962 if (nFixedGroups>=fixArrSize) {fixArrSize*=2; fixGroups.Set(fixArrSize);}
963 fixGroups[nFixedGroups++] = grIDold; // add group to fix
967 // 2) loop over records and add contributions of fixed groups with negative sign
970 for (Long_t i=0;i<ndr;i++) {
971 Long_t iev = i+first;
973 if (!IsRecordAcceptable()) continue;
974 Bool_t suppr = kFALSE;
975 for (int ifx=nFixedGroups;ifx--;)if (fRecord->IsGroupPresent(fixGroups[ifx])) suppr = kTRUE;
976 if (suppr) LocalFit();
981 printf("Suppressed contributions of groups with NPoints<%d :\n",fMinPntValid);
982 for (int i=0;i<nFixedGroups;i++) printf("%d ",fixGroups[i]); printf("\n");
987 // ---------------------- Reject parameters with low statistics ------------<<
989 // add large number to diagonal of fixed params
991 for (int i=fNGloPar;i--;) { // // Reset row and column of fixed params and add 1/sig^2 to free ones
992 // printf("#%3d : Nproc : %5d grp: %d\n",i,fProcPnt[i],fParamGrID[i]);
996 matCGlo.DiagElem(i) = 1.;//float(fNLocEquations*fNLocEquations);
997 // matCGlo.DiagElem(i) = float(fNLocEquations*fNLocEquations);
999 else matCGlo.DiagElem(i) += (fgWeightSigma ? fProcPnt[i] : 1.)/(fSigmaPar[i]*fSigmaPar[i]);
1002 for (int i=fNGloPar;i--;) fDiagCGlo[i] = matCGlo.QueryDiag(i); // save the diagonal elements
1004 // add constraint equations
1005 int nVar = fNGloPar; // Current size of global matrix
1006 for (int i=0; i<fNGloConstraints; i++) {
1007 ReadRecordConstraint(i);
1008 double val = fRecord->GetValue(0);
1009 double sig = fRecord->GetValue(1);
1010 int *indV = fRecord->GetIndex()+2;
1011 double *der = fRecord->GetValue()+2;
1012 int csize = fRecord->GetSize()-2;
1015 for (int jp=csize;jp--;) {
1017 if (fkReGroup[idp]<0) AliFatal(Form("Constain is requested for suppressed parameter #%d",indV[jp]));
1021 // check if after suppression of fixed variables there are non-0 derivatives
1022 // and determine the max statistics of involved params
1023 int nSuppressed = 0;
1025 for (int j=csize;j--;) {
1026 if (fProcPnt[indV[j]]<1) nSuppressed++;
1028 maxStat = TMath::Max(maxStat,fProcPnt[indV[j]]);
1032 if (nSuppressed==csize) {
1033 // AliInfo(Form("Neglecting constraint %d of %d derivatives since no free parameters left",i,csize));
1035 // was this constraint ever created ?
1036 if ( sig==0 && fConstrUsed[i] ) { // this is needed only for constraints with Lagrange multiplier
1037 // to avoid empty row impose dummy constraint on "Lagrange multiplier"
1038 matCGlo.DiagElem(nVar) = 1.;
1039 fVecBGlo[nVar++] = 0;
1044 // account for already accumulated corrections
1045 for (int j=csize; j--;) val -= der[j]*(fInitPar[ indV[j] ]+fDeltaPar[ indV[j] ]);
1047 if (sig > 0) { // this is a gaussian constriant: no Lagrange multipliers are added
1049 double sig2i = (fgWeightSigma ? TMath::Sqrt(maxStat) : 1.)/sig/sig;
1050 for (int ir=0;ir<csize;ir++) {
1052 for (int ic=0;ic<=ir;ic++) { // matrix is symmetric
1054 double vl = der[ir]*der[ic]*sig2i;
1055 if (!IsZero(vl)) matCGlo(iID,jID) += vl;
1057 fVecBGlo[iID] += val*der[ir]*sig2i;
1060 else { // this is exact constriant: Lagrange multipliers must be added
1061 for (int j=csize; j--;) {
1063 if (fProcPnt[jID]<1) continue; // this parameter was fixed, don't put it into constraint
1064 matCGlo(nVar,jID) = float(fNLocEquations)*der[j]; // fMatCGlo is symmetric, only lower triangle is filled
1067 if (matCGlo.QueryDiag(nVar)) matCGlo.DiagElem(nVar) = 0.0;
1068 fVecBGlo[nVar++] = float(fNLocEquations)*val; //RS ? should we use here fNLocFits ?
1069 fConstrUsed[i] = kTRUE;
1073 AliInfo(Form("Obtained %-7ld equations from %-7ld records (%-7ld rejected). Fixed %-4d globals",
1074 fNLocEquations,fNLocFits,fNLocFitsRejected,fNGloFix));
1079 #ifdef _DUMP_EQ_BEFORE_
1080 const char* faildumpB = Form("mp2eq_before%d.dat",fIter);
1081 int defoutB = dup(1);
1082 if (defoutB<0) AliFatal("Failed on dup");
1083 int slvDumpB = open(faildumpB, O_RDWR|O_CREAT, 0666);
1086 printf("Solving%d for %d params\n",fIter,fNGloSize);
1087 matCGlo.Print("10");
1088 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]);
1113 fGloSolveStatus = SolveGlobalMatEq(); // obtain solution for this step
1114 #ifdef _DUMPEQ_AFTER_
1115 const char* faildumpA = Form("mp2eq_after%d.dat",fIter);
1116 int defoutA = dup(1);
1117 int slvDumpA = open(faildumpA, O_RDWR|O_CREAT, 0666);
1120 printf("#Matrix after step %d\n",fIter);
1121 fMatCGlo->Print("10");
1122 printf("#RHS/STAT : NGlo:%d NGloSize:%d\n",fNGloPar,fNGloSize);
1123 for (int i=0;i<fNGloSize;i++) printf("%d %+.10f %d\n",i,fVecBGlo[i],fProcPnt[i]);
1130 printf("Solve %d |",fIter); sws.Print();
1133 AliInfo(Form("Iteration#%2d %s. CPU time: %.1f",fIter,fGloSolveStatus==kFailed ? "Failed":"Converged",sw.CpuTime()));
1134 if (fGloSolveStatus==kFailed) return 0;
1136 for (int i=fNGloPar;i--;) fDeltaPar[i] += fVecBGlo[i]; // Update global parameters values (for iterations)
1138 #ifdef _DUMP_EQ_AFTER_
1139 const char* faildumpA = Form("mp2eq_after%d.dat",fIter);
1140 int defoutA = dup(1);
1141 if (defoutA<0) AliFatal("Failed on dup");
1142 int slvDumpA = open(faildumpA, O_RDWR|O_CREAT, 0666);
1145 printf("Solving%d for %d params\n",fIter,fNGloSize);
1146 matCGlo.Print("10");
1147 for (int i=0;i<fNGloSize;i++) printf("b%2d : %+.10f\n",i,fVecBGlo[i]);
1154 printf("Solved:\n");
1156 for (int i=0;i<fNGloSize;i++) printf("b%2d : %+e (->%+e)\n",i,fVecBGlo[i], fDeltaPar[i]);
1159 PrintGlobalParameters();
1163 //_____________________________________________________________________________________________
1164 Int_t AliMillePede2::SolveGlobalMatEq()
1167 // solve global matrix equation MatCGlob*X=VecBGlo and store the result in the VecBGlo
1170 printf("GlobalMatrix\n");
1171 fMatCGlo->Print("l");
1173 for (int i=0;i<fNGloPar;i++) printf("%d %+e\n",i,fVecBGlo[i]);
1176 if (!fgIsMatGloSparse) {
1178 if (fNLagrangeConstraints==0) { // pos-def systems are faster to solve by Cholesky
1179 if ( ((AliSymMatrix*)fMatCGlo)->SolveChol(fVecBGlo, fgInvChol) ) return fgInvChol ? kInvert:kNoInversion;
1180 else AliInfo("Solution of Global Dense System by Cholesky failed, trying Gaussian Elimiation");
1183 if (((AliSymMatrix*)fMatCGlo)->SolveSpmInv(fVecBGlo, kTRUE)) return kInvert;
1184 else AliInfo("Solution of Global Dense System by Gaussian Elimination failed, trying iterative methods");
1186 // try to solve by minres
1187 TVectorD sol(fNGloSize);
1189 AliMinResSolve *slv = new AliMinResSolve(fMatCGlo,fVecBGlo);
1190 if (!slv) return kFailed;
1192 Bool_t res = kFALSE;
1193 if (fgIterSol == AliMinResSolve::kSolMinRes)
1194 res = slv->SolveMinRes(sol,fgMinResCondType,fgMinResMaxIter,fgMinResTol);
1195 else if (fgIterSol == AliMinResSolve::kSolFGMRes)
1196 res = slv->SolveFGMRES(sol,fgMinResCondType,fgMinResMaxIter,fgMinResTol,fgNKrylovV);
1198 AliInfo(Form("Undefined Iteritive Solver ID=%d, only %d are defined",fgIterSol,AliMinResSolve::kNSolvers));
1201 const char* faildump = "fgmr_failed.dat";
1202 int defout = dup(1);
1204 AliInfo("Failed on dup");
1207 int slvDump = open(faildump, O_RDWR|O_CREAT, 0666);
1211 printf("#Failed to solve using solver %d with PreCond: %d MaxIter: %d Tol: %e NKrylov: %d\n",
1212 fgIterSol,fgMinResCondType,fgMinResMaxIter,fgMinResTol,fgNKrylovV);
1213 printf("#Dump of matrix:\n");
1214 fMatCGlo->Print("10");
1215 printf("#Dump of RHS:\n");
1216 for (int i=0;i<fNGloSize;i++) printf("%d %+.10f\n",i,fVecBGlo[i]);
1220 printf("#Dumped failed matrix and RHS to %s\n",faildump);
1222 else AliInfo("Failed on file open for matrix dumping");
1226 for (int i=fNGloSize;i--;) fVecBGlo[i] = sol[i];
1228 return kNoInversion;
1232 //_____________________________________________________________________________________________
1233 Float_t AliMillePede2::Chi2DoFLim(int nSig, int nDoF) const
1235 /// return the limit in chi^2/nd for n sigmas stdev authorized
1236 // Only n=1, 2, and 3 are expected in input
1238 float sn[3] = {0.47523, 1.690140, 2.782170};
1239 float table[3][30] = {{1.0000, 1.1479, 1.1753, 1.1798, 1.1775, 1.1730, 1.1680, 1.1630,
1240 1.1581, 1.1536, 1.1493, 1.1454, 1.1417, 1.1383, 1.1351, 1.1321,
1241 1.1293, 1.1266, 1.1242, 1.1218, 1.1196, 1.1175, 1.1155, 1.1136,
1242 1.1119, 1.1101, 1.1085, 1.1070, 1.1055, 1.1040},
1243 {4.0000, 3.0900, 2.6750, 2.4290, 2.2628, 2.1415, 2.0481, 1.9736,
1244 1.9124, 1.8610, 1.8171, 1.7791, 1.7457, 1.7161, 1.6897, 1.6658,
1245 1.6442, 1.6246, 1.6065, 1.5899, 1.5745, 1.5603, 1.5470, 1.5346,
1246 1.5230, 1.5120, 1.5017, 1.4920, 1.4829, 1.4742},
1247 {9.0000, 5.9146, 4.7184, 4.0628, 3.6410, 3.3436, 3.1209, 2.9468,
1248 2.8063, 2.6902, 2.5922, 2.5082, 2.4352, 2.3711, 2.3143, 2.2635,
1249 2.2178, 2.1764, 2.1386, 2.1040, 2.0722, 2.0428, 2.0155, 1.9901,
1250 1.9665, 1.9443, 1.9235, 1.9040, 1.8855, 1.8681}};
1256 lNSig = TMath::Max(1,TMath::Min(nSig,3));
1259 return table[lNSig-1][nDoF-1];
1261 else { // approximation
1262 return ((sn[lNSig-1]+TMath::Sqrt(float(2*nDoF-3)))*
1263 (sn[lNSig-1]+TMath::Sqrt(float(2*nDoF-3))))/float(2*nDoF-2);
1268 //_____________________________________________________________________________________________
1269 Int_t AliMillePede2::SetIterations(double lChi2CutFac)
1271 // Number of iterations is calculated from lChi2CutFac
1272 fChi2CutFactor = TMath::Max(1.0, lChi2CutFac);
1274 AliInfo(Form("Initial cut factor is %f",fChi2CutFactor));
1275 fIter = 1; // Initializes the iteration process
1279 //_____________________________________________________________________________________________
1280 Double_t AliMillePede2::GetParError(int iPar) const
1282 // return error for parameter iPar
1283 if (fGloSolveStatus==kInvert) {
1284 if (fkReGroup) iPar = fkReGroup[iPar];
1286 // AliDebug(2,Form("Parameter %d was suppressed in the regrouping",iPar));
1289 double res = fMatCGlo->QueryDiag(iPar);
1290 if (res>=0) return TMath::Sqrt(res);
1295 //_____________________________________________________________________________________________
1296 Double_t AliMillePede2::GetPull(int iPar) const
1298 // return pull for parameter iPar
1299 if (fGloSolveStatus==kInvert) {
1300 if (fkReGroup) iPar = fkReGroup[iPar];
1302 // AliDebug(2,Form("Parameter %d was suppressed in the regrouping",iPar));
1306 return fProcPnt[iPar]>0 && (fSigmaPar[iPar]*fSigmaPar[iPar]-fMatCGlo->QueryDiag(iPar))>0. && fSigmaPar[iPar]>0
1307 ? fDeltaPar[iPar]/TMath::Sqrt(fSigmaPar[iPar]*fSigmaPar[iPar]-fMatCGlo->QueryDiag(iPar)) : 0;
1313 //_____________________________________________________________________________________________
1314 Int_t AliMillePede2::PrintGlobalParameters() const
1316 /// Print the final results into the logfile
1318 double lGlobalCor =0.;
1320 printf("\nMillePede2 output\n");
1321 printf(" Result of fit for global parameters\n");
1322 printf(" ===================================\n");
1323 printf(" I initial final differ lastcor error gcor Npnt\n");
1324 printf("----------------------------------------------------------------------------------------------\n");
1326 int lastPrintedId = -1;
1327 for (int i0=0; i0<fNGloParIni; i0++) {
1328 int i = GetRGId(i0); if (i<0) continue;
1329 if (i!=i0 && lastPrintedId>=0 && i<=lastPrintedId) continue; // grouped param
1331 lError = GetParError(i0);
1335 if (fGloSolveStatus==kInvert && TMath::Abs( (dg=fMatCGlo->QueryDiag(i)) *fDiagCGlo[i]) > 0) {
1336 lGlobalCor = TMath::Sqrt(TMath::Abs(1.0-1.0/(dg*fDiagCGlo[i])));
1337 printf("%4d(%4d)\t %+.6f\t %+.6f\t %+.6f\t %.6f\t %.6f\t %.6f\t %6d\n",
1338 i,i0,fInitPar[i],fInitPar[i]+fDeltaPar[i],fDeltaPar[i],fVecBGlo[i],lError,lGlobalCor,fProcPnt[i]);
1341 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],
1342 fDeltaPar[i],fVecBGlo[i],fProcPnt[i]);
1348 //_____________________________________________________________________________________________
1349 Bool_t AliMillePede2::IsRecordAcceptable()
1351 // validate record according run lists set by the user
1352 static Long_t prevRunID = kMaxInt;
1353 static Bool_t prevAns = kTRUE;
1354 Long_t runID = fRecord->GetRunID();
1355 if (runID!=prevRunID) {
1359 // is run to be rejected?
1360 if (fRejRunList && (n=fRejRunList->GetSize())) {
1362 for (int i=n;i--;) if (runID == (*fRejRunList)[i]) {
1364 printf("New Run to reject: %ld -> %d\n",runID,prevAns);
1368 else if (fAccRunList && (n=fAccRunList->GetSize())) { // is run specifically selected
1370 for (int i=n;i--;) if (runID == (*fAccRunList)[i]) {prevAns = kTRUE; fRunWgh = (*fAccRunListWgh)[i]; break;}
1378 //_____________________________________________________________________________________________
1379 void AliMillePede2::SetRejRunList(const UInt_t *runs, Int_t nruns)
1381 // set the list of runs to be rejected
1382 if (fRejRunList) delete fRejRunList;
1384 if (nruns<1 || !runs) return;
1385 fRejRunList = new TArrayL(nruns);
1386 for (int i=0;i<nruns;i++) (*fRejRunList)[i] = runs[i];
1389 //_____________________________________________________________________________________________
1390 void AliMillePede2::SetAccRunList(const UInt_t *runs, Int_t nruns, const Float_t* wghList)
1392 // set the list of runs to be selected
1393 if (fAccRunList) delete fAccRunList;
1394 if (fAccRunListWgh) delete fAccRunListWgh;
1396 if (nruns<1 || !runs) return;
1397 fAccRunList = new TArrayL(nruns);
1398 fAccRunListWgh = new TArrayF(nruns);
1399 for (int i=0;i<nruns;i++) {
1400 (*fAccRunList)[i] = runs[i];
1401 (*fAccRunListWgh)[i] =wghList ? wghList[i] : 1.0;
1405 //_____________________________________________________________________________________________
1406 void AliMillePede2::SetInitPars(const Double_t* par)
1408 // initialize parameters, account for eventual grouping
1409 for (int i=0;i<fNGloParIni;i++) {
1410 int id = GetRGId(i); if (id<0) continue;
1411 fInitPar[id] = par[i];
1415 //_____________________________________________________________________________________________
1416 void AliMillePede2::SetSigmaPars(const Double_t* par)
1418 // initialize sigmas, account for eventual grouping
1419 for (int i=0;i<fNGloParIni;i++) {
1420 int id = GetRGId(i); if (id<0) continue;
1421 fSigmaPar[id] = par[i];
1425 //_____________________________________________________________________________________________
1426 void AliMillePede2::SetInitPar(Int_t i,Double_t par)
1428 // initialize param, account for eventual grouping
1429 int id = GetRGId(i); if (id<0) return;
1433 //_____________________________________________________________________________________________
1434 void AliMillePede2::SetSigmaPar(Int_t i,Double_t par)
1436 // initialize sigma, account for eventual grouping
1437 int id = GetRGId(i); if (id<0) return;
1438 fSigmaPar[id] = par;