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