forked from alisw/AliRoot
-
Notifications
You must be signed in to change notification settings - Fork 0
/
AliNDLocalRegression.cxx
1166 lines (1070 loc) · 48.2 KB
/
AliNDLocalRegression.cxx
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
/**************************************************************************
* Copyright(c) 2006-07, ALICE Experiment at CERN, All rights reserved. *
* *
* Author: The ALICE Off-line Project. *
* Contributors are mentioned in the code where appropriate. *
* *
* Permission to use, copy, modify and distribute this software and its *
* documentation strictly for non-commercial purposes is hereby granted *
* without fee, provided that the above copyright notice appears in all *
* copies and that both the copyright notice and this permission notice *
* appear in the supporting documentation. The authors make no claims *
* about the suitability of this software for any purpose. It is *
* provided "as is" without express or implied warranty. *
**************************************************************************/
//-------------------------------------------------------------------------
// Implementation of the AliNDLocalRegression class
//-------------------------------------------------------------------------
/*
Related task: https://alice.its.cern.ch/jira/browse/ATO-193
Algorithm secription - see:
Kernel_smoother: Local polynomial regression
http://en.wikipedia.org/w/index.php?title=Kernel_smoother&oldid=627785784
Formally, the local polynomial regression is computed by solving a weighted least square problem.
Weights are provided as a width of the gausian kernel.
Local fit parameters are computed on the grid defined by axis set defiend by THn.
For example use please check UnitTest:
.L $ALICE_ROOT/../src/STAT/test/AliNDLocalRegressionTest.C+
//
Init:
AliNDLocalRegression *pfitNDIdeal=0;
pfitNDIdeal->SetHistogram((THn*)(hN->Clone()));
pfitNDIdeal->SetCuts(3,0.8,1); // outlier rejection setting see /AliNDLocalRegressionTest.C:UnitTestGaussNoisePlusOutliers() for motivation
pfitNDIdeal->MakeFit(treeIn, "val:err", "xyz0:xyz1","Entry$%2==1", "0.05:0.05","2:2",0.001);
Usage:
Double_t xyz[2]={1,2}
pfitNDIdeal->Eval(xyz);
In TFormulas:
pfitNDGaus0->AddVisualCorrection(pfitNDGaus0,2);
pfitNDGaus1->AddVisualCorrection(pfitNDGaus0,3);
treeIn->Draw("(AliNDLocalRegression::GetCorrND(3,xyz0,xyz1)-AliNDLocalRegression::GetCorrND(2,xyz0,xyz1))/sqrt(AliNDLocalRegression::GetCorrNDError(3,xyz0,xyz1)**2+AliNDLocalRegression::GetCorrNDError(2,xyz0,xyz1)**2)>>pullsGaus01(200,-20,20)","","");
To do:
1.) Statistical error of the local interpolation ignores Gaussian kernel weights
errors are overestimated - find a proper mathematical formula to estimate statistical error of estimator
2.) Implent regularization for smoothing - requesting approximate smoothnes in values and derivative
author: [email protected]
*/
#include <TError.h>
#include <TVectorD.h>
#include "AliNDLocalRegression.h"
#include "THn.h"
#include "TObjString.h"
#include "TTreeStream.h"
#include "AliMathBase.h"
#include "TMatrixD.h"
#include "TRobustEstimator.h"
#include "AliMathBase.h"
#include "TStopwatch.h"
ClassImp(AliNDLocalRegression)
TObjArray *AliNDLocalRegression::fgVisualCorrection=0;
// instance of correction for visualization
Int_t AliNDLocalRegression::fgVerboseLevel=1000;
AliNDLocalRegression::AliNDLocalRegression():
TNamed(),
fHistPoints(0), // ND histogram defining regression granularity
fRobustFractionLTS(0), // fraction of data used for the robust mean and robust rms estimator (LTS https://en.wikipedia.org/wiki/Least_trimmed_squares)
fRobustRMSLTSCut(0), // cut on the robust RMS |value-localmean|<fRobustRMSLTSCut*localRMS
fCutType(0), // type of the cut 0- no cut 1-cut localmean=median, 2-cut localmen=rosbut mean
fInputTree(0), // input tree - object is not owner
fStreamer(0), // optional streamer
fFormulaVal(0), // value:err definition formula
fSelection(0), // point selector formula
fFormulaVar(0), //: separated variable definition formula
fKernelWidthFormula(0), //: separated - kernel width for the regression
fPolDimensionFormula(0), //: separated - polynom for the regression
fNParameters(0), // number of local paramters to fit
fLocalFitParam(0), // local fit parameters
fLocalFitQuality(0), // local fit quality
fLocalFitCovar(0), // local fit covariance matrix
fLocalRobustStat(0), // local robust statistic
fBinIndex(0), //[fNParameters] working arrays current bin index
fBinCenter(0), //[fNParameters] working current local variables - bin center
fBinDelta(0), //[fNParameters] working current local variables - bin delta
fBinWidth(0), //[fNParameters] working current local variables - bin delta
fUseBinNorm(kFALSE) // switch make polynom in units of bins (kTRUE) or in natural units (kFALSE)
{
if (!fgVisualCorrection) fgVisualCorrection= new TObjArray;
}
AliNDLocalRegression::AliNDLocalRegression(const char* name, const char* title):
TNamed(name,title),
fHistPoints(0), // ND histogram defining regression granularity
fRobustFractionLTS(0), // fraction of data used for the robust mean and robust rms estimator (LTS https://en.wikipedia.org/wiki/Least_trimmed_squares)
fRobustRMSLTSCut(0), // cut on the robust RMS |value-localmean|<fRobustRMSLTSCut*localRMS
fCutType(0), // type of the cut 0- no cut 1-cut localmean=median, 2-cut localmen=rosbut mean
fInputTree(0), // input tree - object is not owner
fStreamer(0), // optional streamer
fFormulaVal(0), // value:err definition formula
fSelection(0), // point selector formula
fFormulaVar(0), //: separated variable definition formula
fKernelWidthFormula(0), //: separated - kernel width for the regression
fPolDimensionFormula(0), //: separated - polynom for the regression
fNParameters(0), // number of local paramters to fit
fLocalFitParam(0), // local fit parameters
fLocalFitQuality(0), // local fit quality
fLocalFitCovar(0), // local fit covariance matrix
fLocalRobustStat(0), // local robust statistic
fBinIndex(0), //[fNParameters] working arrays current bin index
fBinCenter(0), //[fNParameters] working current local variables - bin center
fBinDelta(0), //[fNParameters] working current local variables - bin delta
fBinWidth(0), //[fNParameters] working current local variables - bin delta
fUseBinNorm(kFALSE) // switch make polynom in units of bins (kTRUE) or in natural units (kFALSE)
{
}
AliNDLocalRegression::~AliNDLocalRegression(){
//
// destructor
//
if (fHistPoints) delete fHistPoints;
if (fStreamer) delete fStreamer;
// fInputTree(0), //! input tree - object is not owner
delete fFormulaVal; // value:err definition formula
delete fSelection; // point selector formula
delete fFormulaVar; //: separated variable definition formula
delete fKernelWidthFormula; //: separated - kernel width for the regression
delete fPolDimensionFormula; //: separated - polynom for the regression
//
if (fLocalFitParam) fLocalFitParam->Delete();
if (fLocalFitQuality) fLocalFitQuality->Delete();
if (fLocalFitCovar) fLocalFitCovar->Delete();
delete fLocalFitParam; // local fit parameters
delete fLocalFitQuality; // local fit quality
delete fLocalFitCovar; // local fit covariance matrix
//
delete fLocalRobustStat;
//
delete[] fBinIndex;
delete[] fBinCenter;
delete[] fBinDelta;
}
void AliNDLocalRegression::SetHistogram(THn* histo ){
//
// Setup the local regression ayout according THn hitogram binning
//
if (fHistPoints!=0){
::Error("AliNDLocalRegression::SetHistogram","Histogram initialized\n");
return ;
}
fHistPoints=histo;
fLocalFitParam = new TObjArray(fHistPoints->GetNbins());
fLocalFitParam->SetOwner(kTRUE);
fLocalFitQuality = new TObjArray(fHistPoints->GetNbins());
fLocalFitQuality->SetOwner(kTRUE);
fLocalFitCovar = new TObjArray(fHistPoints->GetNbins());
fLocalFitCovar->SetOwner(kTRUE);
//
// Check histogram
//
Int_t ndim = histo->GetNdimensions();
Bool_t isOK=kTRUE;
for (Int_t idim=0; idim<ndim; idim++){
TAxis * axis = histo->GetAxis(idim);
if (axis->GetNbins()<2) {
::Error("AliNDLocalRegression::SetHistogram",TString::Format("Invalid binning nbins<2 %d", axis->GetNbins()).Data());
}
if (axis->GetXmin()>=axis->GetXmax()) {
::Error("AliNDLocalRegression::SetHistogram",TString::Format("Invalid range <%f,%f", axis->GetXmin(),axis->GetXmax()).Data());
}
}
}
void AliNDLocalRegression::SetCuts(Double_t nSigma, Double_t robustFraction, Int_t estimator){
//
//
//
fRobustFractionLTS=robustFraction; // fraction of data used for the robust mean and robust rms estimator (LTS https://en.wikipedia.org/wiki/Least_trimmed_squares)
fRobustRMSLTSCut=nSigma; // cut on the robust RMS |value-localmean|<fRobustRMSLTSCut*localRMS
fCutType=estimator; // type of the cut 0- no cut 1-cut localmean=median, 2-cut localmen=rosbut mean
}
Bool_t AliNDLocalRegression::CleanCovariance(){
//
// Clean covariance matrix if not needed anymore
//
if (fLocalFitCovar) delete fLocalFitCovar;
fLocalFitCovar=0;
};
Bool_t AliNDLocalRegression::MakeFit(TTree * tree , const char* formulaVal, const char * formulaVar, const char*selection, const char * formulaKernel, const char * dimensionFormula, Double_t weightCut, Int_t entries, Bool_t useBinNorm){
//
// Make a local fit in grid as specified by the input THn histogram
// Histogram has to be set before invocation of method
//
// Output:
// array of fit parameters and covariance matrices
//
// Input Parameters:
// tree - input tree
// formulaVal - : separated variable:error string
// formulaVar - : separate varaible list
// selection - selection (cut) for TTreeDraw
// kernelWidth - : separated list of width of kernel for local fitting
// dimenstionFormula - dummy for the moment
//
//Algorithm:
// 1.) Check consistency of input data
//
// 2.) Cache input data from tree to the array of vector TVectorD
//
// 3.) Calculate robust local mean and robust local RMS in case outlier removal algorithm specified
//
// 4.) Make local fit
//
// const Double_t kEpsilon=1e-6;
const Int_t kMaxDim=100;
Int_t binRange[kMaxDim]={0};
Double_t nchi2Cut=-2*TMath::Log(weightCut); // transform probability to nsigma cut
if (fHistPoints==NULL){
::Error("AliNDLocalRegression::MakeFit", "ND histogram not initialized\n");
return kFALSE;
}
if (tree==NULL || tree->GetEntries()==0){
::Error("AliNDLocalRegression::MakeFit", "Empty tree\n");
return kFALSE;
}
if (formulaVar==NULL || formulaVar==0) {
::Error("AliNDLocalRegression::MakeFit", "Empty variable list\n");
return kFALSE;
}
if (formulaKernel==NULL) {
::Error("AliNDLocalRegression::MakeFit", "Kernel width not specified\n");
return kFALSE;
}
fUseBinNorm=useBinNorm;
//
fInputTree= tree; // should be better TRef?
fFormulaVal = new TObjString(formulaVal);
fFormulaVar = new TObjString(formulaVar);
fSelection = new TObjString(selection);
fKernelWidthFormula = new TObjString(formulaKernel);
fPolDimensionFormula = new TObjString(dimensionFormula);
TObjArray * arrayFormulaVar=fFormulaVar->String().Tokenize(":");
Int_t nvarFormula = arrayFormulaVar->GetEntries();
if (nvarFormula!=fHistPoints->GetNdimensions()){
::Error("AliNDLocalRegression::MakeFit", "Histogram/points mismatch\n");
return kFALSE;
}
TObjArray * arrayKernel=fKernelWidthFormula->String().Tokenize(":");
Int_t nwidthFormula = arrayKernel->GetEntries();
if (nvarFormula!=nwidthFormula){
delete arrayKernel;
delete arrayFormulaVar;
::Error("AliNDLocalRegression::MakeFit", "Variable/Kernel mismath\n");
return kFALSE;
}
fNParameters=nvarFormula;
//
// 2.) Load input data
//
//
Int_t entriesVal = tree->Draw(formulaVal,selection,"goffpara",entries);
if (entriesVal==0) {
::Error("AliNDLocalRegression::MakeFit", TString::Format("Empty point list\t%s\t%s\n",formulaVal,selection).Data());
return kFALSE;
}
if (tree->GetVal(0)==NULL || (tree->GetVal(1)==NULL)){
::Error("AliNDLocalRegression::MakeFit", TString::Format("Wrong selection\t%s\t%s\n",formulaVar,selection).Data());
return kFALSE;
}
TVectorD values(entriesVal,tree->GetVal(0));
TVectorD errors(entriesVal,tree->GetVal(1));
Double_t *pvalues = values.GetMatrixArray();
Double_t *perrors = errors.GetMatrixArray();
// 2.b) variables
TObjArray pointArray(fNParameters);
Int_t entriesVar = tree->Draw(formulaVar,selection,"goffpara",entries);
if (entriesVal!=entriesVar) {
::Error("AliNDLocalRegression::MakeFit", TString::Format("Wrong selection\t%s\t%s\n",formulaVar,selection).Data());
return kFALSE;
}
for (Int_t ipar=0; ipar<fNParameters; ipar++) pointArray.AddAt(new TVectorD(entriesVar,tree->GetVal(ipar)),ipar);
// 2.c) kernel array
TObjArray kernelArrayI2(fNParameters);
Int_t entriesKernel = tree->Draw(formulaKernel,selection,"goffpara",entries);
for (Int_t ipar=0; ipar<fNParameters; ipar++) {
TVectorD* vdI2 = new TVectorD(entriesVar,tree->GetVal(ipar));
for (int k=entriesVar;k--;) { // to speed up, precalculate inverse squared
double kv = (*vdI2)[k];
if (TMath::Abs(kv)<1e-12) ::Fatal("AliNDLocalRegression::MakeFit", "Kernel width=%f for entry %d of par:%d",kv,k,ipar);
(*vdI2)[k] = 1./(kv*kv);
}
kernelArrayI2.AddAt(vdI2,ipar);
}
//
Double_t *pvecVar[kMaxDim]={0};
Double_t *pvecKernelI2[kMaxDim]={0};
for (Int_t idim=0; idim<pointArray.GetEntries(); idim++){
pvecVar[idim]=((TVectorD*)(pointArray.At(idim)))->GetMatrixArray();
pvecKernelI2[idim]=((TVectorD*)(kernelArrayI2.At(idim)))->GetMatrixArray();
binRange[idim]=fHistPoints->GetAxis(idim)->GetNbins();
}
//
//
//
Int_t nbins = fHistPoints->GetNbins();
fBinIndex = new Int_t[fHistPoints->GetNdimensions()];
fBinCenter = new Double_t[fHistPoints->GetNdimensions()];
fBinDelta = new Double_t[fHistPoints->GetNdimensions()];
fBinWidth = new Double_t[fHistPoints->GetNdimensions()];
//
// 3.)
//
if (fCutType>0 && fRobustRMSLTSCut>0){
MakeRobustStatistic(values, errors, pointArray, kernelArrayI2, weightCut, fRobustFractionLTS);
}
//
// 4.) Make local fits
//
Double_t *binHypFit = new Double_t[2*fHistPoints->GetNdimensions()];
//
TLinearFitter fitter(1+2*fNParameters,TString::Format("hyp%d",2*fNParameters).Data());
for (Int_t ibin=0; ibin<nbins; ibin++) {
fHistPoints->GetBinContent(ibin,fBinIndex);
Bool_t isUnderFlowBin=kFALSE;
Bool_t isOverFlowBin=kFALSE;
for (Int_t idim=0; idim<fNParameters; idim++) {
if (fBinIndex[idim]==0) isUnderFlowBin=kTRUE;
if (fBinIndex[idim]>binRange[idim]) isOverFlowBin=kTRUE;
fBinCenter[idim]=fHistPoints->GetAxis(idim)->GetBinCenter(fBinIndex[idim]);
fBinWidth[idim]=fHistPoints->GetAxis(idim)->GetBinWidth(fBinIndex[idim]);
}
if (isUnderFlowBin || isOverFlowBin) continue;
fitter.ClearPoints();
// add fit points
for (Int_t ipoint=0; ipoint<entriesVal; ipoint++){
Double_t sumChi2=0;
if (fCutType>0 && fRobustRMSLTSCut>0){
Double_t localRMS=(*fLocalRobustStat)(ibin,2);
Double_t localMean=(*fLocalRobustStat)(ibin,1);
Double_t localMedian=(*fLocalRobustStat)(ibin,0);
if (fCutType==1){
if (TMath::Abs(pvalues[ipoint]-localMedian)>fRobustRMSLTSCut*localRMS) continue;
}
if (fCutType==2){
if (TMath::Abs(pvalues[ipoint]-localMean)>fRobustRMSLTSCut*localRMS) continue;
}
}
for (Int_t idim=0; idim<fNParameters; idim++){
//TVectorD &vecVar=*((TVectorD*)(pointArray.UncheckedAt(idim)));
//TVectorD &vecKernel=*((TVectorD*)(kernelArray.UncheckedAt(idim)));
fBinDelta[idim]=pvecVar[idim][ipoint]-fBinCenter[idim];
sumChi2+= (fBinDelta[idim]*fBinDelta[idim]) * pvecKernelI2[idim][ipoint];
if (sumChi2>nchi2Cut) break;//continue;
if (fUseBinNorm){
binHypFit[2*idim]=fBinDelta[idim]/fBinWidth[idim];
binHypFit[2*idim+1]=binHypFit[2*idim]*binHypFit[2*idim];
}else{
binHypFit[2*idim]=fBinDelta[idim];
binHypFit[2*idim+1]=fBinDelta[idim]*fBinDelta[idim];
}
}
if (sumChi2>nchi2Cut) continue;
// Double_t weight=TMath::Exp(-sumChi2*0.5);
// fitter.AddPoint(binHypFit,pvalues[ipoint], perrors[ipoint]/weight);
Double_t weightI=TMath::Exp(sumChi2*0.5);
fitter.AddPoint(binHypFit,pvalues[ipoint], perrors[ipoint]*weightI);
}
TVectorD * fitParam=new TVectorD(fNParameters*2+1);
TVectorD * fitQuality=new TVectorD(3);
TMatrixD * fitCovar=new TMatrixD(fNParameters*2+1,fNParameters*2+1);
Double_t normRMS=0;
Int_t nBinPoints=fitter.GetNpoints();
Bool_t fitOK=kFALSE;
(*fitQuality)[0]=0;
(*fitQuality)[1]=0;
(*fitQuality)[2]=0;
if (fitter.GetNpoints()>fNParameters*2+2){
fitOK = (fitter.Eval()==0);
if (fitOK){
normRMS=fitter.GetChisquare()/(fitter.GetNpoints()-fitter.GetNumberFreeParameters());
fitter.GetParameters(*fitParam);
fitter.GetCovarianceMatrix(*fitCovar);
(*fitQuality)[0]=nBinPoints;
(*fitQuality)[1]=normRMS;
(*fitQuality)[2]=ibin;
fLocalFitParam->AddAt(fitParam,ibin);
fLocalFitQuality->AddAt(fitQuality,ibin);
fLocalFitCovar->AddAt(fitCovar,ibin);
}
}
if (fStreamer){
TVectorD pfBinCenter(fNParameters, fBinCenter);
Double_t median=0,mean=0,rms=0;
if (fLocalRobustStat){
median=(*fLocalRobustStat)(ibin,0);
mean=(*fLocalRobustStat)(ibin,1);
rms=(*fLocalRobustStat)(ibin,2);
}
(*fStreamer)<<"localFit"<<
"ibin="<<ibin<< // bin index
"fitOK="<<fitOK<<
"localMedian="<<median<<
"localMean="<<mean<<
"localRMS="<<rms<<
"nBinPoints="<<nBinPoints<< // center of the bin
"binCenter.="<<&pfBinCenter<< //
"normRMS="<<normRMS<<
"fitParam.="<<fitParam<<
"fitCovar.="<<fitCovar<<
"fitOK="<<fitOK<<
"\n";
}
if (!fitOK) { // avoid memory leak for failed fits
delete fitParam;
delete fitQuality;
delete fitCovar;
}
}
return kTRUE;
}
Bool_t AliNDLocalRegression::MakeRobustStatistic(TVectorD &values,TVectorD &errors, TObjArray &pointArray, TObjArray &kernelArrayI2, Double_t weightCut, Double_t robustFraction){
//
// Calculate robust statistic information
// use raw array to make faster calcualtion
const Int_t kMaxDim=100;
Double_t *pvalues=values.GetMatrixArray();
Double_t *pvecVar[kMaxDim]={0};
Double_t *pvecKernelI2[kMaxDim]={0};
for (Int_t idim=0; idim<pointArray.GetEntries(); idim++){
pvecVar[idim]=((TVectorD*)(pointArray.At(idim)))->GetMatrixArray();
pvecKernelI2[idim]=((TVectorD*)(kernelArrayI2.At(idim)))->GetMatrixArray();
}
Double_t nchi2Cut=-2*TMath::Log(weightCut); // transform probability to nsigma cut
if (robustFraction>1) robustFraction=1;
Int_t nbins = fHistPoints->GetNbins(); //
Int_t npoints= values.GetNrows(); // number of points for fit
if (fLocalRobustStat){
delete fLocalRobustStat;
}
fLocalRobustStat=new TMatrixD(nbins,3);
TVectorD valueLocal(npoints);
for (Int_t ibin=0; ibin<nbins; ibin++){
fHistPoints->GetBinContent(ibin,fBinIndex); //
for (Int_t idim=0; idim<fNParameters; idim++){
fBinCenter[idim]=fHistPoints->GetAxis(idim)->GetBinCenter(fBinIndex[idim]);
fBinWidth[idim]=fHistPoints->GetAxis(idim)->GetBinWidth(fBinIndex[idim]);
}
Int_t indexLocal=0;
for (Int_t ipoint=0; ipoint<npoints; ipoint++){
Double_t sumChi2=0;
for (Int_t idim=0; idim<fNParameters; idim++){
//TVectorD &vecVar=*((TVectorD*)(pointArray.UncheckedAt(idim)));
//TVectorD &vecKernel=*((TVectorD*)(kernelArray.UncheckedAt(idim)));
fBinDelta[idim]=pvecVar[idim][ipoint]-fBinCenter[idim];
sumChi2+= (fBinDelta[idim]*fBinDelta[idim]) * pvecKernelI2[idim][ipoint];
if (sumChi2>nchi2Cut) break; //continue;
}
if (sumChi2>nchi2Cut) continue;
valueLocal[indexLocal]=pvalues[ipoint];
indexLocal++;
}
Double_t median=0,meanX=0, rmsX=0;
if (indexLocal*robustFraction-1>3){
median=TMath::Median(indexLocal,valueLocal.GetMatrixArray());
AliMathBase::EvaluateUni(indexLocal,valueLocal.GetMatrixArray(), meanX,rmsX, indexLocal*robustFraction-1);
}
(*fLocalRobustStat)(ibin,0)=median;
(*fLocalRobustStat)(ibin,1)=meanX;
(*fLocalRobustStat)(ibin,2)=rmsX;
}
}
Double_t AliNDLocalRegression::Eval(Double_t *point ){
//
//
//
const Double_t almost0=0.00000001;
// backward compatibility
if(!fBinWidth){
fBinWidth = new Double_t[fHistPoints->GetNdimensions()];
}
for (Int_t iDim=0; iDim<fNParameters; iDim++){
if (point[iDim]<= fHistPoints->GetAxis(iDim)->GetXmin()) point[iDim]=fHistPoints->GetAxis(iDim)->GetXmin()+almost0*fHistPoints->GetAxis(iDim)->GetBinWidth(0);
if (point[iDim]>= fHistPoints->GetAxis(iDim)->GetXmax()) point[iDim]=fHistPoints->GetAxis(iDim)->GetXmax()-almost0*fHistPoints->GetAxis(iDim)->GetBinWidth(0);
}
Int_t ibin = fHistPoints->GetBin(point);
Bool_t rangeOK=kTRUE;
if (ibin>=fLocalFitParam->GetEntriesFast() ){
rangeOK=kFALSE;
}else{
if (fLocalFitParam->UncheckedAt(ibin)==NULL) {
rangeOK=kFALSE;
}
}
if (!rangeOK) return 0;
fHistPoints->GetBinContent(ibin,fBinIndex);
for (Int_t idim=0; idim<fNParameters; idim++){
fBinCenter[idim]=fHistPoints->GetAxis(idim)->GetBinCenter(fBinIndex[idim]);
fBinWidth[idim]=fHistPoints->GetAxis(idim)->GetBinWidth(fBinIndex[idim]);
}
TVectorD &vecParam = *((TVectorD*)fLocalFitParam->At(ibin));
Double_t value=vecParam[0];
if (!rangeOK) return value;
if (fUseBinNorm) {
for (Int_t ipar=0; ipar<fNParameters; ipar++){
Double_t delta=(point[ipar]-fBinCenter[ipar])/fBinWidth[ipar];
value+=(vecParam[1+2*ipar]+vecParam[1+2*ipar+1]*delta)*delta;
}
} else {
for (Int_t ipar=0; ipar<fNParameters; ipar++){
Double_t delta=(point[ipar]-fBinCenter[ipar]);
value+=(vecParam[1+2*ipar]+vecParam[1+2*ipar+1]*delta)*delta;
}
}
return value;
}
Double_t AliNDLocalRegression::EvalError(Double_t *point ){
//
//
//
if (fLocalFitCovar==NULL) {
::Error("AliNDLocalRegression::EvalError", "Covariance matrix not available");
return 0 ;
}
for (Int_t iDim=0; iDim<fNParameters; iDim++){
if (point[iDim]< fHistPoints->GetAxis(iDim)->GetXmin()) point[iDim]=fHistPoints->GetAxis(iDim)->GetXmin();
if (point[iDim]> fHistPoints->GetAxis(iDim)->GetXmax()) point[iDim]=fHistPoints->GetAxis(iDim)->GetXmax();
}
Int_t ibin = fHistPoints->GetBin(point);
if (fLocalFitParam->At(ibin)==NULL) return 0;
fHistPoints->GetBinContent(ibin,fBinIndex);
for (Int_t idim=0; idim<fNParameters; idim++){
fBinCenter[idim]=fHistPoints->GetAxis(idim)->GetBinCenter(fBinIndex[idim]);
}
TMatrixD &vecCovar = *((TMatrixD*)fLocalFitCovar->At(ibin));
//TVectorD &vecQuality = *((TVectorD*)fLocalFitQuality->At(ibin));
Double_t value=TMath::Sqrt(vecCovar(0,0)); // fill covariance to be used
return value;
}
Bool_t AliNDLocalRegression::Derivative(Double_t *point, Double_t *d)
{
// fill d by partial derivatives
//
const Double_t almost0=0.00000001;
for (Int_t iDim=0; iDim<fNParameters; iDim++){
const TAxis* ax = fHistPoints->GetAxis(iDim);
if (point[iDim]<=ax->GetXmin()) point[iDim] = ax->GetXmin()+almost0*ax->GetBinWidth(0);
else if (point[iDim]>=ax->GetXmax()) point[iDim] = ax->GetXmax()-almost0*ax->GetBinWidth(0);
}
Int_t ibin = fHistPoints->GetBin(point);
if (ibin>=fLocalFitParam->GetEntriesFast() ||
!fLocalFitParam->UncheckedAt(ibin) ) return kFALSE;
//
fHistPoints->GetBinContent(ibin,fBinIndex);
for (Int_t idim=0; idim<fNParameters; idim++){
const TAxis* ax = fHistPoints->GetAxis(idim);
fBinCenter[idim] = ax->GetBinCenter(fBinIndex[idim]);
fBinWidth[idim] = ax->GetBinWidth(fBinIndex[idim]);
}
TVectorD &vecParam = *((TVectorD*)fLocalFitParam->At(ibin));
if (fUseBinNorm){
for (Int_t ipar=0; ipar<fNParameters; ipar++){
Double_t delta=(point[ipar]-fBinCenter[ipar])/fBinWidth[ipar];
d[ipar] = vecParam[1+2*ipar] + 2.0*vecParam[1+2*ipar+1]*delta;
}
}else{
for (Int_t ipar=0; ipar<fNParameters; ipar++){
Double_t delta=(point[ipar]-fBinCenter[ipar]);
d[ipar] = vecParam[1+2*ipar] + 2.0*vecParam[1+2*ipar+1]*delta;
}
}
return kTRUE;
}
Bool_t AliNDLocalRegression::EvalAndDerivative(Double_t *point, Double_t &val, Double_t *d)
{
// fill d by partial derivatives and calculate value in one go
//
const Double_t almost0=0.00000001;
for (Int_t iDim=0; iDim<fNParameters; iDim++){
const TAxis* ax = fHistPoints->GetAxis(iDim);
if (point[iDim]<=ax->GetXmin()) point[iDim] = ax->GetXmin()+almost0*ax->GetBinWidth(0);
else if (point[iDim]>=ax->GetXmax()) point[iDim] = ax->GetXmax()-almost0*ax->GetBinWidth(0);
}
val = 0;
Int_t ibin = fHistPoints->GetBin(point);
if (ibin>=fLocalFitParam->GetEntriesFast() ||
!fLocalFitParam->UncheckedAt(ibin) ) return kFALSE;
//
fHistPoints->GetBinContent(ibin,fBinIndex);
for (Int_t idim=0; idim<fNParameters; idim++){
const TAxis* ax = fHistPoints->GetAxis(idim);
fBinCenter[idim] = ax->GetBinCenter(fBinIndex[idim]);
fBinWidth[idim] = ax->GetBinWidth(fBinIndex[idim]);
}
TVectorD &vecParam = *((TVectorD*)fLocalFitParam->At(ibin));
val = vecParam[0];
if (fUseBinNorm){
for (Int_t ipar=0; ipar<fNParameters; ipar++){
Double_t delta=(point[ipar]-fBinCenter[ipar])/fBinWidth[ipar];
d[ipar] = vecParam[1+2*ipar] + 2.0*vecParam[1+2*ipar+1]*delta;
val += (vecParam[1+2*ipar]+vecParam[1+2*ipar+1]*delta)*delta;
}
} else {
for (Int_t ipar=0; ipar<fNParameters; ipar++){
Double_t delta=(point[ipar]-fBinCenter[ipar]);
d[ipar] = vecParam[1+2*ipar] + 2.0*vecParam[1+2*ipar+1]*delta;
val += (vecParam[1+2*ipar]+vecParam[1+2*ipar+1]*delta)*delta;
}
}
return kTRUE;
}
Int_t AliNDLocalRegression::GetVisualCorrectionIndex(const char *corName){
//
return TMath::Hash(corName)%1000000;
}
void AliNDLocalRegression::AddVisualCorrection(AliNDLocalRegression* corr, Int_t position){
/// make correction available for visualization using
/// TFormula, TFX and TTree::Draw
/// important in order to check corrections and also compute dervied variables
/// e.g correction partial derivatives
///
/// NOTE - class is not owner of correction
if (position==0) {
position=GetVisualCorrectionIndex(corr->GetName());
}
if (!fgVisualCorrection) fgVisualCorrection=new TObjArray(1000000);
if (position>=fgVisualCorrection->GetEntriesFast())
fgVisualCorrection->Expand((position+10)*2);
if (fgVisualCorrection->At(position)!=NULL){
::Error("AliNDLocalRegression::AddVisualCorrection","Correction %d already defined Old: %s New: %s",position,fgVisualCorrection->At(position)->GetName(), corr->GetName());
}
fgVisualCorrection->AddAt(corr, position);
}
AliNDLocalRegression* AliNDLocalRegression::GetVisualCorrection(Int_t position) {
/// Get visula correction registered at index=position
return fgVisualCorrection? (AliNDLocalRegression*)fgVisualCorrection->At(position):0;
}
Double_t AliNDLocalRegression::GetCorrND(Double_t index, Double_t par0){
//
//
AliNDLocalRegression *corr = (AliNDLocalRegression*)fgVisualCorrection->At(index);
if (!corr) return 0;
return corr->Eval(&par0);
}
Double_t AliNDLocalRegression::GetCorrNDError(Double_t index, Double_t par0){
//
//
AliNDLocalRegression *corr = (AliNDLocalRegression*)fgVisualCorrection->At(index);
if (!corr) return 0;
return corr->EvalError(&par0);
}
Double_t AliNDLocalRegression::GetCorrND(Double_t index, Double_t par0, Double_t par1){
//
//
AliNDLocalRegression *corr = (AliNDLocalRegression*)fgVisualCorrection->At(index);
if (!corr) return 0;
Double_t par[2]={par0,par1};
return corr->Eval(par);
}
Double_t AliNDLocalRegression::GetCorrNDError(Double_t index, Double_t par0, Double_t par1){
//
//
AliNDLocalRegression *corr = (AliNDLocalRegression*)fgVisualCorrection->At(index);
if (!corr) return 0;
Double_t par[2]={par0,par1};
return corr->EvalError(par);
}
Double_t AliNDLocalRegression::GetCorrND(Double_t index, Double_t par0, Double_t par1, Double_t par2){
//
//
AliNDLocalRegression *corr = (AliNDLocalRegression*)fgVisualCorrection->At(index);
if (!corr) return 0;
Double_t par[3]={par0,par1,par2};
return corr->Eval(par);
}
Double_t AliNDLocalRegression::GetCorrNDError(Double_t index, Double_t par0, Double_t par1, Double_t par2){
//
//
AliNDLocalRegression *corr = (AliNDLocalRegression*)fgVisualCorrection->At(index);
if (!corr) return 0;
Double_t par[3]={par0,par1,par2};
return corr->EvalError(par);
}
Double_t AliNDLocalRegression::GetCorrND(Double_t index, Double_t par0, Double_t par1, Double_t par2, Double_t par3){
//
//
AliNDLocalRegression *corr = (AliNDLocalRegression*)fgVisualCorrection->At(index);
if (!corr) return 0;
Double_t par[4]={par0,par1,par2,par3};
return corr->Eval(par);
}
Double_t AliNDLocalRegression::GetCorrNDError(Double_t index, Double_t par0, Double_t par1, Double_t par2, Double_t par3){
//
//
AliNDLocalRegression *corr = (AliNDLocalRegression*)fgVisualCorrection->At(index);
if (!corr) return 0;
Double_t par[4]={par0,par1,par2,par3};
return corr->EvalError(par);
}
Bool_t AliNDLocalRegression::AddWeekConstrainsAtBoundaries(Int_t nDims, Int_t *indexes, Double_t *relWeight, TTreeSRedirector* pcstream, Bool_t useCommon){
//
// Adding week constrain AtBoundaries
//
// Technique similar to "Kalman update" of measurement used at boundaries - https://en.wikipedia.org/wiki/Kalman_filter
//
// 1.) Make backup of original parameters
// 2.) Book Kalman matrices
// 3.) Loop over all measurements bins and update mesurements -adding boundary measurements as additional measurement
// relWeight vector specify relative weight of such measurement (err_i=sigma_i*refWeight_i) - not yet implemented
// 4.) replace original parameters with constrained parameters
// procedure can be repeated
/*
Input parameters example:
nDims=2;
Int_t indexes[2]={0,1};
Double_t relWeight0[6]={1,1,1,1,1,1};
Double_t relWeight1[6]={1,1,10,1,1,10};
pcstream=new TTreeSRedirector("constrainStream.root","recreate");
AliNDLocalRegression * regression0 = ( AliNDLocalRegression *)AliNDLocalRegression::GetVisualCorrections()->FindObject("pfitNDGaus0");
AliNDLocalRegression * regression1 = ( AliNDLocalRegression *)AliNDLocalRegression::GetVisualCorrections()->FindObject("pfitNDGaus1");
regressionUpdate0 = (AliNDLocalRegression *)regression0->Clone();
regressionUpdate1 = (AliNDLocalRegression *)regression1->Clone();
AddWeekConstrainsAtBoundaries( regressionUpdate0, nDims, indexes,relWeight0, pcstream);
AddWeekConstrainsAtBoundaries( regressionUpdate0, nDims, indexes,relWeight0, pcstream);
AddWeekConstrainsAtBoundaries( regressionUpdate0, nDims, indexes,relWeight0, pcstream);
AddWeekConstrainsAtBoundaries( regressionUpdate0, nDims, indexes,relWeight0, pcstream);
AddWeekConstrainsAtBoundaries( regressionUpdate1, nDims, indexes,relWeight1, pcstream);
AddWeekConstrainsAtBoundaries( regressionUpdate1, nDims, indexes,relWeight1, pcstream);
AddWeekConstrainsAtBoundaries( regressionUpdate1, nDims, indexes,relWeight1, pcstream);
AddWeekConstrainsAtBoundaries( regressionUpdate1, nDims, indexes,relWeight1, pcstream);
regressionUpdate0->SetName("pfitNDGaus0_Updated");
regressionUpdate1->SetName("pfitNDGaus1_Updated");
AliNDLocalRegression::AddVisualCorrection(regressionUpdate0);
AliNDLocalRegression::AddVisualCorrection(regressionUpdate1);
treeIn->SetAlias( regressionUpdate0->GetName(), TString::Format("AliNDLocalRegression::GetCorrND(%d,xyz0,xyz1+0)", regressionUpdate0->GetVisualCorrectionIndex()).Data());
treeIn->SetAlias( regressionUpdate1->GetName(), TString::Format("AliNDLocalRegression::GetCorrND(%d,xyz0,xyz1+0)", regressionUpdate1->GetVisualCorrectionIndex()).Data());
delete pcstream;
TFile *f = TFile::Open("constrainStream.root")
*/
const Double_t kScale=0.5;
const Double_t singularity_tolerance = 1e-200;
if (fLocalFitCovar==NULL) {
::Error("AliNDLocalRegression::EvalError", "Covariance matrix not available");
return 0 ;
}
//
// 1.) Make backup of original parameters
//
const TObjArray *vecParamOrig = fLocalFitParam;
const TObjArray *vecCovarOrig = fLocalFitCovar;
TObjArray *vecParamUpdated = new TObjArray(fHistPoints->GetNbins());
TObjArray *vecCovarUpdated = new TObjArray(fHistPoints->GetNbins());
//
// 2.) Book local varaibles and Kalman matrices
//
Int_t nParams= 1+2*fNParameters;
Int_t nMeas= nDims*6; // update each dimension specified 2 ends 2 measurements (value and first derivative)
TMatrixD matWeight(nParams,nParams); // weight matrix for side param
TMatrixD matCovarSide(nParams,nParams); // reweighted covariance matrix for side parameters
TMatrixD vecXk(nParams,1); // X vector - parameter of the local fit at bin
TMatrixD covXk(nParams,nParams); // X covariance
TMatrixD matHk(nMeas,nParams); // vector to mesurement (values at boundary of bin)
TMatrixD measR(nMeas,nMeas); // measurement error at boundary as provided by bin in local neigberhood
TMatrixD vecZk(nMeas,1); // measurement at boundary
//
TMatrixD measRBin(nMeas,nMeas); // measurement error bin
TMatrixD vecZkBin(nMeas,1); // measurement bin
TMatrixD matrixTransformBin(nMeas, nParams); // vector to measurement to calculate error matrix current bin
//
TMatrixD vecZkSide(3,1); // measurement side
TMatrixD matrixTransformSide(3,nParams);// vector to measurement to calculate error matrix side bin
//
TMatrixD vecYk(nMeas,1); // Innovation or measurement residual
TMatrixD matHkT(nParams,nMeas);
TMatrixD matSk(nMeas,nMeas); // Innovation (or residual) covariance
TMatrixD matKk(nParams,nMeas); // Optimal Kalman gain
TMatrixD mat1(nParams,nParams); // update covariance matrix
TMatrixD covXk2(nParams,nParams); //
TMatrixD covOut(nParams,nParams); //
mat1.UnitMatrix();
//
// 3.) Loop over all measurements bins and update mesurements -adding boundary measurements as additional measurement
// relWeight vector specify relative weight of such measurement (err_i=sigma_i*refWeight_i
const THn* his = GetHistogram();
Int_t binIndex[999]={0};
Int_t binIndexSide[999]={0};
Int_t nbinsAxis[999]={0};
Double_t binCenter[999]={0};
Double_t binWidth[999]={0};
if (relWeight!=NULL) for (Int_t iParam=0; iParam<nParams; iParam++){
Int_t index=0;
if (iParam<3) index=iParam;
if (iParam>=3) {
Int_t dim=(iParam-3)/2;
Int_t deriv=1+(iParam-3)%2;
index=3*dim+deriv;
}
matWeight(iParam,iParam)=relWeight[index];
}
for (Int_t iDim=0; iDim<nDims; iDim++){nbinsAxis[iDim]=his->GetAxis(iDim)->GetNbins();}
// Int_t nBins=fHistPoints->GetNbins();
Int_t nBins=fLocalFitParam->GetSize();
for (Int_t iBin=0; iBin<nBins; iBin++){ // loop over bins
if (iBin%fgVerboseLevel==0) printf("%d\n",iBin);
//
his->GetBinContent(iBin,binIndex);
for (Int_t iDim=0; iDim<nDims; iDim++) { // fill common info for bin of interest
binCenter[iDim]= his->GetAxis(iDim)->GetBinCenter(binIndex[iDim]);
binWidth[iDim] = his->GetAxis(iDim)->GetBinWidth(binIndex[iDim]);
}
if (fLocalFitParam->UncheckedAt(iBin)==NULL) continue;
Double_t *vecParam0 = ((TVectorD*)(fLocalFitParam->UncheckedAt(iBin)))->GetMatrixArray();
TMatrixD matParam0(nParams,1, vecParam0);
TMatrixD & matCovar0=*(((TMatrixD*)(fLocalFitCovar->UncheckedAt(iBin))));
measR.Zero();
vecZk.Zero();
measRBin.Zero();
vecZkBin.Zero();
matrixTransformBin.Zero();
covXk=matCovar0;
vecXk=matParam0;
//
// neiborhood loop
Int_t constCounter=0;
Int_t constCounterDim[100]={0};
for (Int_t iDim=0; iDim<nDims; iDim++){ // loop in n dim
constCounterDim[iDim]=0; // number of constraints per dimension
for (Int_t iSide=-1; iSide<=1; iSide+=2){ // left right loop
for (Int_t jDim=0; jDim<nDims; jDim++) binIndexSide[jDim]= binIndex[jDim];
vecZkSide.Zero();
matrixTransformSide.Zero();
//
binIndexSide[iDim]+=iSide;
if (binIndexSide[iDim]<0) binIndexSide[iDim]=0;
if (binIndexSide[iDim]>his->GetAxis(iDim)->GetNbins()) binIndexSide[iDim]=his->GetAxis(iDim)->GetNbins();
Bool_t isConst=binIndexSide[iDim]>0 &&binIndexSide[iDim]<=his->GetAxis(iDim)->GetNbins() && (fLocalFitParam)->UncheckedAt(his->GetBin(binIndexSide))!=NULL;
Int_t binSide=his->GetBin(binIndexSide);
if (binSide>=nBins ) binIndexSide[iDim]=binIndex[iDim];
if ((fLocalFitParam)->UncheckedAt(his->GetBin(binIndexSide))==NULL) binIndexSide[iDim]=binIndex[iDim];
if (isConst) {
constCounter++;
constCounterDim[iDim]++;
}
Double_t localCenter=his->GetAxis(iDim)->GetBinCenter(binIndex[iDim]);
Double_t sideCenter= his->GetAxis(iDim)->GetBinCenter(binIndexSide[iDim]);
Double_t position= (iSide<0) ? his->GetAxis(iDim)->GetBinLowEdge(binIndex[iDim]) : his->GetAxis(iDim)->GetBinUpEdge(binIndex[iDim]);
Double_t* vecParamSide = ((TVectorD*)(fLocalFitParam)->UncheckedAt(his->GetBin(binIndexSide)))->GetMatrixArray();
TMatrixD matParamSide(nParams,1, vecParamSide);
if (relWeight==NULL){
matCovarSide=*((TMatrixD*)(fLocalFitCovar->UncheckedAt(his->GetBin(binIndexSide))));
}
if (relWeight!=NULL){
matCovarSide=TMatrixD( matWeight,TMatrixD::kMult,*((TMatrixD*)(fLocalFitCovar->UncheckedAt(his->GetBin(binIndexSide)))));
matCovarSide*=matWeight;
}
if (!isConst) matCovarSide*=1000;
//
Double_t deltaLocal=(position-localCenter);
Double_t deltaSide=(position-sideCenter);
if (fUseBinNorm){
deltaLocal/=fBinWidth[iDim];
deltaSide/=fBinWidth[iDim];
}
//
matrixTransformSide(0,0)=1; matrixTransformSide(0,1+2*iDim)=deltaSide; matrixTransformSide(0,1+2*iDim+1)=deltaSide*deltaSide;
matrixTransformSide(1,1+2*iDim)=1; matrixTransformSide(1,1+2*iDim+1)=2*deltaSide;
matrixTransformSide(2,1+2*iDim+1)=2;
//
Int_t iMeas0=6*iDim+3*(iSide+1)/2;
matrixTransformBin(iMeas0+0,0)=1; matrixTransformBin(iMeas0+0,1+2*iDim)=deltaLocal; matrixTransformBin(iMeas0+0,1+2*iDim+1)=deltaSide*deltaLocal;
matrixTransformBin(iMeas0+1,1+2*iDim)=1; matrixTransformBin(iMeas0+1,1+2*iDim+1)=2*deltaLocal;
matrixTransformBin(iMeas0+2,1+2*iDim+1)=2;
//
for (Int_t iconst=0; iconst<3; iconst++){
Int_t iMeas=iMeas0+iconst;
Double_t localMeasurement=0;
Double_t sideMeasurement=0;
if (iconst==0){ // measurement - derivative 0
localMeasurement=vecParam0[0]+deltaLocal*(vecParam0[1+2*iDim]+vecParam0[2+2*iDim]*deltaLocal);
sideMeasurement=vecParamSide[0]+deltaSide*(vecParamSide[1+2*iDim]+vecParamSide[2+2*iDim]*deltaSide);
}
if (iconst==1){ // measurement -derivative 1
localMeasurement=(vecParam0[1+2*iDim]+2*vecParam0[2+2*iDim]*deltaLocal);
sideMeasurement=(vecParamSide[1+2*iDim]+2*vecParamSide[2+2*iDim]*deltaSide);
}
if (iconst==2){
localMeasurement=2*vecParam0[2+2*iDim];
sideMeasurement=2*vecParamSide[2+2*iDim];
}
vecZkSide(iconst,0)=sideMeasurement;
vecZk(iMeas,0)=sideMeasurement;
if (!isConst) vecZk(iMeas,0)=localMeasurement;
vecZkBin(iMeas,0)=localMeasurement;
}
TMatrixD measRSide0(matrixTransformSide,TMatrixD::kMult,matCovarSide); // (iconst,iconst) = (iconst,nParam)*(nParams,nParams)*(nParams,iconst
TMatrixD matrixTransformSideT(TMatrixD::kTransposed ,matrixTransformSide);
TMatrixD measRSide(measRSide0,TMatrixD::kMult,matrixTransformSideT);
// update measutement Covariance matrix for given side
for (Int_t iconst=0; iconst<3; iconst++)
for (Int_t jconst=0; jconst<3; jconst++){
measR(iMeas0+iconst,iMeas0+jconst)=measRSide(iconst,jconst);
}
if (pcstream){
TMatrixD vecZkSideCheck(matrixTransformSide,TMatrixD::kMult,matParamSide); // (iconst,1) = (iConst,nParam)*(nParams,1)
//
(*pcstream)<<"checkSide"<< // check agreement in 1D
"iBin="<<iBin<<
"iDim="<<iDim<<
"iSide="<<iSide<<
"vecZkSide.="<<&vecZkSide<<
"vecZkSideCheck.="<<&vecZkSideCheck<<
"measRSide.="<<&measRSide<<
"vecZk.="<<&vecZk<<
"vecZkBin.="<<&vecZkBin<<
"\n";
}
}
}
//
//
TMatrixD measRBin0(matrixTransformBin,TMatrixD::kMult,matCovar0); // (iconst,iconst) = (iconst,nParam)*(nParams,nParams)*(nParams,iconst