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