File indexing completed on 2023-03-17 11:25:57
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013 #include "Riostream.h"
0014 #include "SimTransport/TotemRPProtonTransportParametrization/interface/TMultiDimFet.h"
0015 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0016 #include "TMath.h"
0017 #include "TH1.h"
0018 #include "TH2.h"
0019 #include "TROOT.h"
0020 #include "TDecompChol.h"
0021 #include "TDatime.h"
0022 #include <iostream>
0023 #include <map>
0024
0025 #define RADDEG (180. / TMath::Pi())
0026 #define DEGRAD (TMath::Pi() / 180.)
0027 #define HIST_XORIG 0
0028 #define HIST_DORIG 1
0029 #define HIST_XNORM 2
0030 #define HIST_DSHIF 3
0031 #define HIST_RX 4
0032 #define HIST_RD 5
0033 #define HIST_RTRAI 6
0034 #define HIST_RTEST 7
0035 #define PARAM_MAXSTUDY 1
0036 #define PARAM_SEVERAL 2
0037 #define PARAM_RELERR 3
0038 #define PARAM_MAXTERMS 4
0039
0040
0041
0042
0043
0044 ClassImp(TMultiDimFet);
0045
0046
0047
0048
0049
0050
0051 TMultiDimFet::TMultiDimFet() {
0052
0053 fMeanQuantity = 0;
0054 fMaxQuantity = 0;
0055 fMinQuantity = 0;
0056 fSumSqQuantity = 0;
0057 fSumSqAvgQuantity = 0;
0058 fPowerLimit = 1;
0059
0060 fMaxAngle = 0;
0061 fMinAngle = 1;
0062
0063 fNVariables = 0;
0064 fMaxVariables = 0;
0065 fMinVariables = 0;
0066 fSampleSize = 0;
0067
0068 fMaxAngle = 0;
0069 fMinAngle = 0;
0070
0071 fPolyType = kMonomials;
0072 fShowCorrelation = kFALSE;
0073
0074 fIsUserFunction = kFALSE;
0075
0076 fHistograms = nullptr;
0077 fHistogramMask = 0;
0078
0079
0080
0081 }
0082
0083 const TMultiDimFet &TMultiDimFet::operator=(const TMultiDimFet &in) {
0084 if (this == &in) {
0085 return in;
0086 }
0087
0088 fMeanQuantity = in.fMeanQuantity;
0089
0090 fMaxQuantity = 0.0;
0091 fMinQuantity = 0.0;
0092 fSumSqQuantity = 0.0;
0093 fSumSqAvgQuantity = 0.0;
0094
0095 fNVariables = in.fNVariables;
0096
0097 fMaxVariables.ResizeTo(in.fMaxVariables.GetLwb(), in.fMaxVariables.GetUpb());
0098 fMaxVariables = in.fMaxVariables;
0099
0100 fMinVariables.ResizeTo(in.fMinVariables.GetLwb(), in.fMinVariables.GetUpb());
0101 fMinVariables = in.fMinVariables;
0102
0103 fSampleSize = 0;
0104 fTestSampleSize = 0;
0105 fMinAngle = 1;
0106 fMaxAngle = 0.0;
0107
0108 fMaxTerms = in.fMaxTerms;
0109
0110 fMinRelativeError = 0.0;
0111
0112 fMaxPowers.clear();
0113 fPowerLimit = 1;
0114
0115 fMaxFunctions = in.fMaxFunctions;
0116
0117 fFunctionCodes.clear();
0118 fMaxStudy = 0;
0119
0120 fMaxPowersFinal.clear();
0121
0122 fMaxFunctionsTimesNVariables = in.fMaxFunctionsTimesNVariables;
0123 fPowers = in.fPowers;
0124
0125 fPowerIndex = in.fPowerIndex;
0126
0127 fMaxResidual = 0.0;
0128 fMinResidual = 0.0;
0129 fMaxResidualRow = 0;
0130 fMinResidualRow = 0;
0131 fSumSqResidual = 0.0;
0132
0133 fNCoefficients = in.fNCoefficients;
0134
0135 fCoefficients.ResizeTo(in.fCoefficients.GetLwb(), in.fCoefficients.GetUpb());
0136 fCoefficients = in.fCoefficients;
0137
0138 fRMS = 0.0;
0139 fChi2 = 0.0;
0140 fParameterisationCode = 0;
0141 fError = 0.0;
0142 fTestError = 0.0;
0143 fPrecision = 0.0;
0144 fTestPrecision = 0.0;
0145 fCorrelationCoeff = 0.0;
0146 fTestCorrelationCoeff = 0.0;
0147 fHistograms = nullptr;
0148 fHistogramMask = 0;
0149
0150
0151 fPolyType = in.fPolyType;
0152 fShowCorrelation = in.fShowCorrelation;
0153 fIsUserFunction = in.fIsUserFunction;
0154 fIsVerbose = in.fIsVerbose;
0155 return in;
0156 }
0157
0158
0159 TMultiDimFet::TMultiDimFet(Int_t dimension, EMDFPolyType type, Option_t *option)
0160 : TNamed("multidimfit", "Multi-dimensional fit object"),
0161 fQuantity(dimension),
0162 fSqError(dimension),
0163 fVariables(dimension * 100),
0164 fMeanVariables(dimension),
0165 fMaxVariables(dimension),
0166 fMinVariables(dimension) {
0167
0168
0169
0170
0171
0172
0173
0174
0175
0176
0177
0178
0179
0180
0181
0182
0183 fMeanQuantity = 0;
0184 fMaxQuantity = 0;
0185 fMinQuantity = 0;
0186 fSumSqQuantity = 0;
0187 fSumSqAvgQuantity = 0;
0188 fPowerLimit = 1;
0189
0190 fMaxAngle = 0;
0191 fMinAngle = 1;
0192
0193 fNVariables = dimension;
0194 fMaxVariables = 0;
0195 fMaxFunctionsTimesNVariables = fMaxFunctions * fNVariables;
0196 fMinVariables = 0;
0197 fSampleSize = 0;
0198 fTestSampleSize = 0;
0199 fMinRelativeError = 0.01;
0200 fError = 0;
0201 fTestError = 0;
0202 fPrecision = 0;
0203 fTestPrecision = 0;
0204 fParameterisationCode = 0;
0205
0206 fPolyType = type;
0207 fShowCorrelation = kFALSE;
0208 fIsVerbose = kFALSE;
0209
0210 TString opt = option;
0211 opt.ToLower();
0212
0213 if (opt.Contains("k"))
0214 fShowCorrelation = kTRUE;
0215 if (opt.Contains("v"))
0216 fIsVerbose = kTRUE;
0217
0218 fIsUserFunction = kFALSE;
0219
0220 fHistograms = nullptr;
0221 fHistogramMask = 0;
0222
0223 fMaxPowers.resize(dimension);
0224 fMaxPowersFinal.resize(dimension);
0225
0226 }
0227
0228
0229 TMultiDimFet::~TMultiDimFet() {
0230 if (fHistograms)
0231 fHistograms->Clear("nodelete");
0232 delete fHistograms;
0233 }
0234
0235
0236 void TMultiDimFet::AddRow(const Double_t *x, Double_t D, Double_t E) {
0237
0238
0239
0240
0241
0242
0243
0244
0245
0246
0247
0248 if (!x)
0249 return;
0250
0251 if (++fSampleSize == 1) {
0252 fMeanQuantity = D;
0253 fMaxQuantity = D;
0254 fMinQuantity = D;
0255 } else {
0256 fMeanQuantity *= 1 - 1. / Double_t(fSampleSize);
0257 fMeanQuantity += D / Double_t(fSampleSize);
0258 fSumSqQuantity += D * D;
0259
0260 if (D >= fMaxQuantity)
0261 fMaxQuantity = D;
0262 if (D <= fMinQuantity)
0263 fMinQuantity = D;
0264 }
0265
0266
0267
0268 Int_t size = fQuantity.GetNrows();
0269 if (fSampleSize > size) {
0270 fQuantity.ResizeTo(size + size / 2);
0271 fSqError.ResizeTo(size + size / 2);
0272 }
0273
0274
0275 fQuantity(fSampleSize - 1) = D;
0276 fSqError(fSampleSize - 1) = (E == 0 ? D : E);
0277
0278
0279
0280
0281 size = fVariables.GetNrows();
0282 if (fSampleSize * fNVariables > size)
0283 fVariables.ResizeTo(size + size / 2);
0284
0285
0286 Int_t i, j;
0287 for (i = 0; i < fNVariables; i++) {
0288 if (fSampleSize == 1) {
0289 fMeanVariables(i) = x[i];
0290 fMaxVariables(i) = x[i];
0291 fMinVariables(i) = x[i];
0292 } else {
0293 fMeanVariables(i) *= 1 - 1. / Double_t(fSampleSize);
0294 fMeanVariables(i) += x[i] / Double_t(fSampleSize);
0295
0296
0297 if (x[i] >= fMaxVariables(i))
0298 fMaxVariables(i) = x[i];
0299
0300
0301 if (x[i] <= fMinVariables(i))
0302 fMinVariables(i) = x[i];
0303 }
0304
0305
0306 j = (fSampleSize - 1) * fNVariables + i;
0307 fVariables(j) = x[i];
0308 }
0309 }
0310
0311
0312 void TMultiDimFet::AddTestRow(const Double_t *x, Double_t D, Double_t E) {
0313
0314
0315
0316
0317
0318
0319
0320
0321
0322 if (fTestSampleSize++ == 0) {
0323 fTestQuantity.ResizeTo(fNVariables);
0324 fTestSqError.ResizeTo(fNVariables);
0325 fTestVariables.ResizeTo(fNVariables * 100);
0326 }
0327
0328
0329
0330 Int_t size = fTestQuantity.GetNrows();
0331 if (fTestSampleSize > size) {
0332 fTestQuantity.ResizeTo(size + size / 2);
0333 fTestSqError.ResizeTo(size + size / 2);
0334 }
0335
0336
0337 fTestQuantity(fTestSampleSize - 1) = D;
0338 fTestSqError(fTestSampleSize - 1) = (E == 0 ? D : E);
0339
0340
0341
0342
0343 size = fTestVariables.GetNrows();
0344 if (fTestSampleSize * fNVariables > size)
0345 fTestVariables.ResizeTo(size + size / 2);
0346
0347
0348 Int_t i, j;
0349 for (i = 0; i < fNVariables; i++) {
0350 j = fNVariables * (fTestSampleSize - 1) + i;
0351 fTestVariables(j) = x[i];
0352
0353 if (x[i] > fMaxVariables(i))
0354 Warning("AddTestRow", "variable %d (row: %d) too large: %f > %f", i, fTestSampleSize, x[i], fMaxVariables(i));
0355 if (x[i] < fMinVariables(i))
0356 Warning("AddTestRow", "variable %d (row: %d) too small: %f < %f", i, fTestSampleSize, x[i], fMinVariables(i));
0357 }
0358 }
0359
0360
0361 void TMultiDimFet::Clear(Option_t *option) {
0362
0363 Int_t i, j, n = fNVariables, m = fMaxFunctions;
0364
0365
0366 fQuantity.Zero();
0367 fSqError.Zero();
0368 fMeanQuantity = 0;
0369 fMaxQuantity = 0;
0370 fMinQuantity = 0;
0371 fSumSqQuantity = 0;
0372 fSumSqAvgQuantity = 0;
0373
0374
0375 fVariables.Zero();
0376 fNVariables = 0;
0377 fSampleSize = 0;
0378 fMeanVariables.Zero();
0379 fMaxVariables.Zero();
0380 fMinVariables.Zero();
0381
0382
0383 fTestQuantity.Zero();
0384 fTestSqError.Zero();
0385 fTestVariables.Zero();
0386 fTestSampleSize = 0;
0387
0388
0389 fFunctions.Zero();
0390 fMaxFunctions = 0;
0391 fMaxStudy = 0;
0392 fMaxFunctionsTimesNVariables = 0;
0393 fOrthFunctions.Zero();
0394 fOrthFunctionNorms.Zero();
0395
0396
0397 fMinRelativeError = 0;
0398 fMinAngle = 0;
0399 fMaxAngle = 0;
0400 fMaxTerms = 0;
0401
0402
0403 for (i = 0; i < n; i++) {
0404 fMaxPowers[i] = 0;
0405 fMaxPowersFinal[i] = 0;
0406 for (j = 0; j < m; j++)
0407 fPowers[i * n + j] = 0;
0408 }
0409 fPowerLimit = 0;
0410
0411
0412 fMaxResidual = 0;
0413 fMinResidual = 0;
0414 fMaxResidualRow = 0;
0415 fMinResidualRow = 0;
0416 fSumSqResidual = 0;
0417
0418
0419 fNCoefficients = 0;
0420 fOrthCoefficients = 0;
0421 fOrthCurvatureMatrix = 0;
0422 fRMS = 0;
0423 fCorrelationMatrix.Zero();
0424 fError = 0;
0425 fTestError = 0;
0426 fPrecision = 0;
0427 fTestPrecision = 0;
0428
0429
0430 fCoefficients.Zero();
0431 fCoefficientsRMS.Zero();
0432 fResiduals.Zero();
0433 fHistograms->Clear(option);
0434
0435
0436 fPolyType = kMonomials;
0437 fShowCorrelation = kFALSE;
0438 fIsUserFunction = kFALSE;
0439 }
0440
0441
0442 Double_t TMultiDimFet::Eval(const Double_t *x, const Double_t *coeff) const {
0443
0444
0445
0446
0447 Double_t returnValue = fMeanQuantity;
0448 Double_t term = 0;
0449 Int_t i, j;
0450
0451 for (i = 0; i < fNCoefficients; i++) {
0452
0453 term = (coeff ? coeff[i] : fCoefficients(i));
0454 for (j = 0; j < fNVariables; j++) {
0455
0456 Int_t p = fPowers[fPowerIndex[i] * fNVariables + j];
0457 Double_t y = 1 + 2. / (fMaxVariables(j) - fMinVariables(j)) * (x[j] - fMaxVariables(j));
0458 term *= EvalFactor(p, y);
0459 }
0460
0461 returnValue += term;
0462 }
0463 return returnValue;
0464 }
0465
0466 void TMultiDimFet::ReducePolynomial(double error) {
0467 if (error == 0.0)
0468 return;
0469 else {
0470 ZeroDoubiousCoefficients(error);
0471 }
0472 }
0473
0474 void TMultiDimFet::ZeroDoubiousCoefficients(double error) {
0475 typedef std::multimap<double, int> cmt;
0476 cmt m;
0477
0478 for (int i = 0; i < fNCoefficients; i++) {
0479 m.insert(std::pair<double, int>(TMath::Abs(fCoefficients(i)), i));
0480 }
0481
0482 double del_error_abs = 0;
0483 int deleted_terms_count = 0;
0484
0485 for (cmt::iterator it = m.begin(); it != m.end() && del_error_abs < error; ++it) {
0486 if (TMath::Abs(it->first) + del_error_abs < error) {
0487 fCoefficients(it->second) = 0.0;
0488 del_error_abs = TMath::Abs(it->first) + del_error_abs;
0489 deleted_terms_count++;
0490 } else
0491 break;
0492 }
0493
0494 int fNCoefficients_new = fNCoefficients - deleted_terms_count;
0495 TVectorD fCoefficients_new(fNCoefficients_new);
0496 std::vector<Int_t> fPowerIndex_new;
0497
0498 int ind = 0;
0499 for (int i = 0; i < fNCoefficients; i++) {
0500 if (fCoefficients(i) != 0.0) {
0501 fCoefficients_new(ind) = fCoefficients(i);
0502 fPowerIndex_new.push_back(fPowerIndex[i]);
0503 ind++;
0504 }
0505 }
0506 fNCoefficients = fNCoefficients_new;
0507 fCoefficients.ResizeTo(fNCoefficients);
0508 fCoefficients = fCoefficients_new;
0509 fPowerIndex = fPowerIndex_new;
0510 edm::LogInfo("TMultiDimFet") << deleted_terms_count << " terms removed"
0511 << "\n";
0512 }
0513
0514
0515 Double_t TMultiDimFet::EvalControl(const Int_t *iv) {
0516
0517
0518 Double_t s = 0;
0519 Double_t epsilon = 1e-6;
0520 for (Int_t i = 0; i < fNVariables; i++) {
0521 if (fMaxPowers[i] != 1)
0522 s += (epsilon + iv[i] - 1) / (epsilon + fMaxPowers[i] - 1);
0523 }
0524 return s;
0525 }
0526
0527
0528 Double_t TMultiDimFet::EvalFactor(Int_t p, Double_t x) const {
0529
0530
0531 Int_t i = 0;
0532 Double_t p1 = 1;
0533 Double_t p2 = 0;
0534 Double_t p3 = 0;
0535 Double_t r = 0;
0536
0537 switch (p) {
0538 case 1:
0539 r = 1;
0540 break;
0541 case 2:
0542 r = x;
0543 break;
0544 default:
0545 p2 = x;
0546 for (i = 3; i <= p; i++) {
0547 p3 = p2 * x;
0548 if (fPolyType == kLegendre)
0549 p3 = ((2 * i - 3) * p2 * x - (i - 2) * p1) / (i - 1);
0550 else if (fPolyType == kChebyshev)
0551 p3 = 2 * x * p2 - p1;
0552 p1 = p2;
0553 p2 = p3;
0554 }
0555 r = p3;
0556 }
0557
0558 return r;
0559 }
0560
0561
0562 void TMultiDimFet::FindParameterization(double precision) {
0563
0564
0565
0566
0567
0568
0569
0570 MakeNormalized();
0571 MakeCandidates();
0572 MakeParameterization();
0573 MakeCoefficients();
0574 ReducePolynomial(precision);
0575 }
0576
0577
0578
0579
0580
0581
0582
0583
0584
0585
0586
0587
0588
0589
0590
0591
0592
0593
0594
0595
0596
0597
0598
0599
0600
0601
0602
0603
0604
0605
0606
0607
0608
0609
0610
0611
0612
0613
0614
0615
0616
0617
0618
0619
0620
0621
0622
0623
0624
0625
0626
0627
0628
0629
0630
0631
0632
0633
0634
0635
0636
0637
0638
0639
0640
0641
0642
0643
0644
0645
0646
0647
0648
0649
0650
0651
0652
0653
0654
0655
0656
0657
0658 void TMultiDimFet::MakeCandidates() {
0659
0660
0661
0662
0663 Int_t i = 0;
0664 Int_t j = 0;
0665 Int_t k = 0;
0666
0667
0668
0669 Int_t *powers = new Int_t[fNVariables * fMaxFunctions];
0670
0671
0672 Double_t *control = new Double_t[fMaxFunctions];
0673
0674
0675 Int_t *iv = new Int_t[fNVariables];
0676 for (i = 0; i < fNVariables; i++)
0677 iv[i] = 1;
0678
0679 if (!fIsUserFunction) {
0680
0681 Int_t numberFunctions = 0;
0682
0683 while (kTRUE) {
0684
0685 Double_t s = EvalControl(iv);
0686
0687 if (s <= fPowerLimit) {
0688
0689
0690 if (Select(iv)) {
0691 numberFunctions++;
0692
0693
0694
0695 if (numberFunctions > fMaxFunctions)
0696 break;
0697
0698
0699
0700 control[numberFunctions - 1] = Int_t(1.0e+6 * s);
0701
0702
0703 for (i = 0; i < fNVariables; i++) {
0704 j = (numberFunctions - 1) * fNVariables + i;
0705 powers[j] = iv[i];
0706 }
0707 }
0708 }
0709
0710 for (i = 0; i < fNVariables; i++)
0711 if (iv[i] < fMaxPowers[i])
0712 break;
0713
0714
0715
0716 if (i == fNVariables) {
0717 fMaxFunctions = numberFunctions;
0718 fMaxFunctionsTimesNVariables = fMaxFunctions * fNVariables;
0719 break;
0720 }
0721
0722
0723 iv[i]++;
0724
0725 for (j = 0; j < i; j++)
0726 iv[j] = 1;
0727 }
0728 } else {
0729
0730 for (i = 0; i < fMaxFunctions; i++) {
0731
0732 for (j = 0; j < fNVariables; j++) {
0733 powers[i * fNVariables + j] = fPowers[i * fNVariables + j];
0734 iv[j] = fPowers[i * fNVariables + j];
0735 }
0736
0737 control[i] = Int_t(1.0e+6 * EvalControl(iv));
0738 }
0739 }
0740
0741
0742
0743 Int_t *order = new Int_t[fMaxFunctions];
0744 for (i = 0; i < fMaxFunctions; i++)
0745 order[i] = i;
0746 fPowers.resize(fMaxFunctions * fNVariables);
0747
0748 for (i = 0; i < fMaxFunctions; i++) {
0749 Double_t x = control[i];
0750 Int_t l = order[i];
0751 k = i;
0752
0753 for (j = i; j < fMaxFunctions; j++) {
0754 if (control[j] <= x) {
0755 x = control[j];
0756 l = order[j];
0757 k = j;
0758 }
0759 }
0760
0761 if (k != i) {
0762 control[k] = control[i];
0763 control[i] = x;
0764 order[k] = order[i];
0765 order[i] = l;
0766 }
0767 }
0768
0769 for (i = 0; i < fMaxFunctions; i++)
0770 for (j = 0; j < fNVariables; j++)
0771 fPowers[i * fNVariables + j] = powers[order[i] * fNVariables + j];
0772
0773 delete[] control;
0774 delete[] powers;
0775 delete[] order;
0776 delete[] iv;
0777 }
0778
0779 Double_t TMultiDimFet::MakeChi2(const Double_t *coeff) {
0780
0781
0782
0783
0784
0785 fChi2 = 0;
0786 Int_t i, j;
0787 Double_t *x = new Double_t[fNVariables];
0788 for (i = 0; i < fTestSampleSize; i++) {
0789
0790 for (j = 0; j < fNVariables; j++)
0791 x[j] = fTestVariables(i * fNVariables + j);
0792
0793
0794 Double_t f = Eval(x, coeff);
0795
0796
0797 fChi2 += 1. / TMath::Max(fTestSqError(i), 1e-20) * (fTestQuantity(i) - f) * (fTestQuantity(i) - f);
0798 }
0799
0800
0801 delete[] x;
0802
0803 return fChi2;
0804 }
0805
0806
0807 void TMultiDimFet::MakeCode(const char *filename, Option_t *option) {
0808
0809
0810
0811
0812
0813
0814
0815
0816
0817
0818
0819
0820
0821
0822
0823
0824
0825
0826
0827
0828
0829
0830
0831
0832 TString outName(filename);
0833 if (!outName.EndsWith(".C") && !outName.EndsWith(".cxx"))
0834 outName += ".C";
0835
0836 MakeRealCode(outName.Data(), "", option);
0837 }
0838
0839 void TMultiDimFet::MakeCoefficientErrors() {
0840
0841
0842
0843 Int_t i = 0;
0844 Int_t j = 0;
0845 Int_t k = 0;
0846 TVectorD iF(fSampleSize);
0847 TVectorD jF(fSampleSize);
0848 fCoefficientsRMS.ResizeTo(fNCoefficients);
0849
0850 TMatrixDSym curvatureMatrix(fNCoefficients);
0851
0852
0853 for (i = 0; i < fNCoefficients; i++) {
0854 iF = TMatrixDRow(fFunctions, i);
0855 for (j = 0; j <= i; j++) {
0856 jF = TMatrixDRow(fFunctions, j);
0857 for (k = 0; k < fSampleSize; k++)
0858 curvatureMatrix(i, j) += 1 / TMath::Max(fSqError(k), 1e-20) * iF(k) * jF(k);
0859 curvatureMatrix(j, i) = curvatureMatrix(i, j);
0860 }
0861 }
0862
0863
0864 fChi2 = 0;
0865 for (i = 0; i < fSampleSize; i++) {
0866 Double_t f = 0;
0867 for (j = 0; j < fNCoefficients; j++)
0868 f += fCoefficients(j) * fFunctions(j, i);
0869 fChi2 += 1. / TMath::Max(fSqError(i), 1e-20) * (fQuantity(i) - f) * (fQuantity(i) - f);
0870 }
0871
0872
0873 const TVectorD diag = TMatrixDDiag_const(curvatureMatrix);
0874 curvatureMatrix.NormByDiag(diag);
0875
0876 TDecompChol chol(curvatureMatrix);
0877 if (!chol.Decompose())
0878 Error("MakeCoefficientErrors", "curvature matrix is singular");
0879 chol.Invert(curvatureMatrix);
0880
0881 curvatureMatrix.NormByDiag(diag);
0882
0883 for (i = 0; i < fNCoefficients; i++)
0884 fCoefficientsRMS(i) = TMath::Sqrt(curvatureMatrix(i, i));
0885 }
0886
0887
0888 void TMultiDimFet::MakeCoefficients() {
0889
0890
0891
0892
0893
0894
0895
0896
0897 Int_t i = 0, j = 0;
0898 Int_t col = 0, row = 0;
0899
0900
0901 for (col = 1; col < fNCoefficients; col++) {
0902 for (row = col - 1; row > -1; row--) {
0903 fOrthCurvatureMatrix(row, col) = 0;
0904 for (i = row; i <= col; i++)
0905 fOrthCurvatureMatrix(row, col) -= fOrthCurvatureMatrix(i, row) * fOrthCurvatureMatrix(i, col);
0906 }
0907 }
0908
0909
0910 fCoefficients.ResizeTo(fNCoefficients);
0911
0912 for (i = 0; i < fNCoefficients; i++) {
0913 Double_t sum = 0;
0914 for (j = i; j < fNCoefficients; j++)
0915 sum += fOrthCurvatureMatrix(i, j) * fOrthCoefficients(j);
0916 fCoefficients(i) = sum;
0917 }
0918
0919
0920 fResiduals.ResizeTo(fSampleSize);
0921 for (i = 0; i < fSampleSize; i++)
0922 fResiduals(i) = fQuantity(i);
0923
0924 for (i = 0; i < fNCoefficients; i++)
0925 for (j = 0; j < fSampleSize; j++)
0926 fResiduals(j) -= fCoefficients(i) * fFunctions(i, j);
0927
0928
0929
0930 fMinResidual = 10e10;
0931 fMaxResidual = -10e10;
0932 Double_t sqRes = 0;
0933 for (i = 0; i < fSampleSize; i++) {
0934 sqRes += fResiduals(i) * fResiduals(i);
0935 if (fResiduals(i) <= fMinResidual) {
0936 fMinResidual = fResiduals(i);
0937 fMinResidualRow = i;
0938 }
0939 if (fResiduals(i) >= fMaxResidual) {
0940 fMaxResidual = fResiduals(i);
0941 fMaxResidualRow = i;
0942 }
0943 }
0944
0945 fCorrelationCoeff = fSumSqResidual / fSumSqAvgQuantity;
0946 fPrecision = TMath::Sqrt(sqRes / fSumSqQuantity);
0947
0948
0949 if (TESTBIT(fHistogramMask, HIST_RD) || TESTBIT(fHistogramMask, HIST_RTRAI) || TESTBIT(fHistogramMask, HIST_RX)) {
0950 for (i = 0; i < fSampleSize; i++) {
0951 if (TESTBIT(fHistogramMask, HIST_RD))
0952 ((TH2D *)fHistograms->FindObject("res_d"))->Fill(fQuantity(i), fResiduals(i));
0953 if (TESTBIT(fHistogramMask, HIST_RTRAI))
0954 ((TH1D *)fHistograms->FindObject("res_train"))->Fill(fResiduals(i));
0955
0956 if (TESTBIT(fHistogramMask, HIST_RX))
0957 for (j = 0; j < fNVariables; j++)
0958 ((TH2D *)fHistograms->FindObject(Form("res_x_%d", j)))->Fill(fVariables(i * fNVariables + j), fResiduals(i));
0959 }
0960 }
0961 }
0962
0963
0964 void TMultiDimFet::MakeCorrelation() {
0965
0966
0967 if (!fShowCorrelation)
0968 return;
0969
0970 fCorrelationMatrix.ResizeTo(fNVariables, fNVariables + 1);
0971
0972 Double_t d2 = 0;
0973 Double_t ddotXi = 0;
0974 Double_t xiNorm = 0;
0975 Double_t xidotXj = 0;
0976 Double_t xjNorm = 0;
0977
0978 Int_t i, j, k, l, m;
0979 for (i = 0; i < fSampleSize; i++)
0980 d2 += fQuantity(i) * fQuantity(i);
0981
0982 for (i = 0; i < fNVariables; i++) {
0983 ddotXi = 0.;
0984 xiNorm = 0.;
0985 for (j = 0; j < fSampleSize; j++) {
0986
0987 k = j * fNVariables + i;
0988 ddotXi += fQuantity(j) * (fVariables(k) - fMeanVariables(i));
0989 xiNorm += (fVariables(k) - fMeanVariables(i)) * (fVariables(k) - fMeanVariables(i));
0990 }
0991 fCorrelationMatrix(i, 0) = ddotXi / TMath::Sqrt(d2 * xiNorm);
0992
0993 for (j = 0; j < i; j++) {
0994 xidotXj = 0.;
0995 xjNorm = 0.;
0996 for (k = 0; k < fSampleSize; k++) {
0997
0998
0999 l = k * fNVariables + j;
1000 m = k * fNVariables + i;
1001
1002
1003 xidotXj +=
1004 (fVariables(m) - fMeanVariables(i)) * (fVariables(l) - fMeanVariables(j));
1005 xjNorm += (fVariables(l) - fMeanVariables(j)) * (fVariables(l) - fMeanVariables(j));
1006 }
1007 fCorrelationMatrix(i, j + 1) = xidotXj / TMath::Sqrt(xiNorm * xjNorm);
1008 }
1009 }
1010 }
1011
1012
1013 Double_t TMultiDimFet::MakeGramSchmidt(Int_t function) {
1014
1015
1016
1017
1018
1019
1020
1021
1022 Double_t f2 = 0;
1023 fOrthCoefficients(fNCoefficients) = 0;
1024 fOrthFunctionNorms(fNCoefficients) = 0;
1025 Int_t j = 0;
1026 Int_t k = 0;
1027
1028 for (j = 0; j < fSampleSize; j++) {
1029 fFunctions(fNCoefficients, j) = 1;
1030 fOrthFunctions(fNCoefficients, j) = 0;
1031
1032 for (k = 0; k < fNVariables; k++) {
1033 Int_t p = fPowers[function * fNVariables + k];
1034 Double_t x = fVariables(j * fNVariables + k);
1035 fFunctions(fNCoefficients, j) *= EvalFactor(p, x);
1036 }
1037
1038
1039 f2 += fFunctions(fNCoefficients, j) * fFunctions(fNCoefficients, j);
1040
1041 fOrthFunctions(fNCoefficients, j) = fFunctions(fNCoefficients, j);
1042 }
1043
1044
1045 for (j = 0; j < fNCoefficients; j++) {
1046 Double_t fdw = 0;
1047
1048 for (k = 0; k < fSampleSize; k++) {
1049 fdw += fFunctions(fNCoefficients, k) * fOrthFunctions(j, k) / fOrthFunctionNorms(j);
1050 }
1051
1052 fOrthCurvatureMatrix(fNCoefficients, j) = fdw;
1053
1054 for (k = 0; k < fSampleSize; k++)
1055 fOrthFunctions(fNCoefficients, k) -= fdw * fOrthFunctions(j, k);
1056 }
1057
1058 for (j = 0; j < fSampleSize; j++) {
1059
1060 fOrthFunctionNorms(fNCoefficients) += fOrthFunctions(fNCoefficients, j) * fOrthFunctions(fNCoefficients, j);
1061
1062
1063 fOrthCoefficients(fNCoefficients) += fQuantity(j) * fOrthFunctions(fNCoefficients, j);
1064 }
1065
1066
1067 if (!fIsUserFunction)
1068 if (TMath::Sqrt(fOrthFunctionNorms(fNCoefficients) / (f2 + 1e-10)) < TMath::Sin(fMinAngle * DEGRAD))
1069 return 0;
1070
1071
1072
1073
1074
1075 fOrthCurvatureMatrix(fNCoefficients, fNCoefficients) = 1;
1076 Double_t b = fOrthCoefficients(fNCoefficients);
1077 fOrthCoefficients(fNCoefficients) /= fOrthFunctionNorms(fNCoefficients);
1078
1079
1080 Double_t dResidur = fOrthCoefficients(fNCoefficients) * b;
1081
1082 return dResidur;
1083 }
1084
1085
1086 void TMultiDimFet::MakeHistograms(Option_t *option) {
1087
1088
1089
1090
1091
1092
1093
1094
1095
1096
1097
1098
1099
1100
1101
1102
1103
1104 TString opt(option);
1105 opt.ToLower();
1106
1107 if (opt.Length() < 1)
1108 return;
1109
1110 if (!fHistograms)
1111 fHistograms = new TList;
1112
1113
1114 Int_t i = 0;
1115
1116
1117 if (opt.Contains("x") || opt.Contains("a")) {
1118 SETBIT(fHistogramMask, HIST_XORIG);
1119 for (i = 0; i < fNVariables; i++)
1120 if (!fHistograms->FindObject(Form("x_%d_orig", i)))
1121 fHistograms->Add(
1122 new TH1D(Form("x_%d_orig", i), Form("Original variable # %d", i), 100, fMinVariables(i), fMaxVariables(i)));
1123 }
1124
1125
1126 if (opt.Contains("d") || opt.Contains("a")) {
1127 SETBIT(fHistogramMask, HIST_DORIG);
1128 if (!fHistograms->FindObject("d_orig"))
1129 fHistograms->Add(new TH1D("d_orig", "Original Quantity", 100, fMinQuantity, fMaxQuantity));
1130 }
1131
1132
1133 if (opt.Contains("n") || opt.Contains("a")) {
1134 SETBIT(fHistogramMask, HIST_XNORM);
1135 for (i = 0; i < fNVariables; i++)
1136 if (!fHistograms->FindObject(Form("x_%d_norm", i)))
1137 fHistograms->Add(new TH1D(Form("x_%d_norm", i), Form("Normalized variable # %d", i), 100, -1, 1));
1138 }
1139
1140
1141 if (opt.Contains("s") || opt.Contains("a")) {
1142 SETBIT(fHistogramMask, HIST_DSHIF);
1143 if (!fHistograms->FindObject("d_shifted"))
1144 fHistograms->Add(
1145 new TH1D("d_shifted", "Shifted Quantity", 100, fMinQuantity - fMeanQuantity, fMaxQuantity - fMeanQuantity));
1146 }
1147
1148
1149 if (opt.Contains("r1") || opt.Contains("a")) {
1150 SETBIT(fHistogramMask, HIST_RX);
1151 for (i = 0; i < fNVariables; i++)
1152 if (!fHistograms->FindObject(Form("res_x_%d", i)))
1153 fHistograms->Add(new TH2D(Form("res_x_%d", i),
1154 Form("Computed residual versus x_%d", i),
1155 100,
1156 -1,
1157 1,
1158 35,
1159 fMinQuantity - fMeanQuantity,
1160 fMaxQuantity - fMeanQuantity));
1161 }
1162
1163
1164 if (opt.Contains("r2") || opt.Contains("a")) {
1165 SETBIT(fHistogramMask, HIST_RD);
1166 if (!fHistograms->FindObject("res_d"))
1167 fHistograms->Add(new TH2D("res_d",
1168 "Computed residuals vs Quantity",
1169 100,
1170 fMinQuantity - fMeanQuantity,
1171 fMaxQuantity - fMeanQuantity,
1172 35,
1173 fMinQuantity - fMeanQuantity,
1174 fMaxQuantity - fMeanQuantity));
1175 }
1176
1177
1178 if (opt.Contains("r3") || opt.Contains("a")) {
1179 SETBIT(fHistogramMask, HIST_RTRAI);
1180 if (!fHistograms->FindObject("res_train"))
1181 fHistograms->Add(new TH1D("res_train",
1182 "Computed residuals over training sample",
1183 100,
1184 fMinQuantity - fMeanQuantity,
1185 fMaxQuantity - fMeanQuantity));
1186 }
1187 if (opt.Contains("r4") || opt.Contains("a")) {
1188 SETBIT(fHistogramMask, HIST_RTEST);
1189 if (!fHistograms->FindObject("res_test"))
1190 fHistograms->Add(new TH1D("res_test",
1191 "Distribution of residuals from test",
1192 100,
1193 fMinQuantity - fMeanQuantity,
1194 fMaxQuantity - fMeanQuantity));
1195 }
1196 }
1197
1198
1199 void TMultiDimFet::MakeMethod(const Char_t *classname, Option_t *option) {
1200
1201
1202
1203
1204
1205
1206
1207
1208
1209
1210
1211
1212
1213
1214
1215
1216
1217
1218
1219
1220
1221
1222
1223
1224
1225
1226
1227
1228
1229
1230
1231
1232
1233
1234
1235
1236
1237
1238
1239
1240
1241
1242
1243
1244 MakeRealCode(Form("%sMDF.cxx", classname), classname, option);
1245 }
1246
1247
1248 void TMultiDimFet::MakeNormalized() {
1249
1250
1251
1252
1253 Int_t i = 0;
1254 Int_t j = 0;
1255 Int_t k = 0;
1256
1257 for (i = 0; i < fSampleSize; i++) {
1258 if (TESTBIT(fHistogramMask, HIST_DORIG))
1259 ((TH1D *)fHistograms->FindObject("d_orig"))->Fill(fQuantity(i));
1260
1261 fQuantity(i) -= fMeanQuantity;
1262 fSumSqAvgQuantity += fQuantity(i) * fQuantity(i);
1263
1264 if (TESTBIT(fHistogramMask, HIST_DSHIF))
1265 ((TH1D *)fHistograms->FindObject("d_shifted"))->Fill(fQuantity(i));
1266
1267 for (j = 0; j < fNVariables; j++) {
1268 Double_t range = 1. / (fMaxVariables(j) - fMinVariables(j));
1269 k = i * fNVariables + j;
1270
1271
1272 if (TESTBIT(fHistogramMask, HIST_XORIG))
1273 ((TH1D *)fHistograms->FindObject(Form("x_%d_orig", j)))->Fill(fVariables(k));
1274
1275
1276 fVariables(k) = 1 + 2 * range * (fVariables(k) - fMaxVariables(j));
1277
1278
1279 if (TESTBIT(fHistogramMask, HIST_XNORM))
1280 ((TH1D *)fHistograms->FindObject(Form("x_%d_norm", j)))->Fill(fVariables(k));
1281 }
1282 }
1283
1284 fMaxQuantity -= fMeanQuantity;
1285 fMinQuantity -= fMeanQuantity;
1286
1287
1288 for (i = 0; i < fNVariables; i++) {
1289 Double_t range = 1. / (fMaxVariables(i) - fMinVariables(i));
1290 fMeanVariables(i) = 1 + 2 * range * (fMeanVariables(i) - fMaxVariables(i));
1291 }
1292 }
1293
1294
1295 void TMultiDimFet::MakeParameterization() {
1296
1297
1298
1299
1300
1301 Int_t i = -1;
1302 Int_t j = 0;
1303 Int_t k = 0;
1304 Int_t maxPass = 3;
1305 Int_t studied = 0;
1306 Double_t squareResidual = fSumSqAvgQuantity;
1307 fNCoefficients = 0;
1308 fSumSqResidual = fSumSqAvgQuantity;
1309 fFunctions.ResizeTo(fMaxTerms, fSampleSize);
1310 fOrthFunctions.ResizeTo(fMaxTerms, fSampleSize);
1311 fOrthFunctionNorms.ResizeTo(fMaxTerms);
1312 fOrthCoefficients.ResizeTo(fMaxTerms);
1313 fOrthCurvatureMatrix.ResizeTo(fMaxTerms, fMaxTerms);
1314 fFunctions = 1;
1315
1316 fFunctionCodes.resize(fMaxFunctions);
1317 fPowerIndex.resize(fMaxTerms);
1318 Int_t l;
1319 for (l = 0; l < fMaxFunctions; l++)
1320 fFunctionCodes[l] = 0;
1321 for (l = 0; l < fMaxTerms; l++)
1322 fPowerIndex[l] = 0;
1323
1324 if (fMaxAngle != 0)
1325 maxPass = 100;
1326 if (fIsUserFunction)
1327 maxPass = 1;
1328
1329
1330
1331 while (kTRUE) {
1332
1333 if (studied++ >= fMaxStudy) {
1334 fParameterisationCode = PARAM_MAXSTUDY;
1335 break;
1336 }
1337
1338
1339 if (k >= maxPass) {
1340 fParameterisationCode = PARAM_SEVERAL;
1341 break;
1342 }
1343
1344
1345 i++;
1346
1347
1348 if (i == fMaxFunctions) {
1349 if (fMaxAngle != 0)
1350 fMaxAngle += (90 - fMaxAngle) / 2;
1351 i = 0;
1352 studied--;
1353 k++;
1354 continue;
1355 }
1356 if (studied == 1)
1357 fFunctionCodes[i] = 0;
1358 else if (fFunctionCodes[i] >= 2)
1359 continue;
1360
1361
1362 if (fIsVerbose && studied == 1)
1363 edm::LogInfo("TMultiDimFet") << "Coeff SumSqRes Contrib Angle QM Func"
1364 << " Value W^2 Powers"
1365 << "\n";
1366
1367
1368 Double_t dResidur = MakeGramSchmidt(i);
1369
1370 if (dResidur == 0) {
1371
1372
1373 fFunctionCodes[i] = 1;
1374 continue;
1375 }
1376
1377
1378 if (!fIsUserFunction) {
1379
1380 fFunctionCodes[i] = 2;
1381
1382
1383 if (!TestFunction(squareResidual, dResidur)) {
1384 fFunctionCodes[i] = 1;
1385 continue;
1386 }
1387 }
1388
1389
1390
1391
1392
1393 fFunctionCodes[i] = 3;
1394 fPowerIndex[fNCoefficients] = i;
1395 fNCoefficients++;
1396
1397
1398
1399 squareResidual -= dResidur;
1400
1401
1402 for (j = 0; j < fNVariables; j++) {
1403 if (fNCoefficients == 1 || fMaxPowersFinal[j] <= fPowers[i * fNVariables + j] - 1)
1404 fMaxPowersFinal[j] = fPowers[i * fNVariables + j] - 1;
1405 }
1406 Double_t s = EvalControl(&fPowers[i * fNVariables]);
1407
1408
1409 if (fIsVerbose) {
1410 edm::LogVerbatim("TMultiDimFet") << std::setw(5) << fNCoefficients << " " << std::setw(10) << std::setprecision(4)
1411 << squareResidual << " " << std::setw(10) << std::setprecision(4) << dResidur
1412 << " " << std::setw(7) << std::setprecision(3) << fMaxAngle << " "
1413 << std::setw(7) << std::setprecision(3) << s << " " << std::setw(5) << i << " "
1414 << std::setw(10) << std::setprecision(4) << fOrthCoefficients(fNCoefficients - 1)
1415 << " " << std::setw(10) << std::setprecision(4)
1416 << fOrthFunctionNorms(fNCoefficients - 1) << " " << std::flush;
1417 for (j = 0; j < fNVariables; j++)
1418 edm::LogInfo("TMultiDimFet") << " " << fPowers[i * fNVariables + j] - 1 << std::flush;
1419 edm::LogInfo("TMultiDimFet") << "\n";
1420 }
1421
1422 if (fNCoefficients >= fMaxTerms ) {
1423 fParameterisationCode = PARAM_MAXTERMS;
1424 break;
1425 }
1426
1427 Double_t err = TMath::Sqrt(TMath::Max(1e-20, squareResidual) / fSumSqAvgQuantity);
1428 if (err < fMinRelativeError) {
1429 fParameterisationCode = PARAM_RELERR;
1430 break;
1431 }
1432 }
1433
1434 fError = TMath::Max(1e-20, squareResidual);
1435 fSumSqResidual -= fError;
1436 fRMS = TMath::Sqrt(fError / fSampleSize);
1437 }
1438
1439
1440 void TMultiDimFet::MakeRealCode(const char *filename, const char *classname, Option_t *) {
1441
1442
1443
1444
1445
1446
1447 Int_t i, j;
1448
1449 Bool_t isMethod = (classname[0] == '\0' ? kFALSE : kTRUE);
1450 const char *prefix = (isMethod ? Form("%s::", classname) : "");
1451 const char *cv_qual = (isMethod ? "" : "static ");
1452
1453 std::ofstream outFile(filename, std::ios::out | std::ios::trunc);
1454 if (!outFile) {
1455 Error("MakeRealCode", "couldn't open output file '%s'", filename);
1456 return;
1457 }
1458
1459 if (fIsVerbose)
1460 edm::LogInfo("TMultiDimFet") << "Writing on file \"" << filename << "\" ... " << std::flush;
1461
1462
1463
1464
1465 outFile << "// -*- mode: c++ -*-"
1466 << "\n";
1467
1468 outFile << "// "
1469 << "\n"
1470 << "// File " << filename << " generated by TMultiDimFet::MakeRealCode"
1471 << "\n";
1472
1473 TDatime date;
1474 outFile << "// on " << date.AsString() << "\n";
1475
1476 outFile << "// ROOT version " << gROOT->GetVersion() << "\n"
1477 << "//"
1478 << "\n";
1479
1480 outFile << "// This file contains the function "
1481 << "\n"
1482 << "//"
1483 << "\n"
1484 << "// double " << prefix << "MDF(double *x); "
1485 << "\n"
1486 << "//"
1487 << "\n"
1488 << "// For evaluating the parameterization obtained"
1489 << "\n"
1490 << "// from TMultiDimFet and the point x"
1491 << "\n"
1492 << "// "
1493 << "\n"
1494 << "// See TMultiDimFet class documentation for more "
1495 << "information "
1496 << "\n"
1497 << "// "
1498 << "\n";
1499
1500 if (isMethod)
1501
1502 outFile << "#include \"" << classname << ".h\""
1503 << "\n";
1504
1505
1506
1507
1508 outFile << "//"
1509 << "\n"
1510 << "// Static data variables"
1511 << "\n"
1512 << "//"
1513 << "\n";
1514 outFile << cv_qual << "int " << prefix << "gNVariables = " << fNVariables << ";"
1515 << "\n";
1516 outFile << cv_qual << "int " << prefix << "gNCoefficients = " << fNCoefficients << ";"
1517 << "\n";
1518 outFile << cv_qual << "double " << prefix << "gDMean = " << fMeanQuantity << ";"
1519 << "\n";
1520
1521
1522 outFile << "// Assignment to mean vector."
1523 << "\n";
1524 outFile << cv_qual << "double " << prefix << "gXMean[] = {"
1525 << "\n";
1526 for (i = 0; i < fNVariables; i++)
1527 outFile << (i != 0 ? ", " : " ") << fMeanVariables(i) << std::flush;
1528 outFile << " };"
1529 << "\n"
1530 << "\n";
1531
1532
1533 outFile << "// Assignment to minimum vector."
1534 << "\n";
1535 outFile << cv_qual << "double " << prefix << "gXMin[] = {"
1536 << "\n";
1537 for (i = 0; i < fNVariables; i++)
1538 outFile << (i != 0 ? ", " : " ") << fMinVariables(i) << std::flush;
1539 outFile << " };"
1540 << "\n"
1541 << "\n";
1542
1543
1544 outFile << "// Assignment to maximum vector."
1545 << "\n";
1546 outFile << cv_qual << "double " << prefix << "gXMax[] = {"
1547 << "\n";
1548 for (i = 0; i < fNVariables; i++)
1549 outFile << (i != 0 ? ", " : " ") << fMaxVariables(i) << std::flush;
1550 outFile << " };"
1551 << "\n"
1552 << "\n";
1553
1554
1555 outFile << "// Assignment to coefficients vector."
1556 << "\n";
1557 outFile << cv_qual << "double " << prefix << "gCoefficient[] = {" << std::flush;
1558 for (i = 0; i < fNCoefficients; i++)
1559 outFile << (i != 0 ? "," : "") << "\n"
1560 << " " << fCoefficients(i) << std::flush;
1561 outFile << "\n"
1562 << " };"
1563 << "\n"
1564 << "\n";
1565
1566
1567 outFile << "// Assignment to powers vector."
1568 << "\n"
1569 << "// The powers are stored row-wise, that is"
1570 << "\n"
1571 << "// p_ij = " << prefix << "gPower[i * NVariables + j];"
1572 << "\n";
1573 outFile << cv_qual << "int " << prefix << "gPower[] = {" << std::flush;
1574 for (i = 0; i < fNCoefficients; i++) {
1575 for (j = 0; j < fNVariables; j++) {
1576 if (j != 0)
1577 outFile << std::flush << " ";
1578 else
1579 outFile << "\n"
1580 << " ";
1581 outFile << fPowers[fPowerIndex[i] * fNVariables + j]
1582 << (i == fNCoefficients - 1 && j == fNVariables - 1 ? "" : ",") << std::flush;
1583 }
1584 }
1585 outFile << "\n"
1586 << "};"
1587 << "\n"
1588 << "\n";
1589
1590
1591
1592
1593 outFile << "// "
1594 << "\n"
1595 << "// The " << (isMethod ? "method " : "function ") << " double " << prefix << "MDF(double *x)"
1596 << "\n"
1597 << "// "
1598 << "\n";
1599 outFile << "double " << prefix << "MDF(double *x) {"
1600 << "\n"
1601 << " double returnValue = " << prefix << "gDMean;"
1602 << "\n"
1603 << " int i = 0, j = 0, k = 0;"
1604 << "\n"
1605 << " for (i = 0; i < " << prefix << "gNCoefficients ; i++) {"
1606 << "\n"
1607 << " // Evaluate the ith term in the expansion"
1608 << "\n"
1609 << " double term = " << prefix << "gCoefficient[i];"
1610 << "\n"
1611 << " for (j = 0; j < " << prefix << "gNVariables; j++) {"
1612 << "\n"
1613 << " // Evaluate the polynomial in the jth variable."
1614 << "\n"
1615 << " int power = " << prefix << "gPower[" << prefix << "gNVariables * i + j]; "
1616 << "\n"
1617 << " double p1 = 1, p2 = 0, p3 = 0, r = 0;"
1618 << "\n"
1619 << " double v = 1 + 2. / (" << prefix << "gXMax[j] - " << prefix << "gXMin[j]) * (x[j] - " << prefix
1620 << "gXMax[j]);"
1621 << "\n"
1622 << " // what is the power to use!"
1623 << "\n"
1624 << " switch(power) {"
1625 << "\n"
1626 << " case 1: r = 1; break; "
1627 << "\n"
1628 << " case 2: r = v; break; "
1629 << "\n"
1630 << " default: "
1631 << "\n"
1632 << " p2 = v; "
1633 << "\n"
1634 << " for (k = 3; k <= power; k++) { "
1635 << "\n"
1636 << " p3 = p2 * v;"
1637 << "\n";
1638 if (fPolyType == kLegendre)
1639 outFile << " p3 = ((2 * i - 3) * p2 * v - (i - 2) * p1)"
1640 << " / (i - 1);"
1641 << "\n";
1642 if (fPolyType == kChebyshev)
1643 outFile << " p3 = 2 * v * p2 - p1; "
1644 << "\n";
1645 outFile << " p1 = p2; p2 = p3; "
1646 << "\n"
1647 << " }"
1648 << "\n"
1649 << " r = p3;"
1650 << "\n"
1651 << " }"
1652 << "\n"
1653 << " // multiply this term by the poly in the jth var"
1654 << "\n"
1655 << " term *= r; "
1656 << "\n"
1657 << " }"
1658 << "\n"
1659 << " // Add this term to the final result"
1660 << "\n"
1661 << " returnValue += term;"
1662 << "\n"
1663 << " }"
1664 << "\n"
1665 << " return returnValue;"
1666 << "\n"
1667 << "}"
1668 << "\n"
1669 << "\n";
1670
1671
1672 outFile << "// EOF for " << filename << "\n";
1673
1674
1675 outFile.close();
1676
1677 if (fIsVerbose)
1678 edm::LogInfo("TMultiDimFet") << "done"
1679 << "\n";
1680 }
1681
1682
1683 void TMultiDimFet::Print(Option_t *option) const {
1684
1685
1686
1687
1688
1689
1690
1691
1692
1693
1694 Int_t i = 0;
1695 Int_t j = 0;
1696
1697 TString opt(option);
1698 opt.ToLower();
1699
1700 if (opt.Contains("p")) {
1701
1702 edm::LogInfo("TMultiDimFet") << "User parameters:"
1703 << "\n"
1704 << "----------------"
1705 << "\n"
1706 << " Variables: " << fNVariables << "\n"
1707 << " Data points: " << fSampleSize << "\n"
1708 << " Max Terms: " << fMaxTerms << "\n"
1709 << " Power Limit Parameter: " << fPowerLimit << "\n"
1710 << " Max functions: " << fMaxFunctions << "\n"
1711 << " Max functions to study: " << fMaxStudy << "\n"
1712 << " Max angle (optional): " << fMaxAngle << "\n"
1713 << " Min angle: " << fMinAngle << "\n"
1714 << " Relative Error accepted: " << fMinRelativeError << "\n"
1715 << " Maximum Powers: " << std::flush;
1716 for (i = 0; i < fNVariables; i++)
1717 edm::LogInfo("TMultiDimFet") << " " << fMaxPowers[i] - 1 << std::flush;
1718 edm::LogInfo("TMultiDimFet") << "\n"
1719 << "\n"
1720 << " Parameterisation will be done using " << std::flush;
1721 if (fPolyType == kChebyshev)
1722 edm::LogInfo("TMultiDimFet") << "Chebyshev polynomials"
1723 << "\n";
1724 else if (fPolyType == kLegendre)
1725 edm::LogInfo("TMultiDimFet") << "Legendre polynomials"
1726 << "\n";
1727 else
1728 edm::LogInfo("TMultiDimFet") << "Monomials"
1729 << "\n";
1730 edm::LogInfo("TMultiDimFet") << "\n";
1731 }
1732
1733 if (opt.Contains("s")) {
1734
1735 edm::LogInfo("TMultiDimFet") << "Sample statistics:"
1736 << "\n"
1737 << "------------------"
1738 << "\n"
1739 << " D" << std::flush;
1740 for (i = 0; i < fNVariables; i++)
1741 edm::LogInfo("TMultiDimFet") << " " << std::setw(10) << i + 1 << std::flush;
1742 edm::LogInfo("TMultiDimFet") << "\n"
1743 << " Max: " << std::setw(10) << std::setprecision(7) << fMaxQuantity << std::flush;
1744 for (i = 0; i < fNVariables; i++)
1745 edm::LogInfo("TMultiDimFet") << " " << std::setw(10) << std::setprecision(4) << fMaxVariables(i) << std::flush;
1746 edm::LogInfo("TMultiDimFet") << "\n"
1747 << " Min: " << std::setw(10) << std::setprecision(7) << fMinQuantity << std::flush;
1748 for (i = 0; i < fNVariables; i++)
1749 edm::LogInfo("TMultiDimFet") << " " << std::setw(10) << std::setprecision(4) << fMinVariables(i) << std::flush;
1750 edm::LogInfo("TMultiDimFet") << "\n"
1751 << " Mean: " << std::setw(10) << std::setprecision(7) << fMeanQuantity << std::flush;
1752 for (i = 0; i < fNVariables; i++)
1753 edm::LogInfo("TMultiDimFet") << " " << std::setw(10) << std::setprecision(4) << fMeanVariables(i) << std::flush;
1754 edm::LogInfo("TMultiDimFet") << "\n"
1755 << " Function Sum Squares: " << fSumSqQuantity << "\n"
1756 << "\n";
1757 }
1758
1759 if (opt.Contains("r")) {
1760 edm::LogInfo("TMultiDimFet") << "Results of Parameterisation:"
1761 << "\n"
1762 << "----------------------------"
1763 << "\n"
1764 << " Total reduction of square residuals " << fSumSqResidual << "\n"
1765 << " Relative precision obtained: " << fPrecision << "\n"
1766 << " Error obtained: " << fError << "\n"
1767 << " Multiple correlation coefficient: " << fCorrelationCoeff << "\n"
1768 << " Reduced Chi square over sample: " << fChi2 / (fSampleSize - fNCoefficients)
1769 << "\n"
1770 << " Maximum residual value: " << fMaxResidual << "\n"
1771 << " Minimum residual value: " << fMinResidual << "\n"
1772 << " Estimated root mean square: " << fRMS << "\n"
1773 << " Maximum powers used: " << std::flush;
1774 for (j = 0; j < fNVariables; j++)
1775 edm::LogInfo("TMultiDimFet") << fMaxPowersFinal[j] << " " << std::flush;
1776 edm::LogInfo("TMultiDimFet") << "\n"
1777 << " Function codes of candidate functions."
1778 << "\n"
1779 << " 1: considered,"
1780 << " 2: too little contribution,"
1781 << " 3: accepted." << std::flush;
1782 for (i = 0; i < fMaxFunctions; i++) {
1783 if (i % 60 == 0)
1784 edm::LogInfo("TMultiDimFet") << "\n"
1785 << " " << std::flush;
1786 else if (i % 10 == 0)
1787 edm::LogInfo("TMultiDimFet") << " " << std::flush;
1788 edm::LogInfo("TMultiDimFet") << fFunctionCodes[i];
1789 }
1790 edm::LogInfo("TMultiDimFet") << "\n"
1791 << " Loop over candidates stopped because " << std::flush;
1792 switch (fParameterisationCode) {
1793 case PARAM_MAXSTUDY:
1794 edm::LogInfo("TMultiDimFet") << "max allowed studies reached"
1795 << "\n";
1796 break;
1797 case PARAM_SEVERAL:
1798 edm::LogInfo("TMultiDimFet") << "all candidates considered several times"
1799 << "\n";
1800 break;
1801 case PARAM_RELERR:
1802 edm::LogInfo("TMultiDimFet") << "wanted relative error obtained"
1803 << "\n";
1804 break;
1805 case PARAM_MAXTERMS:
1806 edm::LogInfo("TMultiDimFet") << "max number of terms reached"
1807 << "\n";
1808 break;
1809 default:
1810 edm::LogInfo("TMultiDimFet") << "some unknown reason"
1811 << "\n";
1812 break;
1813 }
1814 edm::LogInfo("TMultiDimFet") << "\n";
1815 }
1816
1817 if (opt.Contains("f")) {
1818 edm::LogInfo("TMultiDimFet") << "Results of Fit:"
1819 << "\n"
1820 << "---------------"
1821 << "\n"
1822 << " Test sample size: " << fTestSampleSize << "\n"
1823 << " Multiple correlation coefficient: " << fTestCorrelationCoeff << "\n"
1824 << " Relative precision obtained: " << fTestPrecision << "\n"
1825 << " Error obtained: " << fTestError << "\n"
1826 << " Reduced Chi square over sample: " << fChi2 / (fSampleSize - fNCoefficients)
1827 << "\n"
1828 << "\n";
1829
1830
1831
1832
1833
1834
1835 }
1836
1837 if (opt.Contains("c")) {
1838 edm::LogInfo("TMultiDimFet") << "Coefficients:"
1839 << "\n"
1840 << "-------------"
1841 << "\n"
1842 << " # Value Error Powers"
1843 << "\n"
1844 << " ---------------------------------------"
1845 << "\n";
1846 for (i = 0; i < fNCoefficients; i++) {
1847 edm::LogInfo("TMultiDimFet") << " " << std::setw(3) << i << " " << std::setw(12) << fCoefficients(i) << " "
1848 << std::setw(12) << fCoefficientsRMS(i) << " " << std::flush;
1849 for (j = 0; j < fNVariables; j++)
1850 edm::LogInfo("TMultiDimFet") << " " << std::setw(3) << fPowers[fPowerIndex[i] * fNVariables + j] - 1
1851 << std::flush;
1852 edm::LogInfo("TMultiDimFet") << "\n";
1853 }
1854 edm::LogInfo("TMultiDimFet") << "\n";
1855 }
1856 if (opt.Contains("k") && fCorrelationMatrix.IsValid()) {
1857 edm::LogInfo("TMultiDimFet") << "Correlation Matrix:"
1858 << "\n"
1859 << "-------------------";
1860 fCorrelationMatrix.Print();
1861 }
1862
1863 if (opt.Contains("m")) {
1864 edm::LogInfo("TMultiDimFet") << std::setprecision(25);
1865 edm::LogInfo("TMultiDimFet") << "Parameterization:"
1866 << "\n"
1867 << "-----------------"
1868 << "\n"
1869 << " Normalised variables: "
1870 << "\n";
1871 for (i = 0; i < fNVariables; i++)
1872 edm::LogInfo("TMultiDimFet") << "\ty" << i << "\t:= 1 + 2 * (x" << i << " - " << fMaxVariables(i) << ") / ("
1873 << fMaxVariables(i) << " - " << fMinVariables(i) << ")"
1874 << "\n";
1875 edm::LogInfo("TMultiDimFet") << "\n"
1876 << " f[";
1877 for (i = 0; i < fNVariables; i++) {
1878 edm::LogInfo("TMultiDimFet") << "y" << i;
1879 if (i != fNVariables - 1)
1880 edm::LogInfo("TMultiDimFet") << ", ";
1881 }
1882 edm::LogInfo("TMultiDimFet") << "] := ";
1883 for (Int_t i = 0; i < fNCoefficients; i++) {
1884 if (i != 0)
1885 edm::LogInfo("TMultiDimFet") << " " << (fCoefficients(i) < 0 ? "- " : "+ ") << TMath::Abs(fCoefficients(i));
1886 else
1887 edm::LogInfo("TMultiDimFet") << fCoefficients(i);
1888 for (Int_t j = 0; j < fNVariables; j++) {
1889 Int_t p = fPowers[fPowerIndex[i] * fNVariables + j];
1890 switch (p) {
1891 case 1:
1892 break;
1893 case 2:
1894 edm::LogInfo("TMultiDimFet") << " * y" << j;
1895 break;
1896 default:
1897 switch (fPolyType) {
1898 case kLegendre:
1899 edm::LogInfo("TMultiDimFet") << " * L" << p - 1 << "(y" << j << ")";
1900 break;
1901 case kChebyshev:
1902 edm::LogInfo("TMultiDimFet") << " * C" << p - 1 << "(y" << j << ")";
1903 break;
1904 default:
1905 edm::LogInfo("TMultiDimFet") << " * y" << j << "^" << p - 1;
1906 break;
1907 }
1908 }
1909 }
1910 }
1911 edm::LogInfo("TMultiDimFet") << "\n";
1912 }
1913 }
1914
1915
1916 void TMultiDimFet::PrintPolynomialsSpecial(Option_t *option) const {
1917
1918
1919 Int_t i = 0;
1920
1921
1922 TString opt(option);
1923 opt.ToLower();
1924
1925 if (opt.Contains("m")) {
1926 edm::LogInfo("TMultiDimFet") << std::setprecision(25);
1927 edm::LogInfo("TMultiDimFet") << "Parameterization:"
1928 << "\n"
1929 << "-----------------"
1930 << "\n"
1931 << " Normalised variables: "
1932 << "\n";
1933 for (i = 0; i < fNVariables; i++)
1934 edm::LogInfo("TMultiDimFet") << "\tdouble y" << i << "\t=1+2*(x" << i << "-" << fMaxVariables(i) << ")/("
1935 << fMaxVariables(i) << "-" << fMinVariables(i) << ");"
1936 << "\n";
1937 edm::LogInfo("TMultiDimFet") << "\n"
1938 << " f[";
1939 for (i = 0; i < fNVariables; i++) {
1940 edm::LogInfo("TMultiDimFet") << "y" << i;
1941 if (i != fNVariables - 1)
1942 edm::LogInfo("TMultiDimFet") << ", ";
1943 }
1944 edm::LogInfo("TMultiDimFet") << "] := " << fMeanQuantity << " + ";
1945 for (Int_t i = 0; i < fNCoefficients; i++) {
1946 if (i != 0)
1947 edm::LogInfo("TMultiDimFet") << " " << (fCoefficients(i) < 0 ? "-" : "+") << TMath::Abs(fCoefficients(i));
1948 else
1949 edm::LogInfo("TMultiDimFet") << fCoefficients(i);
1950 for (Int_t j = 0; j < fNVariables; j++) {
1951 Int_t p = fPowers[fPowerIndex[i] * fNVariables + j];
1952 switch (p) {
1953 case 1:
1954 break;
1955 case 2:
1956 edm::LogInfo("TMultiDimFet") << "*y" << j;
1957 break;
1958 default:
1959 switch (fPolyType) {
1960 case kLegendre:
1961 edm::LogInfo("TMultiDimFet") << "*Leg(" << p - 1 << ",y" << j << ")";
1962 break;
1963 case kChebyshev:
1964 edm::LogInfo("TMultiDimFet") << "*C" << p - 1 << "(y" << j << ")";
1965 break;
1966 default:
1967 edm::LogInfo("TMultiDimFet") << "*y" << j << "**" << p - 1;
1968 break;
1969 }
1970 }
1971 }
1972 edm::LogInfo("TMultiDimFet") << "\n";
1973 }
1974 edm::LogInfo("TMultiDimFet") << "\n";
1975 }
1976 }
1977
1978
1979 Bool_t TMultiDimFet::Select(const Int_t *) {
1980
1981
1982
1983
1984
1985
1986
1987
1988
1989 return kTRUE;
1990 }
1991
1992
1993 void TMultiDimFet::SetMaxAngle(Double_t ang) {
1994
1995
1996
1997
1998
1999 if (ang >= 90 || ang < 0) {
2000 Warning("SetMaxAngle", "angle must be in [0,90)");
2001 return;
2002 }
2003
2004 fMaxAngle = ang;
2005 }
2006
2007
2008 void TMultiDimFet::SetMinAngle(Double_t ang) {
2009
2010
2011
2012
2013 if (ang > 90 || ang <= 0) {
2014 Warning("SetMinAngle", "angle must be in [0,90)");
2015 return;
2016 }
2017
2018 fMinAngle = ang;
2019 }
2020
2021
2022 void TMultiDimFet::SetPowers(const Int_t *powers, Int_t terms) {
2023
2024
2025
2026
2027
2028
2029 fIsUserFunction = kTRUE;
2030 fMaxFunctions = terms;
2031 fMaxTerms = terms;
2032 fMaxStudy = terms;
2033 fMaxFunctionsTimesNVariables = fMaxFunctions * fNVariables;
2034 fPowers.resize(fMaxFunctions * fNVariables);
2035 Int_t i, j;
2036 for (i = 0; i < fMaxFunctions; i++)
2037 for (j = 0; j < fNVariables; j++)
2038 fPowers[i * fNVariables + j] = powers[i * fNVariables + j] + 1;
2039 }
2040
2041
2042 void TMultiDimFet::SetPowerLimit(Double_t limit) {
2043
2044
2045
2046
2047 fPowerLimit = limit;
2048 }
2049
2050
2051 void TMultiDimFet::SetMaxPowers(const Int_t *powers) {
2052
2053
2054
2055 if (!powers)
2056 return;
2057
2058 for (Int_t i = 0; i < fNVariables; i++)
2059 fMaxPowers[i] = powers[i] + 1;
2060 }
2061
2062
2063 void TMultiDimFet::SetMinRelativeError(Double_t error) {
2064
2065
2066
2067
2068 fMinRelativeError = error;
2069 }
2070
2071
2072 Bool_t TMultiDimFet::TestFunction(Double_t squareResidual, Double_t dResidur) {
2073
2074
2075
2076
2077
2078 if (fNCoefficients != 0) {
2079
2080 if (fMaxAngle == 0) {
2081
2082 if (dResidur < squareResidual / (fMaxTerms - fNCoefficients + 1 + 1E-10)) {
2083 return kFALSE;
2084 }
2085 } else {
2086
2087
2088 if (TMath::Sqrt(dResidur / fSumSqAvgQuantity) < TMath::Cos(fMaxAngle * DEGRAD)) {
2089 return kFALSE;
2090 }
2091 }
2092 }
2093
2094 return kTRUE;
2095 }
2096
2097
2098
2099
2100
2101
2102
2103
2104
2105
2106
2107