]> git.uio.no Git - u/mrichter/AliRoot.git/blame - STEER/AliMillePede2.cxx
AliACORDEQAChecker::Check fixed
[u/mrichter/AliRoot.git] / STEER / AliMillePede2.cxx
CommitLineData
de34b538 1/**********************************************************************************************/
2/* General class for alignment with large number of degrees of freedom */
3/* Based on the original milliped2 by Volker Blobel */
4/* http://www.desy.de/~blobel/mptalks.html */
5/* */
6/* Author: ruben.shahoyan@cern.ch */
7/* */
8/**********************************************************************************************/
9
8a9ab0eb 10#include "AliMillePede2.h"
11#include "AliLog.h"
12#include <TStopwatch.h>
de34b538 13#include <TFile.h>
14#include <TMath.h>
7c3070ec 15#include <TVectorD.h>
de34b538 16#include "AliMatrixSq.h"
17#include "AliSymMatrix.h"
18#include "AliRectMatrix.h"
19#include "AliMatrixSparse.h"
8a9ab0eb 20
8a9ab0eb 21
de34b538 22using namespace std;
8a9ab0eb 23
8a9ab0eb 24
8a9ab0eb 25ClassImp(AliMillePede2)
26
de34b538 27Bool_t AliMillePede2::fgInvChol = kTRUE; // Invert global matrix with Cholesky solver
28Bool_t AliMillePede2::fgWeightSigma = kTRUE; // weight local constraint by module statistics
29Bool_t AliMillePede2::fgIsMatGloSparse = kFALSE; // use faster dense matrix by default
30Int_t AliMillePede2::fgMinResCondType = 1; // Jacoby preconditioner by default
31Double_t AliMillePede2::fgMinResTol = 1.e-11; // default tolerance
32Int_t AliMillePede2::fgMinResMaxIter = 10000; // default max number of iterations
8a9ab0eb 33Int_t AliMillePede2::fgIterSol = AliMinResSolve::kSolMinRes; // default iterative solver
de34b538 34Int_t AliMillePede2::fgNKrylovV = 240; // default number of Krylov vectors to keep
8a9ab0eb 35
36//_____________________________________________________________________________________________
37AliMillePede2::AliMillePede2()
38: fNLocPar(0),
39 fNGloPar(0),
40 fNGloSize(0),
41//
42 fNLocEquations(0),
43 fIter(0),
44 fMaxIter(10),
45 fNStdDev(3),
46 fNGloConstraints(0),
de34b538 47 fNLagrangeConstraints(0),
8a9ab0eb 48 fNLocFits(0),
49 fNLocFitsRejected(0),
50 fNGloFix(0),
51 fGloSolveStatus(gkFailed),
52//
53 fChi2CutFactor(1.),
54 fChi2CutRef(1.),
55 fResCutInit(100.),
56 fResCut(100.),
de34b538 57 fMinPntValid(1),
8a9ab0eb 58//
de34b538 59 fNGroupsSet(0),
60 fParamGrID(0),
8a9ab0eb 61 fProcPnt(0),
62 fVecBLoc(0),
63 fDiagCGlo(0),
64 fVecBGlo(0),
65 fInitPar(0),
66 fDeltaPar(0),
67 fSigmaPar(0),
68 fIsLinear(0),
69 fConstrUsed(0),
70//
71 fGlo2CGlo(0),
72 fCGlo2Glo(0),
73//
74 fMatCLoc(0),
75 fMatCGlo(0),
76 fMatCGloLoc(0),
de34b538 77 //
78 fFillIndex(0),
79 fFillValue(0),
80 //
8a9ab0eb 81 fDataRecFName("/tmp/mp2_data_records.root"),
82 fRecord(0),
83 fDataRecFile(0),
84 fTreeData(0),
85 //
86 fConstrRecFName("/tmp/mp2_constraints_records.root"),
87 fTreeConstr(0),
88 fConsRecFile(0),
89 fCurrRecDataID(0),
de34b538 90 fCurrRecConstrID(0),
91 fLocFitAdd(kTRUE)
8a9ab0eb 92{}
93
94//_____________________________________________________________________________________________
95AliMillePede2::AliMillePede2(const AliMillePede2& src) :
96 TObject(src),fNLocPar(0),fNGloPar(0),fNGloSize(0),fNLocEquations(0),fIter(0),
de34b538 97 fMaxIter(10),fNStdDev(3),fNGloConstraints(0),fNLagrangeConstraints(0),
98 fNLocFits(0),fNLocFitsRejected(0),
8a9ab0eb 99 fNGloFix(0),fGloSolveStatus(0),fChi2CutFactor(0),fChi2CutRef(0),fResCutInit(0),
de34b538 100 fResCut(0),fMinPntValid(1),fNGroupsSet(0),fParamGrID(0),fProcPnt(0),fVecBLoc(0),fDiagCGlo(0),fVecBGlo(0),
8a9ab0eb 101 fInitPar(0),fDeltaPar(0),fSigmaPar(0),fIsLinear(0),fConstrUsed(0),fGlo2CGlo(0),fCGlo2Glo(0),
de34b538 102 fMatCLoc(0),fMatCGlo(0),fMatCGloLoc(0),fFillIndex(0),fFillValue(0),
103 fDataRecFName(0),fRecord(0),fDataRecFile(0),
104 fTreeData(0),fConstrRecFName(0),fTreeConstr(0),fConsRecFile(0),fCurrRecDataID(0),
105 fCurrRecConstrID(0),fLocFitAdd(kTRUE)
8a9ab0eb 106{printf("Dummy\n");}
107
108//_____________________________________________________________________________________________
109AliMillePede2::~AliMillePede2()
110{
111 CloseDataRecStorage();
112 CloseConsRecStorage();
113 //
de34b538 114 delete[] fParamGrID;
8a9ab0eb 115 delete[] fProcPnt;
116 delete[] fVecBLoc;
117 delete[] fDiagCGlo;
118 delete[] fVecBGlo;
119 delete[] fInitPar;
120 delete[] fDeltaPar;
121 delete[] fSigmaPar;
122 delete[] fGlo2CGlo;
123 delete[] fCGlo2Glo;
124 delete[] fIsLinear;
125 delete[] fConstrUsed;
de34b538 126 delete[] fFillIndex;
127 delete[] fFillValue;
8a9ab0eb 128 //
129 delete fRecord;
130 delete fMatCLoc;
131 delete fMatCGlo;
132 delete fMatCGloLoc;
133}
134
135//_____________________________________________________________________________________________
136Int_t AliMillePede2::InitMille(int nGlo, int nLoc, int lNStdDev,double lResCut, double lResCutInit)
137{
138 //
139 if (nLoc>0) fNLocPar = nLoc;
140 if (nGlo>0) fNGloPar = nGlo;
141 if (lResCutInit>0) fResCutInit = lResCutInit;
142 if (lResCut>0) fResCut = lResCut;
143 if (lNStdDev>0) fNStdDev = lNStdDev;
144 //
145 fNGloSize = fNGloPar;
146 //
147 try {
148 //
149 if (fgIsMatGloSparse) {fMatCGlo = new AliMatrixSparse(fNGloPar); fMatCGlo->SetSymmetric(kTRUE);}
150 else fMatCGlo = new AliSymMatrix(fNGloPar);
151 //
de34b538 152 fFillIndex = new Int_t[fNGloPar];
153 fFillValue = new Double_t[fNGloPar];
154 //
8a9ab0eb 155 fMatCLoc = new AliSymMatrix(fNLocPar);
de34b538 156 fMatCGloLoc = new AliRectMatrix(fNGloPar,fNLocPar);
8a9ab0eb 157 //
de34b538 158 fParamGrID = new Int_t[fNGloPar];
8a9ab0eb 159 fProcPnt = new Int_t[fNGloPar];
160 fVecBLoc = new Double_t[fNLocPar];
161 fDiagCGlo = new Double_t[fNGloPar];
162 //
163 fInitPar = new Double_t[fNGloPar];
164 fDeltaPar = new Double_t[fNGloPar];
165 fSigmaPar = new Double_t[fNGloPar];
166 fIsLinear = new Bool_t[fNGloPar];
167 //
168 fGlo2CGlo = new Int_t[fNGloPar];
169 fCGlo2Glo = new Int_t[fNGloPar];
170 }
171 catch(bad_alloc&) {
172 AliInfo(Form("Failed to allocate the memory for %d global and %d local parameters",fNGloPar,fNLocPar));
173 return 0;
174 }
175 //
176 memset(fVecBLoc ,0,fNLocPar*sizeof(Double_t));
177 memset(fDiagCGlo ,0,fNGloPar*sizeof(Double_t));
178 memset(fInitPar ,0,fNGloPar*sizeof(Double_t));
179 memset(fDeltaPar ,0,fNGloPar*sizeof(Double_t));
180 memset(fSigmaPar ,0,fNGloPar*sizeof(Double_t));
181 memset(fProcPnt ,0,fNGloPar*sizeof(Int_t));
182 //
183 for (int i=fNGloPar;i--;) {
de34b538 184 fGlo2CGlo[i] = fCGlo2Glo[i] = -1;
8a9ab0eb 185 fIsLinear[i] = kTRUE;
de34b538 186 fParamGrID[i] = -1;
8a9ab0eb 187 }
188 //
189 return 1;
190}
191
8a9ab0eb 192//_____________________________________________________________________________________________
193Bool_t AliMillePede2::ImposeDataRecFile(const char* fname)
194{
195 CloseDataRecStorage();
196 SetDataRecFName(fname);
197 return InitDataRecStorage(kTRUE); // open in read mode
198}
199
200//_____________________________________________________________________________________________
201Bool_t AliMillePede2::ImposeConsRecFile(const char* fname)
202{
203 CloseConsRecStorage();
204 SetConsRecFName(fname);
205 return InitConsRecStorage(kTRUE); // open in read mode
206}
207
208//_____________________________________________________________________________________________
209Bool_t AliMillePede2::InitDataRecStorage(Bool_t read)
210{
211 // initialize the buffer for processed measurements records
212 //
213 if (fDataRecFile) {AliInfo("Data Records File is already initialized"); return kFALSE;}
214 //
215 if (!fRecord) fRecord = new AliMillePedeRecord();
216 //
217 fDataRecFile = TFile::Open(GetDataRecFName(),read ? "":"recreate");
218 if (!fDataRecFile) {AliInfo(Form("Failed to initialize data records file %s",GetDataRecFName())); return kFALSE;}
219 //
220 AliInfo(Form("File %s used for derivatives records",GetDataRecFName()));
221 if (read) {
222 fTreeData = (TTree*)fDataRecFile->Get("AliMillePedeRecords_Data");
223 if (!fTreeData) {AliInfo(Form("Did not find data records tree in %s",GetDataRecFName())); return kFALSE;}
224 fTreeData->SetBranchAddress("Record_Data",&fRecord);
225 AliInfo(Form("Found %d derivatives records",fTreeData->GetEntries()));
226 }
227 else {
228 fTreeData = new TTree("AliMillePedeRecords_Data","Data Records for AliMillePede2");
229 fTreeData->Branch("Record_Data","AliMillePedeRecord",&fRecord,32000,99);
230 }
231 fCurrRecDataID = -1;
232 //
233 return kTRUE;
234}
235
236//_____________________________________________________________________________________________
237Bool_t AliMillePede2::InitConsRecStorage(Bool_t read)
238{
239 // initialize the buffer for processed measurements records
240 //
241 if (fConsRecFile) {AliInfo("Constraints Records File is already initialized"); return kFALSE;}
242 //
243 if (!fRecord) fRecord = new AliMillePedeRecord();
244 //
245 fConsRecFile = TFile::Open(GetConsRecFName(),read ? "":"recreate");
246 if (!fConsRecFile) {AliInfo(Form("Failed to initialize constraints records file %s",GetConsRecFName())); return kFALSE;}
247 //
248 AliInfo(Form("File %s used for constraints records",GetConsRecFName()));
249 if (read) {
250 fTreeConstr = (TTree*)fConsRecFile->Get("AliMillePedeRecords_Constraints");
251 if (!fTreeConstr) {AliInfo(Form("Did not find constraints records tree in %s",GetConsRecFName())); return kFALSE;}
252 fTreeConstr->SetBranchAddress("Record_Constraints",&fRecord);
253 AliInfo(Form("Found %d constraints records",fTreeConstr->GetEntries()));
254 //
255 }
256 else {
257 //
258 fTreeConstr = new TTree("AliMillePedeRecords_Constraints","Constraints Records for AliMillePede2");
259 fTreeConstr->Branch("Record_Constraints","AliMillePedeRecord",&fRecord,32000,99);
260 }
261 fCurrRecConstrID = -1;
262 //
263 return kTRUE;
264}
265
266//_____________________________________________________________________________________________
267void AliMillePede2::CloseDataRecStorage()
268{
269 if (fTreeData) {
270 if (fDataRecFile->IsWritable()) {
271 fDataRecFile->cd();
272 fTreeData->Write();
273 }
274 delete fTreeData;
275 fTreeData = 0;
276 fDataRecFile->Close();
277 delete fDataRecFile;
278 fDataRecFile = 0;
279 }
280 //
281}
282
283//_____________________________________________________________________________________________
284void AliMillePede2::CloseConsRecStorage()
285{
286 if (fTreeConstr) {
287 if (fConsRecFile->IsWritable()) {
288 fConsRecFile->cd();
289 fTreeConstr->Write();
290 }
291 delete fTreeConstr;
292 fTreeConstr = 0;
293 fConsRecFile->Close();
294 delete fConsRecFile;
295 fConsRecFile = 0;
296 }
297 //
298}
299
300//_____________________________________________________________________________________________
301Bool_t AliMillePede2::ReadNextRecordData()
302{
303 // read next data record (if any)
304 if (!fTreeData || ++fCurrRecDataID >= fTreeData->GetEntries()) { fCurrRecDataID--; return kFALSE;}
305 fTreeData->GetEntry(fCurrRecDataID);
306 return kTRUE;
307}
308
309//_____________________________________________________________________________________________
310Bool_t AliMillePede2::ReadNextRecordConstraint()
311{
312 // read next constraint record (if any)
313 if (!fTreeConstr || ++fCurrRecConstrID >= fTreeConstr->GetEntries()) { fCurrRecConstrID--; return kFALSE;}
314 fTreeConstr->GetEntry(fCurrRecConstrID);
315 return kTRUE;
316}
317
318//_____________________________________________________________________________________________
319void AliMillePede2::SetLocalEquation(double *dergb, double *derlc, double lMeas, double lSigma)
320{
321 if (!fDataRecFile || !fDataRecFile->IsWritable()) InitDataRecStorage(); // create a buffer to store the data
322 //
323 // write data of single measurement
324 if (lSigma<=0.0) { // If parameter is fixed, then no equation
325 for (int i=fNLocPar; i--;) derlc[i] = 0.0;
326 for (int i=fNGloPar; i--;) dergb[i] = 0.0;
327 return;
328 }
329 //
330 fRecord->AddResidual(lMeas);
331 //
332 // Retrieve local param interesting indices
333 for (int i=0;i<fNLocPar;i++) if (derlc[i]!=0.0) {fRecord->AddIndexValue(i,derlc[i]); derlc[i] = 0.0;}
334 //
335 fRecord->AddWeight( 1.0/lSigma/lSigma );
336 //
337 // Idem for global parameters
de34b538 338 for (int i=0;i<fNGloPar;i++) if (dergb[i]!=0.0) {
339 fRecord->AddIndexValue(i,dergb[i]); dergb[i] = 0.0;
340 fRecord->MarkGroup(fParamGrID[i]);
341 }
8a9ab0eb 342 //
343}
344
345//_____________________________________________________________________________________________
346void AliMillePede2::SetLocalEquation(int *indgb, double *dergb, int ngb, int *indlc,
347 double *derlc,int nlc,double lMeas,double lSigma)
348{
349 // write data of single measurement
350 if (lSigma<=0.0) { // If parameter is fixed, then no equation
351 for (int i=nlc;i--;) derlc[i] = 0.0;
352 for (int i=ngb;i--;) dergb[i] = 0.0;
353 return;
354 }
355 //
356 if (!fDataRecFile || !fDataRecFile->IsWritable()) InitDataRecStorage(); // create a buffer to store the data
357 //
358 fRecord->AddResidual(lMeas);
359 //
360 // Retrieve local param interesting indices
361 for (int i=0;i<nlc;i++) if (derlc[i]!=0.0) {fRecord->AddIndexValue(indlc[i],derlc[i]); derlc[i]=0.; indlc[i]=0;}
362 //
363 fRecord->AddWeight( 1./lSigma/lSigma );
364 //
365 // Idem for global parameters
366 for (int i=0;i<ngb;i++) if (dergb[i]!=0.0) {fRecord->AddIndexValue(indgb[i],dergb[i]); dergb[i]=0.; indgb[i]=0;}
367 //
368}
369
de34b538 370
8a9ab0eb 371//_____________________________________________________________________________________________
de34b538 372void AliMillePede2::SetGlobalConstraint(double *dergb, double val, double sigma)
8a9ab0eb 373{
374 // Define a constraint equation.
375 if (!fConsRecFile || !fConsRecFile->IsWritable()) InitConsRecStorage(); // create a buffer to store the data
376 //
377 fRecord->Reset();
378 fRecord->AddResidual(val);
de34b538 379 fRecord->AddWeight(sigma);
8a9ab0eb 380 for (int i=0; i<fNGloPar; i++) if (dergb[i]!=0) fRecord->AddIndexValue(i,dergb[i]);
381 fNGloConstraints++;
de34b538 382 if (sigma==0) fNLagrangeConstraints++;
383 // printf("NewConstraint:\n"); fRecord->Print(); //RRR
8a9ab0eb 384 SaveRecordConstraint();
385 //
386}
387
388//_____________________________________________________________________________________________
de34b538 389void AliMillePede2::SetGlobalConstraint(const int *indgb, double *dergb, int ngb, double val,double sigma)
8a9ab0eb 390{
391 // Define a constraint equation.
392 if (!fConsRecFile || !fConsRecFile->IsWritable()) InitConsRecStorage(); // create a buffer to store the data
393 fRecord->Reset();
394 fRecord->AddResidual(val);
de34b538 395 fRecord->AddWeight(sigma); // dummy
8a9ab0eb 396 for (int i=0; i<ngb; i++) if (dergb[i]!=0) fRecord->AddIndexValue(indgb[i],dergb[i]);
397 fNGloConstraints++;
de34b538 398 if (sigma==0) fNLagrangeConstraints++;
8a9ab0eb 399 SaveRecordConstraint();
400 //
401}
402
403//_____________________________________________________________________________________________
404Int_t AliMillePede2::LocalFit(double *localParams)
405{
406 /*
407 Perform local parameters fit once all the local equations have been set
408 -----------------------------------------------------------
409 localParams = (if !=0) will contain the fitted track parameters and
410 related errors
411 */
de34b538 412 static int nrefSize = 0;
413 // static TArrayI refLoc,refGlo,nrefLoc,nrefGlo;
414 static Int_t *refLoc=0,*refGlo=0,*nrefLoc=0,*nrefGlo=0;
8a9ab0eb 415 int nPoints = 0;
416 //
de34b538 417 AliSymMatrix &matCLoc = *fMatCLoc;
418 AliMatrixSq &matCGlo = *fMatCGlo;
419 AliRectMatrix &matCGloLoc = *fMatCGloLoc;
8a9ab0eb 420 //
421 memset(fVecBLoc,0,fNLocPar*sizeof(double));
de34b538 422 matCLoc.Reset();
8a9ab0eb 423 //
424 int cnt=0;
425 int recSz = fRecord->GetSize();
426 //
427 while(cnt<recSz) { // Transfer the measurement records to matrices
428 //
429 // extract addresses of residual, weight and pointers on local and global derivatives for each point
de34b538 430 if (nrefSize<=nPoints) {
431 int *tmpA = 0;
432 nrefSize = 2*(nPoints+1);
433 tmpA = refLoc; refLoc = new Int_t[nrefSize]; if (tmpA) memcpy(refLoc,tmpA,nPoints*sizeof(int));
434 tmpA = refGlo; refGlo = new Int_t[nrefSize]; if (tmpA) memcpy(refGlo,tmpA,nPoints*sizeof(int));
435 tmpA = nrefLoc; nrefLoc = new Int_t[nrefSize]; if (tmpA) memcpy(nrefLoc,tmpA,nPoints*sizeof(int));
436 tmpA = nrefGlo; nrefGlo = new Int_t[nrefSize]; if (tmpA) memcpy(nrefGlo,tmpA,nPoints*sizeof(int));
8a9ab0eb 437 }
438 //
de34b538 439 refLoc[nPoints] = ++cnt;
8a9ab0eb 440 int nLoc = 0;
441 while(!fRecord->IsWeight(cnt)) {nLoc++; cnt++;}
de34b538 442 nrefLoc[nPoints] = nLoc;
8a9ab0eb 443 //
de34b538 444 refGlo[nPoints] = ++cnt;
8a9ab0eb 445 int nGlo = 0;
446 while(!fRecord->IsResidual(cnt) && cnt<recSz) {nGlo++; cnt++;}
de34b538 447 nrefGlo[nPoints] = nGlo;
8a9ab0eb 448 //
449 nPoints++;
450 }
451 double vl;
452 //
7c3070ec 453 Int_t maxLocUsed = 0;
454 //
8a9ab0eb 455 for (int ip=nPoints;ip--;) { // Transfer the measurement records to matrices
de34b538 456 double resid = fRecord->GetValue( refLoc[ip]-1 );
457 double weight = fRecord->GetValue( refGlo[ip]-1 );
458 double *derLoc = fRecord->GetValue()+refLoc[ip];
459 double *derGlo = fRecord->GetValue()+refGlo[ip];
460 int *indLoc = fRecord->GetIndex()+refLoc[ip];
461 int *indGlo = fRecord->GetIndex()+refGlo[ip];
462 //
463 for (int i=nrefGlo[ip];i--;) { // suppress the global part (only relevant with iterations)
8a9ab0eb 464 int iID = indGlo[i]; // Global param indice
de34b538 465 if (fSigmaPar[iID]<=0.) continue; // fixed parameter RRRCheck
8a9ab0eb 466 if (fIsLinear[iID]) resid -= derGlo[i]*(fInitPar[iID]+fDeltaPar[iID]); // linear parameter
467 else resid -= derGlo[i]*fDeltaPar[iID]; // nonlinear parameter
468 }
469 //
de34b538 470 // Symmetric matrix, don't bother j>i coeffs
471 for (int i=nrefLoc[ip];i--;) { // Fill local matrix and vector
472 fVecBLoc[ indLoc[i] ] += weight*resid*derLoc[i];
7c3070ec 473 if (indLoc[i]>maxLocUsed) maxLocUsed = indLoc[i];
de34b538 474 for (int j=i+1;j--;) matCLoc(indLoc[i] ,indLoc[j]) += weight*derLoc[i]*derLoc[j];
7c3070ec 475 }
8a9ab0eb 476 //
477 } // end of the transfer of the measurement record to matrices
478 //
7c3070ec 479 matCLoc.SetSizeUsed(++maxLocUsed); // data with B=0 may use less than declared nLocals
480 //
8a9ab0eb 481 // first try to solve by faster Cholesky decomposition, then by Gaussian elimination
de34b538 482 if (!matCLoc.SolveChol(fVecBLoc,kTRUE)) {
8a9ab0eb 483 AliInfo("Failed to solve locals by Cholesky, traying Gaussian Elimination");
de34b538 484 if (!matCLoc.SolveSpmInv(fVecBLoc,kTRUE)) {
8a9ab0eb 485 AliInfo("Failed to solve locals by Gaussian Elimination, skip...");
de34b538 486 matCLoc.Print("d");
8a9ab0eb 487 return 0; // failed to solve
488 }
489 }
490 //
491 // If requested, store the track params and errors
de34b538 492 // printf("locfit: "); for (int i=0;i<fNLocPar;i++) printf("%+e |",fVecBLoc[i]); printf("\n");
7c3070ec 493 if (localParams) for (int i=maxLocUsed; i--;) {
8a9ab0eb 494 localParams[2*i] = fVecBLoc[i];
de34b538 495 localParams[2*i+1] = TMath::Sqrt(TMath::Abs(matCLoc.QueryDiag(i)));
8a9ab0eb 496 }
497 //
498 float lChi2 = 0;
499 int nEq = 0;
500 //
501 for (int ip=nPoints;ip--;) { // Calculate residuals
de34b538 502 double resid = fRecord->GetValue( refLoc[ip]-1 );
503 double weight = fRecord->GetValue( refGlo[ip]-1 );
504 double *derLoc = fRecord->GetValue()+refLoc[ip];
505 double *derGlo = fRecord->GetValue()+refGlo[ip];
506 int *indLoc = fRecord->GetIndex()+refLoc[ip];
507 int *indGlo = fRecord->GetIndex()+refGlo[ip];
8a9ab0eb 508 //
509 // Suppress local and global contribution in residuals;
de34b538 510 for (int i=nrefLoc[ip];i--;) resid -= derLoc[i]*fVecBLoc[ indLoc[i] ]; // local part
8a9ab0eb 511 //
de34b538 512 for (int i=nrefGlo[ip];i--;) { // global part
8a9ab0eb 513 int iID = indGlo[i];
514 if ( fSigmaPar[iID] <= 0.) continue; // fixed parameter RRRCheck
515 if (fIsLinear[iID]) resid -= derGlo[i]*(fInitPar[iID]+fDeltaPar[iID]); // linear parameter
516 else resid -= derGlo[i]*fDeltaPar[iID]; // nonlinear parameter
517 }
518 //
519 // reject the track if the residual is too large (outlier)
de34b538 520 double absres = TMath::Abs(resid);
521 if ( (absres >= fResCutInit && fIter ==1 ) ||
522 (absres >= fResCut && fIter > 1)) {
523 if (fLocFitAdd) fNLocFitsRejected++;
8a9ab0eb 524 return 0;
525 }
526 //
527 lChi2 += weight*resid*resid ; // total chi^2
528 nEq++; // number of equations
529 } // end of Calculate residuals
530 //
7c3070ec 531 int nDoF = nEq-maxLocUsed;
8a9ab0eb 532 lChi2 = (nDoF>0) ? lChi2/nDoF : 0; // Chi^2/dof
533 //
534 if (fNStdDev != 0 && nDoF>0 && lChi2 > Chi2DoFLim(fNStdDev,nDoF)*fChi2CutFactor) { // check final chi2
de34b538 535 if (fLocFitAdd) fNLocFitsRejected++;
8a9ab0eb 536 return 0;
537 }
538 //
de34b538 539 if (fLocFitAdd) {
540 fNLocFits++;
541 fNLocEquations += nEq;
542 }
543 else {
544 fNLocFits--;
545 fNLocEquations -= nEq;
546 }
8a9ab0eb 547 //
548 // local operations are finished, track is accepted
549 // We now update the global parameters (other matrices)
550 //
551 int nGloInFit = 0;
552 //
553 for (int ip=nPoints;ip--;) { // Update matrices
de34b538 554 double resid = fRecord->GetValue( refLoc[ip]-1 );
555 double weight = fRecord->GetValue( refGlo[ip]-1 );
556 double *derLoc = fRecord->GetValue()+refLoc[ip];
557 double *derGlo = fRecord->GetValue()+refGlo[ip];
558 int *indLoc = fRecord->GetIndex()+refLoc[ip];
559 int *indGlo = fRecord->GetIndex()+refGlo[ip];
560 //
561 for (int i=nrefGlo[ip];i--;) { // suppress the global part
8a9ab0eb 562 int iID = indGlo[i]; // Global param indice
563 if ( fSigmaPar[iID] <= 0.) continue; // fixed parameter RRRCheck
de34b538 564 if (fIsLinear[iID]) resid -= derGlo[i]*(fInitPar[iID]+fDeltaPar[iID]); // linear parameter
565 else resid -= derGlo[i]*fDeltaPar[iID]; // nonlinear parameter
8a9ab0eb 566 }
567 //
de34b538 568 for (int ig=nrefGlo[ip];ig--;) {
8a9ab0eb 569 int iIDg = indGlo[ig]; // Global param indice (the matrix line)
570 if ( fSigmaPar[iIDg] <= 0.) continue; // fixed parameter RRRCheck
de34b538 571 if (fLocFitAdd) fVecBGlo[ iIDg ] += weight*resid*derGlo[ig]; //!!!
572 else fVecBGlo[ iIDg ] -= weight*resid*derGlo[ig]; //!!!
8a9ab0eb 573 //
574 // First of all, the global/global terms (exactly like local matrix)
de34b538 575 int nfill = 0;
576 for (int jg=ig+1;jg--;) { // matCGlo is symmetric by construction
8a9ab0eb 577 int jIDg = indGlo[jg];
578 if ( fSigmaPar[jIDg] <= 0.) continue; // fixed parameter RRRCheck
de34b538 579 if ( (vl = weight*derGlo[ig]*derGlo[jg])!=0 ) {
580 fFillIndex[nfill] = jIDg;
581 fFillValue[nfill++] = fLocFitAdd ? vl:-vl;
582 }
583 }
584 if (nfill) matCGlo.AddToRow(iIDg,fFillValue,fFillIndex,nfill);
8a9ab0eb 585 //
586 // Now we have also rectangular matrices containing global/local terms.
587 int iCIDg = fGlo2CGlo[iIDg]; // compressed Index of index
588 if (iCIDg == -1) {
de34b538 589 Double_t *rowGL = matCGloLoc(nGloInFit);
7c3070ec 590 for (int k=maxLocUsed;k--;) rowGL[k] = 0.0; // reset the row
8a9ab0eb 591 iCIDg = fGlo2CGlo[iIDg] = nGloInFit;
592 fCGlo2Glo[nGloInFit++] = iIDg;
593 }
594 //
de34b538 595 Double_t *rowGLIDg = matCGloLoc(iCIDg);
596 for (int il=nrefLoc[ip];il--;) rowGLIDg[ indLoc[il] ] += weight*derGlo[ig]*derLoc[il];
597 fProcPnt[iIDg] += fLocFitAdd ? 1:-1; // update counter
8a9ab0eb 598 //
599 }
600 } // end of Update matrices
601 //
602 // calculate fMatCGlo -= fMatCGloLoc * fMatCLoc * fMatCGloLoc^T
603 // and fVecBGlo -= fMatCGloLoc * fVecBLoc
604 //
de34b538 605 //-------------------------------------------------------------- >>>
8a9ab0eb 606 double vll;
607 for (int iCIDg=0; iCIDg<nGloInFit; iCIDg++) {
608 int iIDg = fCGlo2Glo[iCIDg];
609 //
610 vl = 0;
de34b538 611 Double_t *rowGLIDg = matCGloLoc(iCIDg);
7c3070ec 612 for (int kl=0;kl<maxLocUsed;kl++) if (rowGLIDg[kl]) vl += rowGLIDg[kl]*fVecBLoc[kl];
de34b538 613 if (vl!=0) fVecBGlo[iIDg] -= fLocFitAdd ? vl : -vl;
8a9ab0eb 614 //
de34b538 615 int nfill = 0;
8a9ab0eb 616 for (int jCIDg=0;jCIDg<=iCIDg; jCIDg++) {
617 int jIDg = fCGlo2Glo[jCIDg];
618 //
619 vl = 0;
de34b538 620 Double_t *rowGLJDg = matCGloLoc(jCIDg);
7c3070ec 621 for (int kl=0;kl<maxLocUsed;kl++) {
8a9ab0eb 622 // diag terms
de34b538 623 if ( (vll=rowGLIDg[kl]*rowGLJDg[kl])!=0 ) vl += matCLoc.QueryDiag(kl)*vll;
8a9ab0eb 624 //
625 // off-diag terms
626 for (int ll=0;ll<kl;ll++) {
de34b538 627 if ( (vll=rowGLIDg[kl]*rowGLJDg[ll])!=0 ) vl += matCLoc(kl,ll)*vll;
628 if ( (vll=rowGLIDg[ll]*rowGLJDg[kl])!=0 ) vl += matCLoc(kl,ll)*vll;
8a9ab0eb 629 }
630 }
de34b538 631 if (vl!=0) {
632 fFillIndex[nfill] = jIDg;
633 fFillValue[nfill++] = fLocFitAdd ? -vl : vl;
634 }
8a9ab0eb 635 }
de34b538 636 if (nfill) matCGlo.AddToRow(iIDg,fFillValue,fFillIndex,nfill);
8a9ab0eb 637 }
638 //
639 // reset compressed index array
640 //
de34b538 641 for (int i=nGloInFit;i--;) {
642 fGlo2CGlo[ fCGlo2Glo[i] ] = -1;
643 fCGlo2Glo[i] = -1;
644 }
8a9ab0eb 645 //
de34b538 646 //---------------------------------------------------- <<<
8a9ab0eb 647 return 1;
648}
649
650//_____________________________________________________________________________________________
651Int_t AliMillePede2::GlobalFit(Double_t *par, Double_t *error, Double_t *pull)
652{
653 // performs a requested number of global iterations
654 fIter = 1;
655 //
656 TStopwatch sw; sw.Start();
657 //
658 int res = 0;
659 AliInfo("Starting Global fit.");
660 while (fIter<=fMaxIter) {
661 //
662 res = GlobalFitIteration();
663 if (!res) break;
664 //
665 if (fChi2CutFactor != fChi2CutRef) {
666 fChi2CutFactor = TMath::Sqrt(fChi2CutFactor);
667 if (fChi2CutFactor < 1.2*fChi2CutRef) {
668 fChi2CutFactor = fChi2CutRef;
669 fIter = fMaxIter - 1; // Last iteration
670 }
671 }
672 fIter++;
673 }
674 //
675 sw.Stop();
676 AliInfo(Form("Global fit %s, CPU time: %.1f",res ? "Converged":"Failed",sw.CpuTime()));
677 if (!res) return 0;
678 //
679 if (par) for (int i=fNGloPar;i--;) par[i] = fInitPar[i]+fDeltaPar[i];
680 //
681 if (fGloSolveStatus==gkInvert) { // errors on params are available
de34b538 682 if (error) for (int i=fNGloPar;i--;) error[i] = fProcPnt[i]>0 ? TMath::Sqrt(TMath::Abs(fMatCGlo->QueryDiag(i))) : 0.;
683 if (pull) for (int i=fNGloPar;i--;) pull[i] = fProcPnt[i]>0 && (fSigmaPar[i]*fSigmaPar[i]-fMatCGlo->QueryDiag(i))>0. && fSigmaPar[i]>0
684 ? fDeltaPar[i]/TMath::Sqrt(fSigmaPar[i]*fSigmaPar[i]-fMatCGlo->QueryDiag(i)) : 0;
8a9ab0eb 685 }
686 //
687 return 1;
688}
689
690//_____________________________________________________________________________________________
691Int_t AliMillePede2::GlobalFitIteration()
692{
693 // perform global parameters fit once all the local equations have been fitted
694 //
695 AliInfo(Form("Global Fit Iteration#%2d (Local Fit Chi^2 cut factor: %.2f)",fIter,fChi2CutFactor));
696 //
697 if (!fNGloPar || !fTreeData) {
698 AliInfo("No data was stored, aborting iteration");
699 return 0;
700 }
de34b538 701 TStopwatch sw,sws;
702 sw.Start();
703 sws.Stop();
8a9ab0eb 704 //
705 if (!fConstrUsed) {
706 fConstrUsed = new Bool_t[fNGloConstraints];
707 memset(fConstrUsed,0,fNGloConstraints*sizeof(Bool_t));
708 }
709 // Reset all info specific for this step
de34b538 710 AliMatrixSq& matCGlo = *fMatCGlo;
711 matCGlo.Reset();
8a9ab0eb 712 memset(fProcPnt,0,fNGloPar*sizeof(Int_t));
713 //
714 fNGloConstraints = fTreeConstr ? fTreeConstr->GetEntries() : 0;
715 //
de34b538 716 // count number of Lagrange constraints: they need new row/cols to be added
717 fNLagrangeConstraints = 0;
718 for (int i=0; i<fNGloConstraints; i++) {
719 ReadRecordConstraint(i);
720 if ( fRecord->GetValue(1)==0 ) fNLagrangeConstraints++; // exact constraint (no error) -> Lagrange multiplier
721 }
722 //
8a9ab0eb 723 // if needed, readjust the size of the global vector (for matrices this is done automatically)
de34b538 724 if (!fVecBGlo || fNGloSize!=fNGloPar+fNLagrangeConstraints) {
8a9ab0eb 725 delete[] fVecBGlo; // in case some constraint was added between the two manual iterations
de34b538 726 fNGloSize = fNGloPar+fNLagrangeConstraints;
8a9ab0eb 727 fVecBGlo = new Double_t[fNGloSize];
728 }
729 memset(fVecBGlo,0,fNGloSize*sizeof(double));
730 //
731 fNLocFits = 0;
732 fNLocFitsRejected = 0;
de34b538 733 fNLocEquations = 0;
8a9ab0eb 734 //
735 // Process data records and build the matrices
736 Long_t ndr = fTreeData->GetEntries();
737 AliInfo(Form("Building the Global matrix from %d data records",ndr));
de34b538 738 if (!ndr) return 0;
739 //
740 TStopwatch swt; swt.Start();
741 fLocFitAdd = kTRUE; // add contributions of matching tracks
8a9ab0eb 742 for (Long_t i=0;i<ndr;i++) {
743 ReadRecordData(i);
744 LocalFit();
7c3070ec 745 if ( (i%int(0.2*ndr)) == 0) printf("%.1f%% of local fits done\n", double(100.*i)/ndr);
8a9ab0eb 746 }
de34b538 747 swt.Stop();
748 printf("%ld local fits done: ", ndr);
749 swt.Print();
750 sw.Start(kFALSE);
751 //
8a9ab0eb 752 //
de34b538 753 // ---------------------- Reject parameters with low statistics ------------>>
8a9ab0eb 754 fNGloFix = 0;
de34b538 755 if (fMinPntValid>1 && fNGroupsSet) {
756 //
757 printf("Checking parameters with statistics < %d\n",fMinPntValid);
758 TStopwatch swsup;
759 swsup.Start();
760 // 1) build the list of parameters to fix
761 Int_t fixArrSize = 10;
762 Int_t nFixedGroups = 0;
763 TArrayI fixGroups(fixArrSize);
764 //
765 for (int i=fNGloPar;i--;) { // // Reset row and column of fixed params and add 1/sig^2 to free ones
766 int grID = fParamGrID[i];
767 if (grID<0) continue; // not in the group
768 if (fProcPnt[i]>=fMinPntValid) continue;
769 fProcPnt[i] = 0;
770 // check if the group is already accounted
771 Bool_t fnd = kFALSE;
772 for (int j=nFixedGroups;j--;) if (fixGroups[j]==grID) {fnd=kTRUE; break;}
773 if (fnd) continue;
774 if (nFixedGroups>=fixArrSize) {fixArrSize*=2; fixGroups.Set(fixArrSize);}
775 fixGroups[nFixedGroups++] = grID; // add group to fix
776 //
777 }
778 // 2) loop over records and add contributions of fixed groups with negative sign
779 fLocFitAdd = kFALSE;
780 //
781 for (Long_t i=0;i<ndr;i++) {
782 ReadRecordData(i);
783 Bool_t suppr = kFALSE;
784 for (int ifx=nFixedGroups;ifx--;)if (fRecord->IsGroupPresent(fixGroups[ifx])) suppr = kTRUE;
785 if (suppr) LocalFit();
786 }
787 fLocFitAdd = kTRUE;
788 //
789 if (nFixedGroups) {
790 printf("Suppressed contributions of groups with NPoints<%d :\n",fMinPntValid);
791 for (int i=0;i<nFixedGroups;i++) printf("%d ",fixGroups[i]); printf("\n");
792 }
793 swsup.Stop();
794 swsup.Print();
795 }
796 // ---------------------- Reject parameters with low statistics ------------<<
797 //
798 // add large number to diagonal of fixed params
8a9ab0eb 799 //
8a9ab0eb 800 for (int i=fNGloPar;i--;) { // // Reset row and column of fixed params and add 1/sig^2 to free ones
7c3070ec 801 // printf("#%3d : Nproc : %5d grp: %d\n",i,fProcPnt[i],fParamGrID[i]);
de34b538 802 if (fProcPnt[i]<1) {
8a9ab0eb 803 fNGloFix++;
8a9ab0eb 804 fVecBGlo[i] = 0.;
de34b538 805 matCGlo.DiagElem(i) = float(fNLocEquations*fNLocEquations);
8a9ab0eb 806 }
de34b538 807 else matCGlo.DiagElem(i) += (fgWeightSigma ? fProcPnt[i] : 1.)/(fSigmaPar[i]*fSigmaPar[i]);
8a9ab0eb 808 }
809 //
de34b538 810 for (int i=fNGloPar;i--;) fDiagCGlo[i] = matCGlo.QueryDiag(i); // save the diagonal elements
811 //
8a9ab0eb 812 // add constraint equations
813 int nVar = fNGloPar; // Current size of global matrix
814 for (int i=0; i<fNGloConstraints; i++) {
815 ReadRecordConstraint(i);
816 double val = fRecord->GetValue(0);
de34b538 817 double sig = fRecord->GetValue(1);
8a9ab0eb 818 int *indV = fRecord->GetIndex()+2;
819 double *der = fRecord->GetValue()+2;
820 int csize = fRecord->GetSize()-2;
821 //
de34b538 822 // check if after suppression of fixed variables there are non-0 derivatives
823 // and determine the max statistics of involved params
824 int nSuppressed = 0;
825 int maxStat = 1;
826 for (int j=csize;j--;) {
827 if (fProcPnt[indV[j]]<1) nSuppressed++;
828 else {
829 maxStat = TMath::Max(maxStat,fProcPnt[indV[j]]);
830 }
831 }
8a9ab0eb 832 //
de34b538 833 if (nSuppressed==csize) {
834 // AliInfo(Form("Neglecting constraint %d of %d derivatives since no free parameters left",i,csize));
8a9ab0eb 835 //
836 // was this constraint ever created ?
de34b538 837 if ( sig==0 && fConstrUsed[i] ) { // this is needed only for constraints with Lagrange multiplier
8a9ab0eb 838 // to avoid empty row impose dummy constraint on "Lagrange multiplier"
de34b538 839 matCGlo.DiagElem(nVar) = 1.;
8a9ab0eb 840 fVecBGlo[nVar++] = 0;
841 }
842 continue;
843 }
844 //
de34b538 845 // account for already accumulated corrections
846 for (int j=csize; j--;) val -= der[j]*(fInitPar[ indV[j] ]+fDeltaPar[ indV[j] ]);
8a9ab0eb 847 //
de34b538 848 if (sig > 0) { // this is a gaussian constriant: no Lagrange multipliers are added
849 //
850 double sig2i = (fgWeightSigma ? TMath::Sqrt(maxStat) : 1.)/sig/sig;
851 for (int ir=0;ir<csize;ir++) {
852 int iID = indV[ir];
853 for (int ic=0;ic<=ir;ic++) { // matrix is symmetric
854 int jID = indV[ic];
855 double vl = der[ir]*der[ic]*sig2i;
856 if (vl!=0) matCGlo(iID,jID) += vl;
857 }
858 fVecBGlo[iID] += val*der[ir]*sig2i;
859 }
860 }
861 else { // this is exact constriant: Lagrange multipliers must be added
862 for (int j=csize; j--;) {
863 int jID = indV[j];
864 if (fProcPnt[jID]<1) continue; // this parameter was fixed, don't put it into constraint
865 matCGlo(nVar,jID) = float(fNLocEquations)*der[j]; // fMatCGlo is symmetric, only lower triangle is filled
866 }
867 //
868 if (matCGlo.QueryDiag(nVar)) matCGlo.DiagElem(nVar) = 0.0;
869 fVecBGlo[nVar++] = float(fNLocEquations)*val; //RS ? should we use here fNLocFits ?
870 fConstrUsed[i] = kTRUE;
871 }
8a9ab0eb 872 }
873 //
874 AliInfo(Form("Obtained %-7ld equations from %-7ld records (%-7ld rejected). Fixed %-4d globals",
875 fNLocEquations,fNLocFits,fNLocFitsRejected,fNGloFix));
876
877 //
de34b538 878 sws.Start();
8a9ab0eb 879 fGloSolveStatus = SolveGlobalMatEq(); // obtain solution for this step
de34b538 880 sws.Stop();
881 printf("Solve %d |",fIter); sws.Print();
8a9ab0eb 882 //
883 sw.Stop();
884 AliInfo(Form("Iteration#%2d %s. CPU time: %.1f",fIter,fGloSolveStatus==gkFailed ? "Failed":"Converged",sw.CpuTime()));
885 if (fGloSolveStatus==gkFailed) return 0;
886 //
887 for (int i=fNGloPar;i--;) fDeltaPar[i] += fVecBGlo[i]; // Update global parameters values (for iterations)
888 //
de34b538 889 // PrintGlobalParameters();
8a9ab0eb 890 return 1;
891}
892
893//_____________________________________________________________________________________________
894Int_t AliMillePede2::SolveGlobalMatEq()
895{
896 //
897 // solve global matrix equation MatCGlob*X=VecBGlo and store the result in the VecBGlo
898 //
de34b538 899 /*
900 printf("GlobalMatrix\n");
901 fMatCGlo->Print();
902 printf("RHS\n");
903 for (int i=0;i<fNGloPar;i++) printf("%d %+e\n",i,fVecBGlo[i]);
904 */
905 //
8a9ab0eb 906 if (!fgIsMatGloSparse) {
907 //
de34b538 908 if (fNLagrangeConstraints==0) { // pos-def systems are faster to solve by Cholesky
909 if ( ((AliSymMatrix*)fMatCGlo)->SolveChol(fVecBGlo, fgInvChol) ) return fgInvChol ? gkInvert:gkNoInversion;
8a9ab0eb 910 else AliInfo("Solution of Global Dense System by Cholesky failed, trying Gaussian Elimiation");
911 }
912 //
913 if (((AliSymMatrix*)fMatCGlo)->SolveSpmInv(fVecBGlo, kTRUE)) return gkInvert;
914 else AliInfo("Solution of Global Dense System by Gaussian Elimiation failed, trying iterative methods");
915 }
916 // try to solve by minres
de34b538 917 TVectorD sol(fNGloSize);
8a9ab0eb 918 //
919 AliMinResSolve *slv = new AliMinResSolve(fMatCGlo,fVecBGlo);
920 if (!slv) return gkFailed;
921 //
922 Bool_t res = kFALSE;
923 if (fgIterSol == AliMinResSolve::kSolMinRes)
de34b538 924 res = slv->SolveMinRes(sol,fgMinResCondType,fgMinResMaxIter,fgMinResTol);
8a9ab0eb 925 else if (fgIterSol == AliMinResSolve::kSolFGMRes)
de34b538 926 res = slv->SolveFGMRES(sol,fgMinResCondType,fgMinResMaxIter,fgMinResTol,fgNKrylovV);
8a9ab0eb 927 else
928 AliInfo(Form("Undefined Iteritive Solver ID=%d, only %d are defined",fgIterSol,AliMinResSolve::kNSolvers));
929 //
930 if (!res) return gkFailed;
de34b538 931 for (int i=fNGloSize;i--;) fVecBGlo[i] = sol[i];
932 return gkNoInversion;
8a9ab0eb 933 //
934}
935
936//_____________________________________________________________________________________________
937Float_t AliMillePede2::Chi2DoFLim(int nSig, int nDoF) const
938{
939 /// return the limit in chi^2/nd for n sigmas stdev authorized
940 // Only n=1, 2, and 3 are expected in input
941 int lNSig;
942 float sn[3] = {0.47523, 1.690140, 2.782170};
943 float table[3][30] = {{1.0000, 1.1479, 1.1753, 1.1798, 1.1775, 1.1730, 1.1680, 1.1630,
944 1.1581, 1.1536, 1.1493, 1.1454, 1.1417, 1.1383, 1.1351, 1.1321,
945 1.1293, 1.1266, 1.1242, 1.1218, 1.1196, 1.1175, 1.1155, 1.1136,
946 1.1119, 1.1101, 1.1085, 1.1070, 1.1055, 1.1040},
947 {4.0000, 3.0900, 2.6750, 2.4290, 2.2628, 2.1415, 2.0481, 1.9736,
948 1.9124, 1.8610, 1.8171, 1.7791, 1.7457, 1.7161, 1.6897, 1.6658,
949 1.6442, 1.6246, 1.6065, 1.5899, 1.5745, 1.5603, 1.5470, 1.5346,
950 1.5230, 1.5120, 1.5017, 1.4920, 1.4829, 1.4742},
951 {9.0000, 5.9146, 4.7184, 4.0628, 3.6410, 3.3436, 3.1209, 2.9468,
952 2.8063, 2.6902, 2.5922, 2.5082, 2.4352, 2.3711, 2.3143, 2.2635,
953 2.2178, 2.1764, 2.1386, 2.1040, 2.0722, 2.0428, 2.0155, 1.9901,
954 1.9665, 1.9443, 1.9235, 1.9040, 1.8855, 1.8681}};
955
956 if (nDoF < 1) {
957 return 0.0;
958 }
959 else {
960 lNSig = TMath::Max(1,TMath::Min(nSig,3));
961
962 if (nDoF <= 30) {
963 return table[lNSig-1][nDoF-1];
964 }
965 else { // approximation
966 return ((sn[lNSig-1]+TMath::Sqrt(float(2*nDoF-3)))*
967 (sn[lNSig-1]+TMath::Sqrt(float(2*nDoF-3))))/float(2*nDoF-2);
968 }
969 }
970}
971
972//_____________________________________________________________________________________________
973Int_t AliMillePede2::SetIterations(double lChi2CutFac)
974{
975 // Number of iterations is calculated from lChi2CutFac
976 fChi2CutFactor = TMath::Max(1.0, lChi2CutFac);
977 //
978 AliInfo(Form("Initial cut factor is %f",fChi2CutFactor));
979 fIter = 1; // Initializes the iteration process
980 return 1;
981}
982
983//_____________________________________________________________________________________________
984Double_t AliMillePede2::GetParError(int iPar) const
985{
986 // return error for parameter iPar
987 if (fGloSolveStatus==gkInvert) {
de34b538 988 double res = fMatCGlo->QueryDiag(iPar);
8a9ab0eb 989 if (res>=0) return TMath::Sqrt(res);
990 }
991 return -1.;
992}
993
994
995//_____________________________________________________________________________________________
996Int_t AliMillePede2::PrintGlobalParameters() const
997{
998 /// Print the final results into the logfile
999 double lError = 0.;
1000 double lGlobalCor =0.;
1001
1002 AliInfo("");
1003 AliInfo(" Result of fit for global parameters");
1004 AliInfo(" ===================================");
1005 AliInfo(" I initial final differ lastcor error gcor Npnt");
1006 AliInfo("----------------------------------------------------------------------------------------------");
1007 //
1008 for (int i=0; i<fNGloPar; i++) {
1009 lError = GetParError(i);
1010 lGlobalCor = 0.0;
1011 //
1012 double dg;
de34b538 1013 if (fGloSolveStatus==gkInvert && TMath::Abs( (dg=fMatCGlo->QueryDiag(i)) *fDiagCGlo[i]) > 0) {
8a9ab0eb 1014 lGlobalCor = TMath::Sqrt(TMath::Abs(1.0-1.0/(dg*fDiagCGlo[i])));
1015 AliInfo(Form("%d\t %.6f\t %.6f\t %.6f\t %.6f\t %.6f\t %.6f\t %6d",
1016 i,fInitPar[i],fInitPar[i]+fDeltaPar[i],fDeltaPar[i],fVecBGlo[i],lError,lGlobalCor,fProcPnt[i]));
1017 }
1018 else {
1019 AliInfo(Form("%d\t %.6f\t %.6f\t %.6f\t %.6f\t OFF\t OFF\t %6d",i,fInitPar[i],fInitPar[i]+fDeltaPar[i],
1020 fDeltaPar[i],fVecBGlo[i],fProcPnt[i]));
1021 }
1022 }
1023 return 1;
1024}