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