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>
21 #include "AliMatrixSq.h"
22 #include "AliSymMatrix.h"
23 #include "AliRectMatrix.h"
24 #include "AliMatrixSparse.h"
28 #include <sys/types.h>
36 ClassImp(AliMillePede2)
38 Bool_t AliMillePede2::fgInvChol = kTRUE; // Invert global matrix with Cholesky solver
39 Bool_t AliMillePede2::fgWeightSigma = kTRUE; // weight local constraint by module statistics
40 Bool_t AliMillePede2::fgIsMatGloSparse = kFALSE; // use faster dense matrix by default
41 Int_t AliMillePede2::fgMinResCondType = 1; // Jacoby preconditioner by default
42 Double_t AliMillePede2::fgMinResTol = 1.e-11; // default tolerance
43 Int_t AliMillePede2::fgMinResMaxIter = 10000; // default max number of iterations
44 Int_t AliMillePede2::fgIterSol = AliMinResSolve::kSolMinRes; // default iterative solver
45 Int_t AliMillePede2::fgNKrylovV = 240; // default number of Krylov vectors to keep
47 //_____________________________________________________________________________________________
48 AliMillePede2::AliMillePede2()
59 fNLagrangeConstraints(0),
63 fGloSolveStatus(kFailed),
93 fRecDataTreeName("AliMillePedeRecords_Data"),
94 fRecConsTreeName("AliMillePedeRecords_Consaints"),
95 fRecDataBranchName("Record_Data"),
96 fRecConsBranchName("Record_Consaints"),
98 fDataRecFName("/tmp/mp2_data_records.root"),
104 fConstrRecFName("/tmp/mp2_constraints_records.root"),
117 //_____________________________________________________________________________________________
118 AliMillePede2::AliMillePede2(const AliMillePede2& src) :
119 TObject(src),fNLocPar(0),fNGloPar(0),fNGloParIni(0),fNGloSize(0),fNLocEquations(0),fIter(0),
120 fMaxIter(10),fNStdDev(3),fNGloConstraints(0),fNLagrangeConstraints(0),
121 fNLocFits(0),fNLocFitsRejected(0),
122 fNGloFix(0),fGloSolveStatus(0),fChi2CutFactor(0),fChi2CutRef(0),fResCutInit(0),
123 fResCut(0),fMinPntValid(1),fNGroupsSet(0),fParamGrID(0),fProcPnt(0),fVecBLoc(0),fDiagCGlo(0),fVecBGlo(0),
124 fInitPar(0),fDeltaPar(0),fSigmaPar(0),fIsLinear(0),fConstrUsed(0),fGlo2CGlo(0),fCGlo2Glo(0),
125 fMatCLoc(0),fMatCGlo(0),fMatCGloLoc(0),fFillIndex(0),fFillValue(0),
126 fRecDataTreeName(0),fRecConsTreeName(0),fRecDataBranchName(0),fRecConsBranchName(0),
127 fDataRecFName(0),fRecord(0),fDataRecFile(0),
128 fTreeData(0),fRecFileStatus(0),fConstrRecFName(0),fTreeConstr(0),fConsRecFile(0),fCurrRecDataID(0),
129 fCurrRecConstrID(0),fLocFitAdd(kTRUE),
137 //_____________________________________________________________________________________________
138 AliMillePede2::~AliMillePede2()
141 CloseDataRecStorage();
142 CloseConsRecStorage();
155 delete[] fConstrUsed;
167 //_____________________________________________________________________________________________
168 Int_t AliMillePede2::InitMille(int nGlo, int nLoc, int lNStdDev,double lResCut, double lResCutInit, const Int_t* regroup)
173 if (regroup) { // regrouping is requested
175 int ng = 0; // recalculate N globals
177 for (int i=0;i<nGlo;i++) if (regroup[i]>=0) {ng++; if (regroup[i]>maxPID) maxPID = regroup[i];}
178 AliInfo(Form("Gegrouping is requested: from %d raw to %d free globals, max global ID=%d",nGlo,ng,maxPID));
179 if (ng != maxPID-1) AliFatal("Wrong regrouping requested");
182 if (nLoc>0) fNLocPar = nLoc;
183 if (nGlo>0) fNGloPar = nGlo;
184 if (lResCutInit>0) fResCutInit = lResCutInit;
185 if (lResCut>0) fResCut = lResCut;
186 if (lNStdDev>0) fNStdDev = lNStdDev;
188 AliInfo(Form("NLoc: %d NGlo: %d",fNLocPar,fNGloPar));
190 fNGloSize = fNGloPar;
192 if (fgIsMatGloSparse) {fMatCGlo = new AliMatrixSparse(fNGloPar); fMatCGlo->SetSymmetric(kTRUE);}
193 else fMatCGlo = new AliSymMatrix(fNGloPar);
195 fFillIndex = new Int_t[fNGloPar];
196 fFillValue = new Double_t[fNGloPar];
198 fMatCLoc = new AliSymMatrix(fNLocPar);
199 fMatCGloLoc = new AliRectMatrix(fNGloPar,fNLocPar);
201 fParamGrID = new Int_t[fNGloPar];
202 fProcPnt = new Int_t[fNGloPar];
203 fVecBLoc = new Double_t[fNLocPar];
204 fDiagCGlo = new Double_t[fNGloPar];
206 fInitPar = new Double_t[fNGloPar];
207 fDeltaPar = new Double_t[fNGloPar];
208 fSigmaPar = new Double_t[fNGloPar];
209 fIsLinear = new Bool_t[fNGloPar];
211 fGlo2CGlo = new Int_t[fNGloPar];
212 fCGlo2Glo = new Int_t[fNGloPar];
214 memset(fVecBLoc ,0,fNLocPar*sizeof(Double_t));
215 memset(fDiagCGlo ,0,fNGloPar*sizeof(Double_t));
216 memset(fInitPar ,0,fNGloPar*sizeof(Double_t));
217 memset(fDeltaPar ,0,fNGloPar*sizeof(Double_t));
218 memset(fSigmaPar ,0,fNGloPar*sizeof(Double_t));
219 memset(fProcPnt ,0,fNGloPar*sizeof(Int_t));
221 for (int i=fNGloPar;i--;) {
222 fGlo2CGlo[i] = fCGlo2Glo[i] = -1;
223 fIsLinear[i] = kTRUE;
230 //_____________________________________________________________________________________________
231 Bool_t AliMillePede2::ImposeDataRecFile(const char* fname)
233 // set filename for records
234 CloseDataRecStorage();
235 SetDataRecFName(fname);
236 return InitDataRecStorage(kTRUE); // open in read mode
239 //_____________________________________________________________________________________________
240 Bool_t AliMillePede2::ImposeConsRecFile(const char* fname)
242 // set filename for constraints
243 CloseConsRecStorage();
244 SetConsRecFName(fname);
245 return InitConsRecStorage(kTRUE); // open in read mode
248 //_____________________________________________________________________________________________
249 Bool_t AliMillePede2::InitDataRecStorage(Bool_t read)
251 // initialize the buffer for processed measurements records
253 if (fTreeData) {AliInfo("Data Records File is already initialized"); return kFALSE;}
255 if (!fRecord) fRecord = new AliMillePedeRecord();
257 if (!read) { // write mode: cannot use chain
258 fDataRecFile = TFile::Open(GetDataRecFName(),"recreate");
259 if (!fDataRecFile) {AliFatal(Form("Failed to initialize data records file %s",GetDataRecFName())); return kFALSE;}
260 AliInfo(Form("File %s used for derivatives records",GetDataRecFName()));
261 fTreeData = new TTree(GetRecDataTreeName(),"Data Records for AliMillePede2");
262 fTreeData->Branch(GetRecDataBranchName(),"AliMillePedeRecord",&fRecord,32000,99);
265 TChain* ch = new TChain(GetRecDataTreeName());
267 if (fDataRecFName.EndsWith(".root")) ch->AddFile(fDataRecFName);
268 else { // assume text file with list of filenames
270 ifstream inpf(fDataRecFName.Data());
271 if (!inpf.good()) {AliInfo(Form("Failed on input records list %s\n",fDataRecFName.Data())); return kFALSE;}
274 while ( !(recfName.ReadLine(inpf)).eof() ) {
275 recfName = recfName.Strip(TString::kBoth,' ');
276 if (recfName.BeginsWith("//") || recfName.BeginsWith("#") || !recfName.EndsWith(".root")) { // comment
277 AliInfo(Form("Skip %s\n",recfName.Data()));
281 recfName = recfName.Strip(TString::kBoth,',');
282 recfName = recfName.Strip(TString::kBoth,'"');
283 gSystem->ExpandPathName(recfName);
284 printf("Adding %s\n",recfName.Data());
285 ch->AddFile(recfName.Data());
289 Long64_t nent = ch->GetEntries();
290 if (nent<1) { AliInfo("Obtained chain is empty"); return kFALSE;}
292 fTreeData->SetBranchAddress(GetRecDataBranchName(),&fRecord);
293 AliInfo(Form("Found %lld derivatives records",nent));
296 fRecFileStatus = read ? 1:2;
301 //_____________________________________________________________________________________________
302 Bool_t AliMillePede2::InitConsRecStorage(Bool_t read)
304 // initialize the buffer for processed measurements records
306 if (fConsRecFile) {AliInfo("Constraints Records File is already initialized"); return kFALSE;}
308 if (!fRecord) fRecord = new AliMillePedeRecord();
310 fConsRecFile = TFile::Open(GetConsRecFName(),read ? "":"recreate");
311 if (!fConsRecFile) {AliInfo(Form("Failed to initialize constraints records file %s",GetConsRecFName())); return kFALSE;}
313 AliInfo(Form("File %s used for constraints records",GetConsRecFName()));
315 fTreeConstr = (TTree*)fConsRecFile->Get(GetRecConsTreeName());
316 if (!fTreeConstr) {AliInfo(Form("Did not find constraints records tree in %s",GetConsRecFName())); return kFALSE;}
317 fTreeConstr->SetBranchAddress(GetRecConsBranchName(),&fRecord);
318 AliInfo(Form("Found %lld constraints records",fTreeConstr->GetEntries()));
323 fTreeConstr = new TTree(GetRecConsTreeName(),"Constraints Records for AliMillePede2");
324 fTreeConstr->Branch(GetRecConsBranchName(),"AliMillePedeRecord",&fRecord,32000,99);
326 fCurrRecConstrID = -1;
331 //_____________________________________________________________________________________________
332 void AliMillePede2::CloseDataRecStorage()
334 // close records file
336 if (fDataRecFile && fDataRecFile->IsWritable()) {
343 fDataRecFile->Close();
352 //_____________________________________________________________________________________________
353 void AliMillePede2::CloseConsRecStorage()
355 // close constraints file
357 if (fConsRecFile->IsWritable()) {
359 fTreeConstr->Write();
363 fConsRecFile->Close();
370 //_____________________________________________________________________________________________
371 Bool_t AliMillePede2::ReadNextRecordData()
373 // read next data record (if any)
374 if (!fTreeData || ++fCurrRecDataID >= fTreeData->GetEntries()) { fCurrRecDataID--; return kFALSE;}
375 fTreeData->GetEntry(fCurrRecDataID);
379 //_____________________________________________________________________________________________
380 Bool_t AliMillePede2::ReadNextRecordConstraint()
382 // read next constraint record (if any)
383 if (!fTreeConstr || ++fCurrRecConstrID >= fTreeConstr->GetEntries()) { fCurrRecConstrID--; return kFALSE;}
384 fTreeConstr->GetEntry(fCurrRecConstrID);
388 //_____________________________________________________________________________________________
389 void AliMillePede2::SetRecordWeight(double wgh)
392 if (fRecFileStatus<2) InitDataRecStorage(); // create a buffer to store the data
393 fRecord->SetWeight(wgh);
396 //_____________________________________________________________________________________________
397 void AliMillePede2::SetRecordRun(Int_t run)
400 if (fRecFileStatus<2) InitDataRecStorage(); // create a buffer to store the data
401 fRecord->SetRunID(run);
404 //_____________________________________________________________________________________________
405 void AliMillePede2::SetLocalEquation(double *dergb, double *derlc, double lMeas, double lSigma)
407 // assing derivs of loc.eq.
408 if (fRecFileStatus<2) InitDataRecStorage(); // create a buffer to store the data
410 // write data of single measurement
411 if (lSigma<=0.0) { // If parameter is fixed, then no equation
412 for (int i=fNLocPar; i--;) derlc[i] = 0.0;
413 for (int i=fNGloParIni; i--;) dergb[i] = 0.0;
417 fRecord->AddResidual(lMeas);
419 // Retrieve local param interesting indices
420 for (int i=0;i<fNLocPar;i++) if (!IsZero(derlc[i])) {fRecord->AddIndexValue(i,derlc[i]); derlc[i] = 0.0;}
422 fRecord->AddWeight( 1.0/lSigma/lSigma );
424 // Idem for global parameters
425 for (int i=0;i<fNGloParIni;i++) if (!IsZero(dergb[i])) {
426 fRecord->AddIndexValue(i,dergb[i]); dergb[i] = 0.0;
427 int idrg = GetRGId(i);
428 fRecord->MarkGroup(idrg<0 ? -1 : fParamGrID[i]);
434 //_____________________________________________________________________________________________
435 void AliMillePede2::SetLocalEquation(int *indgb, double *dergb, int ngb, int *indlc,
436 double *derlc,int nlc,double lMeas,double lSigma)
438 // write data of single measurement. Note: the records ignore regrouping, store direct parameters
439 if (lSigma<=0.0) { // If parameter is fixed, then no equation
440 for (int i=nlc;i--;) derlc[i] = 0.0;
441 for (int i=ngb;i--;) dergb[i] = 0.0;
445 if (fRecFileStatus<2) InitDataRecStorage(); // create a buffer to store the data
447 fRecord->AddResidual(lMeas);
449 // Retrieve local param interesting indices
450 for (int i=0;i<nlc;i++) if (!IsZero(derlc[i])) {fRecord->AddIndexValue(indlc[i],derlc[i]); derlc[i]=0.; indlc[i]=0;}
452 fRecord->AddWeight( 1./lSigma/lSigma );
454 // Idem for global parameters
455 for (int i=0;i<ngb;i++) if (!IsZero(dergb[i])) {fRecord->AddIndexValue(indgb[i],dergb[i]); dergb[i]=0.; indgb[i]=0;}
460 //_____________________________________________________________________________________________
461 void AliMillePede2::SetGlobalConstraint(const double *dergb, double val, double sigma)
463 // Define a constraint equation.
464 if (!fConsRecFile || !fConsRecFile->IsWritable()) InitConsRecStorage(); // create a buffer to store the data
467 fRecord->AddResidual(val);
468 fRecord->AddWeight(sigma);
469 for (int i=0; i<fNGloParIni; i++) if (!IsZero(dergb[i])) fRecord->AddIndexValue(i,dergb[i]);
471 if (IsZero(sigma)) fNLagrangeConstraints++;
472 // printf("NewConstraint:\n"); fRecord->Print(); //RRR
473 SaveRecordConstraint();
477 //_____________________________________________________________________________________________
478 void AliMillePede2::SetGlobalConstraint(const int *indgb, const double *dergb, int ngb, double val,double sigma)
480 // Define a constraint equation.
481 if (!fConsRecFile || !fConsRecFile->IsWritable()) InitConsRecStorage(); // create a buffer to store the data
483 fRecord->AddResidual(val);
484 fRecord->AddWeight(sigma); // dummy
485 for (int i=0; i<ngb; i++) if (!IsZero(dergb[i])) fRecord->AddIndexValue(indgb[i],dergb[i]);
487 if (IsZero(sigma)) fNLagrangeConstraints++;
488 SaveRecordConstraint();
492 //_____________________________________________________________________________________________
493 Int_t AliMillePede2::LocalFit(double *localParams)
496 Perform local parameters fit once all the local equations have been set
497 -----------------------------------------------------------
498 localParams = (if !=0) will contain the fitted track parameters and
501 static int nrefSize = 0;
502 // static TArrayI refLoc,refGlo,nrefLoc,nrefGlo;
503 static Int_t *refLoc=0,*refGlo=0,*nrefLoc=0,*nrefGlo=0;
506 AliSymMatrix &matCLoc = *fMatCLoc;
507 AliMatrixSq &matCGlo = *fMatCGlo;
508 AliRectMatrix &matCGloLoc = *fMatCGloLoc;
510 memset(fVecBLoc,0,fNLocPar*sizeof(double));
514 int recSz = fRecord->GetSize();
516 while(cnt<recSz) { // Transfer the measurement records to matrices
518 // extract addresses of residual, weight and pointers on local and global derivatives for each point
519 if (nrefSize<=nPoints) {
521 nrefSize = 2*(nPoints+1);
522 tmpA = refLoc; refLoc = new Int_t[nrefSize]; if (tmpA) memcpy(refLoc,tmpA,nPoints*sizeof(int));
523 tmpA = refGlo; refGlo = new Int_t[nrefSize]; if (tmpA) memcpy(refGlo,tmpA,nPoints*sizeof(int));
524 tmpA = nrefLoc; nrefLoc = new Int_t[nrefSize]; if (tmpA) memcpy(nrefLoc,tmpA,nPoints*sizeof(int));
525 tmpA = nrefGlo; nrefGlo = new Int_t[nrefSize]; if (tmpA) memcpy(nrefGlo,tmpA,nPoints*sizeof(int));
528 refLoc[nPoints] = ++cnt;
530 while(!fRecord->IsWeight(cnt)) {nLoc++; cnt++;}
531 nrefLoc[nPoints] = nLoc;
533 refGlo[nPoints] = ++cnt;
535 while(!fRecord->IsResidual(cnt) && cnt<recSz) {nGlo++; cnt++;}
536 nrefGlo[nPoints] = nGlo;
542 double gloWgh = fRecord->GetWeight(); // global weight for this set
543 Int_t maxLocUsed = 0;
545 for (int ip=nPoints;ip--;) { // Transfer the measurement records to matrices
546 double resid = fRecord->GetValue( refLoc[ip]-1 );
547 double weight = fRecord->GetValue( refGlo[ip]-1 )*gloWgh;
548 double *derLoc = fRecord->GetValue()+refLoc[ip];
549 double *derGlo = fRecord->GetValue()+refGlo[ip];
550 int *indLoc = fRecord->GetIndex()+refLoc[ip];
551 int *indGlo = fRecord->GetIndex()+refGlo[ip];
553 for (int i=nrefGlo[ip];i--;) { // suppress the global part (only relevant with iterations)
555 // if regrouping was requested, do it here
557 int idtmp = fkReGroup[ indGlo[i] ];
558 if (idtmp == kFixParID) indGlo[i] = kFixParID; // fixed param in regrouping
561 int iID = indGlo[i]; // Global param indice
562 if (iID<0 || fSigmaPar[iID]<=0.) continue; // fixed parameter RRRCheck
563 if (fIsLinear[iID]) resid -= derGlo[i]*(fInitPar[iID]+fDeltaPar[iID]); // linear parameter
564 else resid -= derGlo[i]*fDeltaPar[iID]; // nonlinear parameter
567 // Symmetric matrix, don't bother j>i coeffs
568 for (int i=nrefLoc[ip];i--;) { // Fill local matrix and vector
569 fVecBLoc[ indLoc[i] ] += weight*resid*derLoc[i];
570 if (indLoc[i]>maxLocUsed) maxLocUsed = indLoc[i];
571 for (int j=i+1;j--;) matCLoc(indLoc[i] ,indLoc[j]) += weight*derLoc[i]*derLoc[j];
574 } // end of the transfer of the measurement record to matrices
576 matCLoc.SetSizeUsed(++maxLocUsed); // data with B=0 may use less than declared nLocals
580 printf("\nBefore\nLocalMatrix: "); matCLoc.Print("l");
581 printf("RHSLoc: "); for (int i=0;i<fNLocPar;i++) printf("%+e |",fVecBLoc[i]); printf("\n");
583 // first try to solve by faster Cholesky decomposition, then by Gaussian elimination
584 if (!matCLoc.SolveChol(fVecBLoc,kTRUE)) {
585 AliInfo("Failed to solve locals by Cholesky, trying Gaussian Elimination");
586 if (!matCLoc.SolveSpmInv(fVecBLoc,kTRUE)) {
587 AliInfo("Failed to solve locals by Gaussian Elimination, skip...");
589 return 0; // failed to solve
593 // If requested, store the track params and errors
594 //RRR printf("locfit: "); for (int i=0;i<fNLocPar;i++) printf("%+e |",fVecBLoc[i]); printf("\n");
596 if (localParams) for (int i=maxLocUsed; i--;) {
597 localParams[2*i] = fVecBLoc[i];
598 localParams[2*i+1] = TMath::Sqrt(TMath::Abs(matCLoc.QueryDiag(i)));
604 for (int ip=nPoints;ip--;) { // Calculate residuals
605 double resid = fRecord->GetValue( refLoc[ip]-1 );
606 double weight = fRecord->GetValue( refGlo[ip]-1 )*gloWgh;
607 double *derLoc = fRecord->GetValue()+refLoc[ip];
608 double *derGlo = fRecord->GetValue()+refGlo[ip];
609 int *indLoc = fRecord->GetIndex()+refLoc[ip];
610 int *indGlo = fRecord->GetIndex()+refGlo[ip];
612 // Suppress local and global contribution in residuals;
613 for (int i=nrefLoc[ip];i--;) resid -= derLoc[i]*fVecBLoc[ indLoc[i] ]; // local part
615 for (int i=nrefGlo[ip];i--;) { // global part
617 if ( iID<0 || fSigmaPar[iID] <= 0.) continue; // fixed parameter RRRCheck
618 if (fIsLinear[iID]) resid -= derGlo[i]*(fInitPar[iID]+fDeltaPar[iID]); // linear parameter
619 else resid -= derGlo[i]*fDeltaPar[iID]; // nonlinear parameter
622 // reject the track if the residual is too large (outlier)
623 double absres = TMath::Abs(resid);
624 if ( (absres >= fResCutInit && fIter ==1 ) ||
625 (absres >= fResCut && fIter > 1)) {
626 if (fLocFitAdd) fNLocFitsRejected++;
627 // printf("reject res %5ld %+e\n",fCurrRecDataID,resid);
631 lChi2 += weight*resid*resid ; // total chi^2
632 nEq++; // number of equations
633 } // end of Calculate residuals
636 int nDoF = nEq-maxLocUsed;
637 lChi2 = (nDoF>0) ? lChi2/nDoF : 0; // Chi^2/dof
639 if (fNStdDev != 0 && nDoF>0 && lChi2 > Chi2DoFLim(fNStdDev,nDoF)*fChi2CutFactor) { // check final chi2
640 if (fLocFitAdd) fNLocFitsRejected++;
641 // printf("reject chi2 %5ld: %+e\n",fCurrRecDataID, lChi2);
647 fNLocEquations += nEq;
651 fNLocEquations -= nEq;
654 // local operations are finished, track is accepted
655 // We now update the global parameters (other matrices)
659 for (int ip=nPoints;ip--;) { // Update matrices
660 double resid = fRecord->GetValue( refLoc[ip]-1 );
661 double weight = fRecord->GetValue( refGlo[ip]-1 )*gloWgh;
662 double *derLoc = fRecord->GetValue()+refLoc[ip];
663 double *derGlo = fRecord->GetValue()+refGlo[ip];
664 int *indLoc = fRecord->GetIndex()+refLoc[ip];
665 int *indGlo = fRecord->GetIndex()+refGlo[ip];
667 for (int i=nrefGlo[ip];i--;) { // suppress the global part
668 int iID = indGlo[i]; // Global param indice
669 if ( iID<0 || fSigmaPar[iID] <= 0.) continue; // fixed parameter RRRCheck
670 if (fIsLinear[iID]) resid -= derGlo[i]*(fInitPar[iID]+fDeltaPar[iID]); // linear parameter
671 else resid -= derGlo[i]*fDeltaPar[iID]; // nonlinear parameter
674 for (int ig=nrefGlo[ip];ig--;) {
675 int iIDg = indGlo[ig]; // Global param indice (the matrix line)
676 if ( iIDg<0 || fSigmaPar[iIDg] <= 0.) continue; // fixed parameter RRRCheck
677 if (fLocFitAdd) fVecBGlo[ iIDg ] += weight*resid*derGlo[ig]; //!!!
678 else fVecBGlo[ iIDg ] -= weight*resid*derGlo[ig]; //!!!
680 // First of all, the global/global terms (exactly like local matrix)
682 for (int jg=ig+1;jg--;) { // matCGlo is symmetric by construction
683 int jIDg = indGlo[jg];
684 if ( jIDg<0 || fSigmaPar[jIDg] <= 0.) continue; // fixed parameter RRRCheck
685 if ( !IsZero(vl = weight*derGlo[ig]*derGlo[jg]) ) {
686 fFillIndex[nfill] = jIDg;
687 fFillValue[nfill++] = fLocFitAdd ? vl:-vl;
690 if (nfill) matCGlo.AddToRow(iIDg,fFillValue,fFillIndex,nfill);
692 // Now we have also rectangular matrices containing global/local terms.
693 int iCIDg = fGlo2CGlo[iIDg]; // compressed Index of index
695 Double_t *rowGL = matCGloLoc(nGloInFit);
696 for (int k=maxLocUsed;k--;) rowGL[k] = 0.0; // reset the row
697 iCIDg = fGlo2CGlo[iIDg] = nGloInFit;
698 fCGlo2Glo[nGloInFit++] = iIDg;
701 Double_t *rowGLIDg = matCGloLoc(iCIDg);
702 for (int il=nrefLoc[ip];il--;) rowGLIDg[ indLoc[il] ] += weight*derGlo[ig]*derLoc[il];
703 fProcPnt[iIDg] += fLocFitAdd ? 1:-1; // update counter
706 } // end of Update matrices
709 printf("After GLO\n");
710 printf("MatCLoc: "); fMatCLoc->Print("l");
711 printf("MatCGlo: "); fMatCGlo->Print("l");
712 printf("MatCGlLc:"); fMatCGloLoc->Print("l");
713 printf("BGlo: "); for (int i=0; i<fNGloPar; i++) printf("%+e |",fVecBGlo[i]); printf("\n");
715 // calculate fMatCGlo -= fMatCGloLoc * fMatCLoc * fMatCGloLoc^T
716 // and fVecBGlo -= fMatCGloLoc * fVecBLoc
718 //-------------------------------------------------------------- >>>
720 for (int iCIDg=0; iCIDg<nGloInFit; iCIDg++) {
721 int iIDg = fCGlo2Glo[iCIDg];
724 Double_t *rowGLIDg = matCGloLoc(iCIDg);
725 for (int kl=0;kl<maxLocUsed;kl++) if (rowGLIDg[kl]) vl += rowGLIDg[kl]*fVecBLoc[kl];
726 if (!IsZero(vl)) fVecBGlo[iIDg] -= fLocFitAdd ? vl : -vl;
729 for (int jCIDg=0;jCIDg<=iCIDg; jCIDg++) {
730 int jIDg = fCGlo2Glo[jCIDg];
733 Double_t *rowGLJDg = matCGloLoc(jCIDg);
734 for (int kl=0;kl<maxLocUsed;kl++) {
736 if ( (!IsZero(vll=rowGLIDg[kl]*rowGLJDg[kl])) ) vl += matCLoc.QueryDiag(kl)*vll;
739 for (int ll=0;ll<kl;ll++) {
740 if ( !IsZero(vll=rowGLIDg[kl]*rowGLJDg[ll]) ) vl += matCLoc(kl,ll)*vll;
741 if ( !IsZero(vll=rowGLIDg[ll]*rowGLJDg[kl]) ) vl += matCLoc(kl,ll)*vll;
745 fFillIndex[nfill] = jIDg;
746 fFillValue[nfill++] = fLocFitAdd ? -vl : vl;
749 if (nfill) matCGlo.AddToRow(iIDg,fFillValue,fFillIndex,nfill);
752 // reset compressed index array
755 printf("After GLOLoc\n");
756 printf("MatCGlo: "); fMatCGlo->Print("");
757 printf("BGlo: "); for (int i=0; i<fNGloPar; i++) printf("%+e |",fVecBGlo[i]); printf("\n");
759 for (int i=nGloInFit;i--;) {
760 fGlo2CGlo[ fCGlo2Glo[i] ] = -1;
764 //---------------------------------------------------- <<<
768 //_____________________________________________________________________________________________
769 Int_t AliMillePede2::GlobalFit(Double_t *par, Double_t *error, Double_t *pull)
771 // performs a requested number of global iterations
774 TStopwatch sw; sw.Start();
777 AliInfo("Starting Global fit.");
778 while (fIter<=fMaxIter) {
780 res = GlobalFitIteration();
783 if (!IsZero(fChi2CutFactor-fChi2CutRef)) {
784 fChi2CutFactor = TMath::Sqrt(fChi2CutFactor);
785 if (fChi2CutFactor < 1.2*fChi2CutRef) {
786 fChi2CutFactor = fChi2CutRef;
787 //RRR fIter = fMaxIter - 1; // Last iteration
794 AliInfo(Form("Global fit %s, CPU time: %.1f",res ? "Converged":"Failed",sw.CpuTime()));
797 if (par) for (int i=fNGloParIni;i--;) par[i] = GetFinalParam(i);
799 if (fGloSolveStatus==kInvert) { // errors on params are available
800 if (error) for (int i=fNGloParIni;i--;) error[i] = GetFinalError(i);
801 if (pull) for (int i=fNGloParIni;i--;) pull[i] = GetPull(i);
807 //_____________________________________________________________________________________________
808 Int_t AliMillePede2::GlobalFitIteration()
810 // perform global parameters fit once all the local equations have been fitted
812 AliInfo(Form("Global Fit Iteration#%2d (Local Fit Chi^2 cut factor: %.2f)",fIter,fChi2CutFactor));
814 if (!fNGloPar || !fTreeData) {
815 AliInfo("No data was stored, stopping iteration");
823 fConstrUsed = new Bool_t[fNGloConstraints];
824 memset(fConstrUsed,0,fNGloConstraints*sizeof(Bool_t));
826 // Reset all info specific for this step
827 AliMatrixSq& matCGlo = *fMatCGlo;
829 memset(fProcPnt,0,fNGloPar*sizeof(Int_t));
831 fNGloConstraints = fTreeConstr ? fTreeConstr->GetEntries() : 0;
833 // count number of Lagrange constraints: they need new row/cols to be added
834 fNLagrangeConstraints = 0;
835 for (int i=0; i<fNGloConstraints; i++) {
836 ReadRecordConstraint(i);
837 if ( IsZero(fRecord->GetValue(1)) ) fNLagrangeConstraints++; // exact constraint (no error) -> Lagrange multiplier
840 // if needed, readjust the size of the global vector (for matrices this is done automatically)
841 if (!fVecBGlo || fNGloSize!=fNGloPar+fNLagrangeConstraints) {
842 delete[] fVecBGlo; // in case some constraint was added between the two manual iterations
843 fNGloSize = fNGloPar+fNLagrangeConstraints;
844 fVecBGlo = new Double_t[fNGloSize];
846 memset(fVecBGlo,0,fNGloSize*sizeof(double));
849 fNLocFitsRejected = 0;
852 // Process data records and build the matrices
853 Long_t ndr = fTreeData->GetEntries();
854 Long_t first = fSelFirst>0 ? fSelFirst : 0;
855 Long_t last = fSelLast<1 ? ndr : (fSelLast>=ndr ? ndr : fSelLast+Long_t(1));
858 AliInfo(Form("Building the Global matrix from data records %ld : %ld",first,last));
861 TStopwatch swt; swt.Start();
862 fLocFitAdd = kTRUE; // add contributions of matching tracks
863 for (Long_t i=0;i<ndr;i++) {
864 Long_t iev = i+first;
866 if (!IsRecordAcceptable()) continue;
868 if ( (i%int(0.2*ndr)) == 0) printf("%.1f%% of local fits done\n", double(100.*i)/ndr);
871 printf("%ld local fits done: ", ndr);
873 printf("MatCGlo: "); fMatCGlo->Print("l");
874 printf("BGlo: "); for (int i=0; i<fNGloPar; i++) printf("%+e |",fVecBGlo[i]); printf("\n");
880 // ---------------------- Reject parameters with low statistics ------------>>
882 if (fMinPntValid>1 && fNGroupsSet) {
884 printf("Checking parameters with statistics < %d\n",fMinPntValid);
887 // 1) build the list of parameters to fix
888 Int_t fixArrSize = 10;
889 Int_t nFixedGroups = 0;
890 TArrayI fixGroups(fixArrSize);
894 double oldMin = 1.e20;
895 double oldMax =-1.e20;
897 for (int i=fNGloPar;i--;) { // // Reset row and column of fixed params and add 1/sig^2 to free ones
898 int grID = fParamGrID[i];
899 if (grID<0) continue; // not in the group
901 if (grID!=grIDold) { // starting new group
902 if (grIDold>=0) { // decide if the group has enough statistics
903 if (oldMin<fMinPntValid && oldMax<2*fMinPntValid) { // suppress group
904 for (int iold=oldStart;iold>i;iold--) fProcPnt[iold] = 0;
905 Bool_t fnd = kFALSE; // check if the group is already accounted
906 for (int j=nFixedGroups;j--;) if (fixGroups[j]==grIDold) {fnd=kTRUE; break;}
908 if (nFixedGroups>=fixArrSize) {fixArrSize*=2; fixGroups.Set(fixArrSize);}
909 fixGroups[nFixedGroups++] = grIDold; // add group to fix
913 grIDold = grID; // mark the start of the new group
918 if (oldMin>fProcPnt[i]) oldMin = fProcPnt[i];
919 if (oldMax<fProcPnt[i]) oldMax = fProcPnt[i];
922 // extra check for the last group
923 if (grIDold>=0 && oldMin<fMinPntValid && oldMax<2*fMinPntValid) { // suppress group
924 for (int iold=oldStart;iold--;) fProcPnt[iold] = 0;
925 Bool_t fnd = kFALSE; // check if the group is already accounted
926 for (int j=nFixedGroups;j--;) if (fixGroups[j]==grIDold) {fnd=kTRUE; break;}
928 if (nFixedGroups>=fixArrSize) {fixArrSize*=2; fixGroups.Set(fixArrSize);}
929 fixGroups[nFixedGroups++] = grIDold; // add group to fix
933 // 2) loop over records and add contributions of fixed groups with negative sign
936 for (Long_t i=0;i<ndr;i++) {
937 Long_t iev = i+first;
939 if (!IsRecordAcceptable()) continue;
940 Bool_t suppr = kFALSE;
941 for (int ifx=nFixedGroups;ifx--;)if (fRecord->IsGroupPresent(fixGroups[ifx])) suppr = kTRUE;
942 if (suppr) LocalFit();
947 printf("Suppressed contributions of groups with NPoints<%d :\n",fMinPntValid);
948 for (int i=0;i<nFixedGroups;i++) printf("%d ",fixGroups[i]); printf("\n");
953 // ---------------------- Reject parameters with low statistics ------------<<
955 // add large number to diagonal of fixed params
957 for (int i=fNGloPar;i--;) { // // Reset row and column of fixed params and add 1/sig^2 to free ones
958 // printf("#%3d : Nproc : %5d grp: %d\n",i,fProcPnt[i],fParamGrID[i]);
962 matCGlo.DiagElem(i) = 1.;//float(fNLocEquations*fNLocEquations);
963 // matCGlo.DiagElem(i) = float(fNLocEquations*fNLocEquations);
965 else matCGlo.DiagElem(i) += (fgWeightSigma ? fProcPnt[i] : 1.)/(fSigmaPar[i]*fSigmaPar[i]);
968 for (int i=fNGloPar;i--;) fDiagCGlo[i] = matCGlo.QueryDiag(i); // save the diagonal elements
970 // add constraint equations
971 int nVar = fNGloPar; // Current size of global matrix
972 for (int i=0; i<fNGloConstraints; i++) {
973 ReadRecordConstraint(i);
974 double val = fRecord->GetValue(0);
975 double sig = fRecord->GetValue(1);
976 int *indV = fRecord->GetIndex()+2;
977 double *der = fRecord->GetValue()+2;
978 int csize = fRecord->GetSize()-2;
981 for (int jp=csize;jp--;) {
983 if (fkReGroup[idp]<0) AliFatal(Form("Constain is requested for suppressed parameter #%d",indV[jp]));
987 // check if after suppression of fixed variables there are non-0 derivatives
988 // and determine the max statistics of involved params
991 for (int j=csize;j--;) {
992 if (fProcPnt[indV[j]]<1) nSuppressed++;
994 maxStat = TMath::Max(maxStat,fProcPnt[indV[j]]);
998 if (nSuppressed==csize) {
999 // AliInfo(Form("Neglecting constraint %d of %d derivatives since no free parameters left",i,csize));
1001 // was this constraint ever created ?
1002 if ( sig==0 && fConstrUsed[i] ) { // this is needed only for constraints with Lagrange multiplier
1003 // to avoid empty row impose dummy constraint on "Lagrange multiplier"
1004 matCGlo.DiagElem(nVar) = 1.;
1005 fVecBGlo[nVar++] = 0;
1010 // account for already accumulated corrections
1011 for (int j=csize; j--;) val -= der[j]*(fInitPar[ indV[j] ]+fDeltaPar[ indV[j] ]);
1013 if (sig > 0) { // this is a gaussian constriant: no Lagrange multipliers are added
1015 double sig2i = (fgWeightSigma ? TMath::Sqrt(maxStat) : 1.)/sig/sig;
1016 for (int ir=0;ir<csize;ir++) {
1018 for (int ic=0;ic<=ir;ic++) { // matrix is symmetric
1020 double vl = der[ir]*der[ic]*sig2i;
1021 if (!IsZero(vl)) matCGlo(iID,jID) += vl;
1023 fVecBGlo[iID] += val*der[ir]*sig2i;
1026 else { // this is exact constriant: Lagrange multipliers must be added
1027 for (int j=csize; j--;) {
1029 if (fProcPnt[jID]<1) continue; // this parameter was fixed, don't put it into constraint
1030 matCGlo(nVar,jID) = float(fNLocEquations)*der[j]; // fMatCGlo is symmetric, only lower triangle is filled
1033 if (matCGlo.QueryDiag(nVar)) matCGlo.DiagElem(nVar) = 0.0;
1034 fVecBGlo[nVar++] = float(fNLocEquations)*val; //RS ? should we use here fNLocFits ?
1035 fConstrUsed[i] = kTRUE;
1039 AliInfo(Form("Obtained %-7ld equations from %-7ld records (%-7ld rejected). Fixed %-4d globals",
1040 fNLocEquations,fNLocFits,fNLocFitsRejected,fNGloFix));
1046 printf("Solving:\n");
1048 for (int i=0;i<fNGloSize;i++) printf("b%2d : %+e\n",i,fVecBGlo[i]);
1050 fGloSolveStatus = SolveGlobalMatEq(); // obtain solution for this step
1052 printf("Solve %d |",fIter); sws.Print();
1055 AliInfo(Form("Iteration#%2d %s. CPU time: %.1f",fIter,fGloSolveStatus==kFailed ? "Failed":"Converged",sw.CpuTime()));
1056 if (fGloSolveStatus==kFailed) return 0;
1058 for (int i=fNGloPar;i--;) fDeltaPar[i] += fVecBGlo[i]; // Update global parameters values (for iterations)
1060 // PrintGlobalParameters();
1064 //_____________________________________________________________________________________________
1065 Int_t AliMillePede2::SolveGlobalMatEq()
1068 // solve global matrix equation MatCGlob*X=VecBGlo and store the result in the VecBGlo
1071 printf("GlobalMatrix\n");
1074 for (int i=0;i<fNGloPar;i++) printf("%d %+e\n",i,fVecBGlo[i]);
1077 if (!fgIsMatGloSparse) {
1079 if (fNLagrangeConstraints==0) { // pos-def systems are faster to solve by Cholesky
1080 if ( ((AliSymMatrix*)fMatCGlo)->SolveChol(fVecBGlo, fgInvChol) ) return fgInvChol ? kInvert:kNoInversion;
1081 else AliInfo("Solution of Global Dense System by Cholesky failed, trying Gaussian Elimiation");
1084 if (((AliSymMatrix*)fMatCGlo)->SolveSpmInv(fVecBGlo, kTRUE)) return kInvert;
1085 else AliInfo("Solution of Global Dense System by Gaussian Elimination failed, trying iterative methods");
1087 // try to solve by minres
1088 TVectorD sol(fNGloSize);
1090 AliMinResSolve *slv = new AliMinResSolve(fMatCGlo,fVecBGlo);
1091 if (!slv) return kFailed;
1093 Bool_t res = kFALSE;
1094 if (fgIterSol == AliMinResSolve::kSolMinRes)
1095 res = slv->SolveMinRes(sol,fgMinResCondType,fgMinResMaxIter,fgMinResTol);
1096 else if (fgIterSol == AliMinResSolve::kSolFGMRes)
1097 res = slv->SolveFGMRES(sol,fgMinResCondType,fgMinResMaxIter,fgMinResTol,fgNKrylovV);
1099 AliInfo(Form("Undefined Iteritive Solver ID=%d, only %d are defined",fgIterSol,AliMinResSolve::kNSolvers));
1102 const char* faildump = "fgmr_failed.dat";
1103 int defout = dup(1);
1105 AliInfo("Failed on dup");
1108 int slvDump = open(faildump, O_RDWR|O_CREAT, 0666);
1112 printf("#Failed to solve using solver %d with PreCond: %d MaxIter: %d Tol: %e NKrylov: %d\n",
1113 fgIterSol,fgMinResCondType,fgMinResMaxIter,fgMinResTol,fgNKrylovV);
1114 printf("#Dump of matrix:\n");
1115 fMatCGlo->Print("10");
1116 printf("#Dump of RHS:\n");
1117 for (int i=0;i<fNGloSize;i++) printf("%d %+.10f\n",i,fVecBGlo[i]);
1121 printf("#Dumped failed matrix and RHS to %s\n",faildump);
1123 else AliInfo("Failed on file open for matrix dumping");
1127 for (int i=fNGloSize;i--;) fVecBGlo[i] = sol[i];
1128 return kNoInversion;
1132 //_____________________________________________________________________________________________
1133 Float_t AliMillePede2::Chi2DoFLim(int nSig, int nDoF) const
1135 /// return the limit in chi^2/nd for n sigmas stdev authorized
1136 // Only n=1, 2, and 3 are expected in input
1138 float sn[3] = {0.47523, 1.690140, 2.782170};
1139 float table[3][30] = {{1.0000, 1.1479, 1.1753, 1.1798, 1.1775, 1.1730, 1.1680, 1.1630,
1140 1.1581, 1.1536, 1.1493, 1.1454, 1.1417, 1.1383, 1.1351, 1.1321,
1141 1.1293, 1.1266, 1.1242, 1.1218, 1.1196, 1.1175, 1.1155, 1.1136,
1142 1.1119, 1.1101, 1.1085, 1.1070, 1.1055, 1.1040},
1143 {4.0000, 3.0900, 2.6750, 2.4290, 2.2628, 2.1415, 2.0481, 1.9736,
1144 1.9124, 1.8610, 1.8171, 1.7791, 1.7457, 1.7161, 1.6897, 1.6658,
1145 1.6442, 1.6246, 1.6065, 1.5899, 1.5745, 1.5603, 1.5470, 1.5346,
1146 1.5230, 1.5120, 1.5017, 1.4920, 1.4829, 1.4742},
1147 {9.0000, 5.9146, 4.7184, 4.0628, 3.6410, 3.3436, 3.1209, 2.9468,
1148 2.8063, 2.6902, 2.5922, 2.5082, 2.4352, 2.3711, 2.3143, 2.2635,
1149 2.2178, 2.1764, 2.1386, 2.1040, 2.0722, 2.0428, 2.0155, 1.9901,
1150 1.9665, 1.9443, 1.9235, 1.9040, 1.8855, 1.8681}};
1156 lNSig = TMath::Max(1,TMath::Min(nSig,3));
1159 return table[lNSig-1][nDoF-1];
1161 else { // approximation
1162 return ((sn[lNSig-1]+TMath::Sqrt(float(2*nDoF-3)))*
1163 (sn[lNSig-1]+TMath::Sqrt(float(2*nDoF-3))))/float(2*nDoF-2);
1168 //_____________________________________________________________________________________________
1169 Int_t AliMillePede2::SetIterations(double lChi2CutFac)
1171 // Number of iterations is calculated from lChi2CutFac
1172 fChi2CutFactor = TMath::Max(1.0, lChi2CutFac);
1174 AliInfo(Form("Initial cut factor is %f",fChi2CutFactor));
1175 fIter = 1; // Initializes the iteration process
1179 //_____________________________________________________________________________________________
1180 Double_t AliMillePede2::GetParError(int iPar) const
1182 // return error for parameter iPar
1183 if (fGloSolveStatus==kInvert) {
1184 if (fkReGroup) iPar = fkReGroup[iPar];
1185 if (iPar<0) {AliInfo(Form("Parameter %d was suppressed in the regrouping",iPar)); return 0;}
1186 double res = fMatCGlo->QueryDiag(iPar);
1187 if (res>=0) return TMath::Sqrt(res);
1192 //_____________________________________________________________________________________________
1193 Double_t AliMillePede2::GetPull(int iPar) const
1195 // return pull for parameter iPar
1196 if (fGloSolveStatus==kInvert) {
1197 if (fkReGroup) iPar = fkReGroup[iPar];
1198 if (iPar<0) {AliInfo(Form("Parameter %d was suppressed in the regrouping",iPar)); return 0;}
1200 return fProcPnt[iPar]>0 && (fSigmaPar[iPar]*fSigmaPar[iPar]-fMatCGlo->QueryDiag(iPar))>0. && fSigmaPar[iPar]>0
1201 ? fDeltaPar[iPar]/TMath::Sqrt(fSigmaPar[iPar]*fSigmaPar[iPar]-fMatCGlo->QueryDiag(iPar)) : 0;
1207 //_____________________________________________________________________________________________
1208 Int_t AliMillePede2::PrintGlobalParameters() const
1210 /// Print the final results into the logfile
1212 double lGlobalCor =0.;
1215 AliInfo(" Result of fit for global parameters");
1216 AliInfo(" ===================================");
1217 AliInfo(" I initial final differ lastcor error gcor Npnt");
1218 AliInfo("----------------------------------------------------------------------------------------------");
1220 for (int i0=0; i0<fNGloParIni; i0++) {
1221 int i = GetRGId(i0); if (i<0) continue;
1222 lError = GetParError(i);
1226 if (fGloSolveStatus==kInvert && TMath::Abs( (dg=fMatCGlo->QueryDiag(i)) *fDiagCGlo[i]) > 0) {
1227 lGlobalCor = TMath::Sqrt(TMath::Abs(1.0-1.0/(dg*fDiagCGlo[i])));
1228 AliInfo(Form("%d\t %.6f\t %.6f\t %.6f\t %.6f\t %.6f\t %.6f\t %6d",
1229 i,fInitPar[i],fInitPar[i]+fDeltaPar[i],fDeltaPar[i],fVecBGlo[i],lError,lGlobalCor,fProcPnt[i]));
1232 AliInfo(Form("%d\t %.6f\t %.6f\t %.6f\t %.6f\t OFF\t OFF\t %6d",i,fInitPar[i],fInitPar[i]+fDeltaPar[i],
1233 fDeltaPar[i],fVecBGlo[i],fProcPnt[i]));
1239 //_____________________________________________________________________________________________
1240 Bool_t AliMillePede2::IsRecordAcceptable() const
1242 // validate record according run lists set by the user
1243 static Long_t prevRunID = kMaxInt;
1244 static Bool_t prevAns = kTRUE;
1245 Long_t runID = fRecord->GetRunID();
1246 if (runID!=prevRunID) {
1249 // is run to be rejected?
1250 if (fRejRunList && (n=fRejRunList->GetSize())) {
1252 for (int i=n;i--;) if (runID == (*fRejRunList)[i]) {
1254 printf("New Run to reject: %ld -> %d\n",runID,prevAns);
1258 else if (fAccRunList && (n=fAccRunList->GetSize())) { // is run specifically selected
1260 for (int i=n;i--;) if (runID == (*fAccRunList)[i]) {prevAns = kTRUE; break;}
1268 //_____________________________________________________________________________________________
1269 void AliMillePede2::SetRejRunList(const UInt_t *runs, Int_t nruns)
1271 // set the list of runs to be rejected
1272 if (fRejRunList) delete fRejRunList;
1274 if (nruns<1 || !runs) return;
1275 fRejRunList = new TArrayL(nruns);
1276 for (int i=0;i<nruns;i++) (*fRejRunList)[i] = runs[i];
1279 //_____________________________________________________________________________________________
1280 void AliMillePede2::SetAccRunList(const UInt_t *runs, Int_t nruns)
1282 // set the list of runs to be selected
1283 if (fAccRunList) delete fAccRunList;
1285 if (nruns<1 || !runs) return;
1286 fAccRunList = new TArrayL(nruns);
1287 for (int i=0;i<nruns;i++) (*fAccRunList)[i] = runs[i];
1290 //_____________________________________________________________________________________________
1291 void AliMillePede2::SetInitPars(const Double_t* par)
1293 // initialize parameters, account for eventual grouping
1294 for (int i=0;i<fNGloParIni;i++) {
1295 int id = GetRGId(i); if (id<0) continue;
1296 fInitPar[id] = par[i];
1300 //_____________________________________________________________________________________________
1301 void AliMillePede2::SetSigmaPars(const Double_t* par)
1303 // initialize sigmas, account for eventual grouping
1304 for (int i=0;i<fNGloParIni;i++) {
1305 int id = GetRGId(i); if (id<0) continue;
1306 fSigmaPar[id] = par[i];
1310 //_____________________________________________________________________________________________
1311 void AliMillePede2::SetInitPar(Int_t i,Double_t par)
1313 // initialize param, account for eventual grouping
1314 int id = GetRGId(i); if (id<0) return;
1318 //_____________________________________________________________________________________________
1319 void AliMillePede2::SetSigmaPar(Int_t i,Double_t par)
1321 // initialize sigma, account for eventual grouping
1322 int id = GetRGId(i); if (id<0) return;
1323 fSigmaPar[id] = par;