File indexing completed on 2024-10-04 22:55:07
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 *this;
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 *this;
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 assert(fNVariables > 0);
0670 Int_t *powers = new Int_t[fNVariables * fMaxFunctions];
0671
0672
0673 Double_t *control = new Double_t[fMaxFunctions];
0674
0675
0676 Int_t *iv = new Int_t[fNVariables];
0677 for (i = 0; i < fNVariables; i++)
0678 iv[i] = 1;
0679
0680 if (!fIsUserFunction) {
0681
0682 Int_t numberFunctions = 0;
0683
0684 while (kTRUE) {
0685
0686 Double_t s = EvalControl(iv);
0687
0688 if (s <= fPowerLimit) {
0689
0690
0691 if (Select(iv)) {
0692 numberFunctions++;
0693
0694
0695
0696 if (numberFunctions > fMaxFunctions)
0697 break;
0698
0699
0700
0701 control[numberFunctions - 1] = Int_t(1.0e+6 * s);
0702
0703
0704 for (i = 0; i < fNVariables; i++) {
0705 j = (numberFunctions - 1) * fNVariables + i;
0706 powers[j] = iv[i];
0707 }
0708 }
0709 }
0710
0711 for (i = 0; i < fNVariables; i++)
0712 if (iv[i] < fMaxPowers[i])
0713 break;
0714
0715
0716
0717 if (i == fNVariables) {
0718 fMaxFunctions = numberFunctions;
0719 fMaxFunctionsTimesNVariables = fMaxFunctions * fNVariables;
0720 break;
0721 }
0722
0723
0724 iv[i]++;
0725
0726 for (j = 0; j < i; j++)
0727 iv[j] = 1;
0728 }
0729 } else {
0730
0731 for (i = 0; i < fMaxFunctions; i++) {
0732
0733 for (j = 0; j < fNVariables; j++) {
0734 powers[i * fNVariables + j] = fPowers[i * fNVariables + j];
0735 iv[j] = fPowers[i * fNVariables + j];
0736 }
0737
0738 control[i] = Int_t(1.0e+6 * EvalControl(iv));
0739 }
0740 }
0741
0742
0743
0744 Int_t *order = new Int_t[fMaxFunctions];
0745 for (i = 0; i < fMaxFunctions; i++)
0746 order[i] = i;
0747 fPowers.resize(fMaxFunctions * fNVariables);
0748
0749 for (i = 0; i < fMaxFunctions; i++) {
0750 Double_t x = control[i];
0751 Int_t l = order[i];
0752 k = i;
0753
0754 for (j = i; j < fMaxFunctions; j++) {
0755 if (control[j] <= x) {
0756 x = control[j];
0757 l = order[j];
0758 k = j;
0759 }
0760 }
0761
0762 if (k != i) {
0763 control[k] = control[i];
0764 control[i] = x;
0765 order[k] = order[i];
0766 order[i] = l;
0767 }
0768 }
0769
0770 for (i = 0; i < fMaxFunctions; i++)
0771 for (j = 0; j < fNVariables; j++)
0772 fPowers[i * fNVariables + j] = powers[order[i] * fNVariables + j];
0773
0774 delete[] control;
0775 delete[] powers;
0776 delete[] order;
0777 delete[] iv;
0778 }
0779
0780 Double_t TMultiDimFet::MakeChi2(const Double_t *coeff) {
0781
0782
0783
0784
0785
0786 fChi2 = 0;
0787 Int_t i, j;
0788 Double_t *x = new Double_t[fNVariables];
0789 for (i = 0; i < fTestSampleSize; i++) {
0790
0791 for (j = 0; j < fNVariables; j++)
0792 x[j] = fTestVariables(i * fNVariables + j);
0793
0794
0795 Double_t f = Eval(x, coeff);
0796
0797
0798 fChi2 += 1. / TMath::Max(fTestSqError(i), 1e-20) * (fTestQuantity(i) - f) * (fTestQuantity(i) - f);
0799 }
0800
0801
0802 delete[] x;
0803
0804 return fChi2;
0805 }
0806
0807
0808 void TMultiDimFet::MakeCode(const char *filename, Option_t *option) {
0809
0810
0811
0812
0813
0814
0815
0816
0817
0818
0819
0820
0821
0822
0823
0824
0825
0826
0827
0828
0829
0830
0831
0832
0833 TString outName(filename);
0834 if (!outName.EndsWith(".C") && !outName.EndsWith(".cxx"))
0835 outName += ".C";
0836
0837 MakeRealCode(outName.Data(), "", option);
0838 }
0839
0840 void TMultiDimFet::MakeCoefficientErrors() {
0841
0842
0843
0844 Int_t i = 0;
0845 Int_t j = 0;
0846 Int_t k = 0;
0847 TVectorD iF(fSampleSize);
0848 TVectorD jF(fSampleSize);
0849 fCoefficientsRMS.ResizeTo(fNCoefficients);
0850
0851 TMatrixDSym curvatureMatrix(fNCoefficients);
0852
0853
0854 for (i = 0; i < fNCoefficients; i++) {
0855 iF = TMatrixDRow(fFunctions, i);
0856 for (j = 0; j <= i; j++) {
0857 jF = TMatrixDRow(fFunctions, j);
0858 for (k = 0; k < fSampleSize; k++)
0859 curvatureMatrix(i, j) += 1 / TMath::Max(fSqError(k), 1e-20) * iF(k) * jF(k);
0860 curvatureMatrix(j, i) = curvatureMatrix(i, j);
0861 }
0862 }
0863
0864
0865 fChi2 = 0;
0866 for (i = 0; i < fSampleSize; i++) {
0867 Double_t f = 0;
0868 for (j = 0; j < fNCoefficients; j++)
0869 f += fCoefficients(j) * fFunctions(j, i);
0870 fChi2 += 1. / TMath::Max(fSqError(i), 1e-20) * (fQuantity(i) - f) * (fQuantity(i) - f);
0871 }
0872
0873
0874 const TVectorD diag = TMatrixDDiag_const(curvatureMatrix);
0875 curvatureMatrix.NormByDiag(diag);
0876
0877 TDecompChol chol(curvatureMatrix);
0878 if (!chol.Decompose())
0879 Error("MakeCoefficientErrors", "curvature matrix is singular");
0880 chol.Invert(curvatureMatrix);
0881
0882 curvatureMatrix.NormByDiag(diag);
0883
0884 for (i = 0; i < fNCoefficients; i++)
0885 fCoefficientsRMS(i) = TMath::Sqrt(curvatureMatrix(i, i));
0886 }
0887
0888
0889 void TMultiDimFet::MakeCoefficients() {
0890
0891
0892
0893
0894
0895
0896
0897
0898 Int_t i = 0, j = 0;
0899 Int_t col = 0, row = 0;
0900
0901
0902 for (col = 1; col < fNCoefficients; col++) {
0903 for (row = col - 1; row > -1; row--) {
0904 fOrthCurvatureMatrix(row, col) = 0;
0905 for (i = row; i <= col; i++)
0906 fOrthCurvatureMatrix(row, col) -= fOrthCurvatureMatrix(i, row) * fOrthCurvatureMatrix(i, col);
0907 }
0908 }
0909
0910
0911 fCoefficients.ResizeTo(fNCoefficients);
0912
0913 for (i = 0; i < fNCoefficients; i++) {
0914 Double_t sum = 0;
0915 for (j = i; j < fNCoefficients; j++)
0916 sum += fOrthCurvatureMatrix(i, j) * fOrthCoefficients(j);
0917 fCoefficients(i) = sum;
0918 }
0919
0920
0921 fResiduals.ResizeTo(fSampleSize);
0922 for (i = 0; i < fSampleSize; i++)
0923 fResiduals(i) = fQuantity(i);
0924
0925 for (i = 0; i < fNCoefficients; i++)
0926 for (j = 0; j < fSampleSize; j++)
0927 fResiduals(j) -= fCoefficients(i) * fFunctions(i, j);
0928
0929
0930
0931 fMinResidual = 10e10;
0932 fMaxResidual = -10e10;
0933 Double_t sqRes = 0;
0934 for (i = 0; i < fSampleSize; i++) {
0935 sqRes += fResiduals(i) * fResiduals(i);
0936 if (fResiduals(i) <= fMinResidual) {
0937 fMinResidual = fResiduals(i);
0938 fMinResidualRow = i;
0939 }
0940 if (fResiduals(i) >= fMaxResidual) {
0941 fMaxResidual = fResiduals(i);
0942 fMaxResidualRow = i;
0943 }
0944 }
0945
0946 fCorrelationCoeff = fSumSqResidual / fSumSqAvgQuantity;
0947 fPrecision = TMath::Sqrt(sqRes / fSumSqQuantity);
0948
0949
0950 if (TESTBIT(fHistogramMask, HIST_RD) || TESTBIT(fHistogramMask, HIST_RTRAI) || TESTBIT(fHistogramMask, HIST_RX)) {
0951 for (i = 0; i < fSampleSize; i++) {
0952 if (TESTBIT(fHistogramMask, HIST_RD))
0953 ((TH2D *)fHistograms->FindObject("res_d"))->Fill(fQuantity(i), fResiduals(i));
0954 if (TESTBIT(fHistogramMask, HIST_RTRAI))
0955 ((TH1D *)fHistograms->FindObject("res_train"))->Fill(fResiduals(i));
0956
0957 if (TESTBIT(fHistogramMask, HIST_RX))
0958 for (j = 0; j < fNVariables; j++)
0959 ((TH2D *)fHistograms->FindObject(Form("res_x_%d", j)))->Fill(fVariables(i * fNVariables + j), fResiduals(i));
0960 }
0961 }
0962 }
0963
0964
0965 void TMultiDimFet::MakeCorrelation() {
0966
0967
0968 if (!fShowCorrelation)
0969 return;
0970
0971 fCorrelationMatrix.ResizeTo(fNVariables, fNVariables + 1);
0972
0973 Double_t d2 = 0;
0974 Double_t ddotXi = 0;
0975 Double_t xiNorm = 0;
0976 Double_t xidotXj = 0;
0977 Double_t xjNorm = 0;
0978
0979 Int_t i, j, k, l, m;
0980 for (i = 0; i < fSampleSize; i++)
0981 d2 += fQuantity(i) * fQuantity(i);
0982
0983 for (i = 0; i < fNVariables; i++) {
0984 ddotXi = 0.;
0985 xiNorm = 0.;
0986 for (j = 0; j < fSampleSize; j++) {
0987
0988 k = j * fNVariables + i;
0989 ddotXi += fQuantity(j) * (fVariables(k) - fMeanVariables(i));
0990 xiNorm += (fVariables(k) - fMeanVariables(i)) * (fVariables(k) - fMeanVariables(i));
0991 }
0992 fCorrelationMatrix(i, 0) = ddotXi / TMath::Sqrt(d2 * xiNorm);
0993
0994 for (j = 0; j < i; j++) {
0995 xidotXj = 0.;
0996 xjNorm = 0.;
0997 for (k = 0; k < fSampleSize; k++) {
0998
0999
1000 l = k * fNVariables + j;
1001 m = k * fNVariables + i;
1002
1003
1004 xidotXj +=
1005 (fVariables(m) - fMeanVariables(i)) * (fVariables(l) - fMeanVariables(j));
1006 xjNorm += (fVariables(l) - fMeanVariables(j)) * (fVariables(l) - fMeanVariables(j));
1007 }
1008 fCorrelationMatrix(i, j + 1) = xidotXj / TMath::Sqrt(xiNorm * xjNorm);
1009 }
1010 }
1011 }
1012
1013
1014 Double_t TMultiDimFet::MakeGramSchmidt(Int_t function) {
1015
1016
1017
1018
1019
1020
1021
1022
1023 Double_t f2 = 0;
1024 fOrthCoefficients(fNCoefficients) = 0;
1025 fOrthFunctionNorms(fNCoefficients) = 0;
1026 Int_t j = 0;
1027 Int_t k = 0;
1028
1029 for (j = 0; j < fSampleSize; j++) {
1030 fFunctions(fNCoefficients, j) = 1;
1031 fOrthFunctions(fNCoefficients, j) = 0;
1032
1033 for (k = 0; k < fNVariables; k++) {
1034 Int_t p = fPowers[function * fNVariables + k];
1035 Double_t x = fVariables(j * fNVariables + k);
1036 fFunctions(fNCoefficients, j) *= EvalFactor(p, x);
1037 }
1038
1039
1040 f2 += fFunctions(fNCoefficients, j) * fFunctions(fNCoefficients, j);
1041
1042 fOrthFunctions(fNCoefficients, j) = fFunctions(fNCoefficients, j);
1043 }
1044
1045
1046 for (j = 0; j < fNCoefficients; j++) {
1047 Double_t fdw = 0;
1048
1049 for (k = 0; k < fSampleSize; k++) {
1050 fdw += fFunctions(fNCoefficients, k) * fOrthFunctions(j, k) / fOrthFunctionNorms(j);
1051 }
1052
1053 fOrthCurvatureMatrix(fNCoefficients, j) = fdw;
1054
1055 for (k = 0; k < fSampleSize; k++)
1056 fOrthFunctions(fNCoefficients, k) -= fdw * fOrthFunctions(j, k);
1057 }
1058
1059 for (j = 0; j < fSampleSize; j++) {
1060
1061 fOrthFunctionNorms(fNCoefficients) += fOrthFunctions(fNCoefficients, j) * fOrthFunctions(fNCoefficients, j);
1062
1063
1064 fOrthCoefficients(fNCoefficients) += fQuantity(j) * fOrthFunctions(fNCoefficients, j);
1065 }
1066
1067
1068 if (!fIsUserFunction)
1069 if (TMath::Sqrt(fOrthFunctionNorms(fNCoefficients) / (f2 + 1e-10)) < TMath::Sin(fMinAngle * DEGRAD))
1070 return 0;
1071
1072
1073
1074
1075
1076 fOrthCurvatureMatrix(fNCoefficients, fNCoefficients) = 1;
1077 Double_t b = fOrthCoefficients(fNCoefficients);
1078 fOrthCoefficients(fNCoefficients) /= fOrthFunctionNorms(fNCoefficients);
1079
1080
1081 Double_t dResidur = fOrthCoefficients(fNCoefficients) * b;
1082
1083 return dResidur;
1084 }
1085
1086
1087 void TMultiDimFet::MakeHistograms(Option_t *option) {
1088
1089
1090
1091
1092
1093
1094
1095
1096
1097
1098
1099
1100
1101
1102
1103
1104
1105 TString opt(option);
1106 opt.ToLower();
1107
1108 if (opt.Length() < 1)
1109 return;
1110
1111 if (!fHistograms)
1112 fHistograms = new TList;
1113
1114
1115 Int_t i = 0;
1116
1117
1118 if (opt.Contains("x") || opt.Contains("a")) {
1119 SETBIT(fHistogramMask, HIST_XORIG);
1120 for (i = 0; i < fNVariables; i++)
1121 if (!fHistograms->FindObject(Form("x_%d_orig", i)))
1122 fHistograms->Add(
1123 new TH1D(Form("x_%d_orig", i), Form("Original variable # %d", i), 100, fMinVariables(i), fMaxVariables(i)));
1124 }
1125
1126
1127 if (opt.Contains("d") || opt.Contains("a")) {
1128 SETBIT(fHistogramMask, HIST_DORIG);
1129 if (!fHistograms->FindObject("d_orig"))
1130 fHistograms->Add(new TH1D("d_orig", "Original Quantity", 100, fMinQuantity, fMaxQuantity));
1131 }
1132
1133
1134 if (opt.Contains("n") || opt.Contains("a")) {
1135 SETBIT(fHistogramMask, HIST_XNORM);
1136 for (i = 0; i < fNVariables; i++)
1137 if (!fHistograms->FindObject(Form("x_%d_norm", i)))
1138 fHistograms->Add(new TH1D(Form("x_%d_norm", i), Form("Normalized variable # %d", i), 100, -1, 1));
1139 }
1140
1141
1142 if (opt.Contains("s") || opt.Contains("a")) {
1143 SETBIT(fHistogramMask, HIST_DSHIF);
1144 if (!fHistograms->FindObject("d_shifted"))
1145 fHistograms->Add(
1146 new TH1D("d_shifted", "Shifted Quantity", 100, fMinQuantity - fMeanQuantity, fMaxQuantity - fMeanQuantity));
1147 }
1148
1149
1150 if (opt.Contains("r1") || opt.Contains("a")) {
1151 SETBIT(fHistogramMask, HIST_RX);
1152 for (i = 0; i < fNVariables; i++)
1153 if (!fHistograms->FindObject(Form("res_x_%d", i)))
1154 fHistograms->Add(new TH2D(Form("res_x_%d", i),
1155 Form("Computed residual versus x_%d", i),
1156 100,
1157 -1,
1158 1,
1159 35,
1160 fMinQuantity - fMeanQuantity,
1161 fMaxQuantity - fMeanQuantity));
1162 }
1163
1164
1165 if (opt.Contains("r2") || opt.Contains("a")) {
1166 SETBIT(fHistogramMask, HIST_RD);
1167 if (!fHistograms->FindObject("res_d"))
1168 fHistograms->Add(new TH2D("res_d",
1169 "Computed residuals vs Quantity",
1170 100,
1171 fMinQuantity - fMeanQuantity,
1172 fMaxQuantity - fMeanQuantity,
1173 35,
1174 fMinQuantity - fMeanQuantity,
1175 fMaxQuantity - fMeanQuantity));
1176 }
1177
1178
1179 if (opt.Contains("r3") || opt.Contains("a")) {
1180 SETBIT(fHistogramMask, HIST_RTRAI);
1181 if (!fHistograms->FindObject("res_train"))
1182 fHistograms->Add(new TH1D("res_train",
1183 "Computed residuals over training sample",
1184 100,
1185 fMinQuantity - fMeanQuantity,
1186 fMaxQuantity - fMeanQuantity));
1187 }
1188 if (opt.Contains("r4") || opt.Contains("a")) {
1189 SETBIT(fHistogramMask, HIST_RTEST);
1190 if (!fHistograms->FindObject("res_test"))
1191 fHistograms->Add(new TH1D("res_test",
1192 "Distribution of residuals from test",
1193 100,
1194 fMinQuantity - fMeanQuantity,
1195 fMaxQuantity - fMeanQuantity));
1196 }
1197 }
1198
1199
1200 void TMultiDimFet::MakeMethod(const Char_t *classname, Option_t *option) {
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
1245 MakeRealCode(Form("%sMDF.cxx", classname), classname, option);
1246 }
1247
1248
1249 void TMultiDimFet::MakeNormalized() {
1250
1251
1252
1253
1254 Int_t i = 0;
1255 Int_t j = 0;
1256 Int_t k = 0;
1257
1258 for (i = 0; i < fSampleSize; i++) {
1259 if (TESTBIT(fHistogramMask, HIST_DORIG))
1260 ((TH1D *)fHistograms->FindObject("d_orig"))->Fill(fQuantity(i));
1261
1262 fQuantity(i) -= fMeanQuantity;
1263 fSumSqAvgQuantity += fQuantity(i) * fQuantity(i);
1264
1265 if (TESTBIT(fHistogramMask, HIST_DSHIF))
1266 ((TH1D *)fHistograms->FindObject("d_shifted"))->Fill(fQuantity(i));
1267
1268 for (j = 0; j < fNVariables; j++) {
1269 Double_t range = 1. / (fMaxVariables(j) - fMinVariables(j));
1270 k = i * fNVariables + j;
1271
1272
1273 if (TESTBIT(fHistogramMask, HIST_XORIG))
1274 ((TH1D *)fHistograms->FindObject(Form("x_%d_orig", j)))->Fill(fVariables(k));
1275
1276
1277 fVariables(k) = 1 + 2 * range * (fVariables(k) - fMaxVariables(j));
1278
1279
1280 if (TESTBIT(fHistogramMask, HIST_XNORM))
1281 ((TH1D *)fHistograms->FindObject(Form("x_%d_norm", j)))->Fill(fVariables(k));
1282 }
1283 }
1284
1285 fMaxQuantity -= fMeanQuantity;
1286 fMinQuantity -= fMeanQuantity;
1287
1288
1289 for (i = 0; i < fNVariables; i++) {
1290 Double_t range = 1. / (fMaxVariables(i) - fMinVariables(i));
1291 fMeanVariables(i) = 1 + 2 * range * (fMeanVariables(i) - fMaxVariables(i));
1292 }
1293 }
1294
1295
1296 void TMultiDimFet::MakeParameterization() {
1297
1298
1299
1300
1301
1302 Int_t i = -1;
1303 Int_t j = 0;
1304 Int_t k = 0;
1305 Int_t maxPass = 3;
1306 Int_t studied = 0;
1307 Double_t squareResidual = fSumSqAvgQuantity;
1308 fNCoefficients = 0;
1309 fSumSqResidual = fSumSqAvgQuantity;
1310 fFunctions.ResizeTo(fMaxTerms, fSampleSize);
1311 fOrthFunctions.ResizeTo(fMaxTerms, fSampleSize);
1312 fOrthFunctionNorms.ResizeTo(fMaxTerms);
1313 fOrthCoefficients.ResizeTo(fMaxTerms);
1314 fOrthCurvatureMatrix.ResizeTo(fMaxTerms, fMaxTerms);
1315 fFunctions = 1;
1316
1317 fFunctionCodes.resize(fMaxFunctions);
1318 fPowerIndex.resize(fMaxTerms);
1319 Int_t l;
1320 for (l = 0; l < fMaxFunctions; l++)
1321 fFunctionCodes[l] = 0;
1322 for (l = 0; l < fMaxTerms; l++)
1323 fPowerIndex[l] = 0;
1324
1325 if (fMaxAngle != 0)
1326 maxPass = 100;
1327 if (fIsUserFunction)
1328 maxPass = 1;
1329
1330
1331
1332 while (kTRUE) {
1333
1334 if (studied++ >= fMaxStudy) {
1335 fParameterisationCode = PARAM_MAXSTUDY;
1336 break;
1337 }
1338
1339
1340 if (k >= maxPass) {
1341 fParameterisationCode = PARAM_SEVERAL;
1342 break;
1343 }
1344
1345
1346 i++;
1347
1348
1349 if (i == fMaxFunctions) {
1350 if (fMaxAngle != 0)
1351 fMaxAngle += (90 - fMaxAngle) / 2;
1352 i = 0;
1353 studied--;
1354 k++;
1355 continue;
1356 }
1357 if (studied == 1)
1358 fFunctionCodes[i] = 0;
1359 else if (fFunctionCodes[i] >= 2)
1360 continue;
1361
1362
1363 if (fIsVerbose && studied == 1)
1364 edm::LogInfo("TMultiDimFet") << "Coeff SumSqRes Contrib Angle QM Func"
1365 << " Value W^2 Powers"
1366 << "\n";
1367
1368
1369 Double_t dResidur = MakeGramSchmidt(i);
1370
1371 if (dResidur == 0) {
1372
1373
1374 fFunctionCodes[i] = 1;
1375 continue;
1376 }
1377
1378
1379 if (!fIsUserFunction) {
1380
1381 fFunctionCodes[i] = 2;
1382
1383
1384 if (!TestFunction(squareResidual, dResidur)) {
1385 fFunctionCodes[i] = 1;
1386 continue;
1387 }
1388 }
1389
1390
1391
1392
1393
1394 fFunctionCodes[i] = 3;
1395 fPowerIndex[fNCoefficients] = i;
1396 fNCoefficients++;
1397
1398
1399
1400 squareResidual -= dResidur;
1401
1402
1403 for (j = 0; j < fNVariables; j++) {
1404 if (fNCoefficients == 1 || fMaxPowersFinal[j] <= fPowers[i * fNVariables + j] - 1)
1405 fMaxPowersFinal[j] = fPowers[i * fNVariables + j] - 1;
1406 }
1407 Double_t s = EvalControl(&fPowers[i * fNVariables]);
1408
1409
1410 if (fIsVerbose) {
1411 edm::LogVerbatim("TMultiDimFet") << std::setw(5) << fNCoefficients << " " << std::setw(10) << std::setprecision(4)
1412 << squareResidual << " " << std::setw(10) << std::setprecision(4) << dResidur
1413 << " " << std::setw(7) << std::setprecision(3) << fMaxAngle << " "
1414 << std::setw(7) << std::setprecision(3) << s << " " << std::setw(5) << i << " "
1415 << std::setw(10) << std::setprecision(4) << fOrthCoefficients(fNCoefficients - 1)
1416 << " " << std::setw(10) << std::setprecision(4)
1417 << fOrthFunctionNorms(fNCoefficients - 1) << " " << std::flush;
1418 for (j = 0; j < fNVariables; j++)
1419 edm::LogInfo("TMultiDimFet") << " " << fPowers[i * fNVariables + j] - 1 << std::flush;
1420 edm::LogInfo("TMultiDimFet") << "\n";
1421 }
1422
1423 if (fNCoefficients >= fMaxTerms ) {
1424 fParameterisationCode = PARAM_MAXTERMS;
1425 break;
1426 }
1427
1428 Double_t err = TMath::Sqrt(TMath::Max(1e-20, squareResidual) / fSumSqAvgQuantity);
1429 if (err < fMinRelativeError) {
1430 fParameterisationCode = PARAM_RELERR;
1431 break;
1432 }
1433 }
1434
1435 fError = TMath::Max(1e-20, squareResidual);
1436 fSumSqResidual -= fError;
1437 fRMS = TMath::Sqrt(fError / fSampleSize);
1438 }
1439
1440
1441 void TMultiDimFet::MakeRealCode(const char *filename, const char *classname, Option_t *) {
1442
1443
1444
1445
1446
1447
1448 Int_t i, j;
1449
1450 Bool_t isMethod = (classname[0] == '\0' ? kFALSE : kTRUE);
1451 const char *prefix = (isMethod ? Form("%s::", classname) : "");
1452 const char *cv_qual = (isMethod ? "" : "static ");
1453
1454 std::ofstream outFile(filename, std::ios::out | std::ios::trunc);
1455 if (!outFile) {
1456 Error("MakeRealCode", "couldn't open output file '%s'", filename);
1457 return;
1458 }
1459
1460 if (fIsVerbose)
1461 edm::LogInfo("TMultiDimFet") << "Writing on file \"" << filename << "\" ... " << std::flush;
1462
1463
1464
1465
1466 outFile << "// -*- mode: c++ -*-"
1467 << "\n";
1468
1469 outFile << "// "
1470 << "\n"
1471 << "// File " << filename << " generated by TMultiDimFet::MakeRealCode"
1472 << "\n";
1473
1474 TDatime date;
1475 outFile << "// on " << date.AsString() << "\n";
1476
1477 outFile << "// ROOT version " << gROOT->GetVersion() << "\n"
1478 << "//"
1479 << "\n";
1480
1481 outFile << "// This file contains the function "
1482 << "\n"
1483 << "//"
1484 << "\n"
1485 << "// double " << prefix << "MDF(double *x); "
1486 << "\n"
1487 << "//"
1488 << "\n"
1489 << "// For evaluating the parameterization obtained"
1490 << "\n"
1491 << "// from TMultiDimFet and the point x"
1492 << "\n"
1493 << "// "
1494 << "\n"
1495 << "// See TMultiDimFet class documentation for more "
1496 << "information "
1497 << "\n"
1498 << "// "
1499 << "\n";
1500
1501 if (isMethod)
1502
1503 outFile << "#include \"" << classname << ".h\""
1504 << "\n";
1505
1506
1507
1508
1509 outFile << "//"
1510 << "\n"
1511 << "// Static data variables"
1512 << "\n"
1513 << "//"
1514 << "\n";
1515 outFile << cv_qual << "int " << prefix << "gNVariables = " << fNVariables << ";"
1516 << "\n";
1517 outFile << cv_qual << "int " << prefix << "gNCoefficients = " << fNCoefficients << ";"
1518 << "\n";
1519 outFile << cv_qual << "double " << prefix << "gDMean = " << fMeanQuantity << ";"
1520 << "\n";
1521
1522
1523 outFile << "// Assignment to mean vector."
1524 << "\n";
1525 outFile << cv_qual << "double " << prefix << "gXMean[] = {"
1526 << "\n";
1527 for (i = 0; i < fNVariables; i++)
1528 outFile << (i != 0 ? ", " : " ") << fMeanVariables(i) << std::flush;
1529 outFile << " };"
1530 << "\n"
1531 << "\n";
1532
1533
1534 outFile << "// Assignment to minimum vector."
1535 << "\n";
1536 outFile << cv_qual << "double " << prefix << "gXMin[] = {"
1537 << "\n";
1538 for (i = 0; i < fNVariables; i++)
1539 outFile << (i != 0 ? ", " : " ") << fMinVariables(i) << std::flush;
1540 outFile << " };"
1541 << "\n"
1542 << "\n";
1543
1544
1545 outFile << "// Assignment to maximum vector."
1546 << "\n";
1547 outFile << cv_qual << "double " << prefix << "gXMax[] = {"
1548 << "\n";
1549 for (i = 0; i < fNVariables; i++)
1550 outFile << (i != 0 ? ", " : " ") << fMaxVariables(i) << std::flush;
1551 outFile << " };"
1552 << "\n"
1553 << "\n";
1554
1555
1556 outFile << "// Assignment to coefficients vector."
1557 << "\n";
1558 outFile << cv_qual << "double " << prefix << "gCoefficient[] = {" << std::flush;
1559 for (i = 0; i < fNCoefficients; i++)
1560 outFile << (i != 0 ? "," : "") << "\n"
1561 << " " << fCoefficients(i) << std::flush;
1562 outFile << "\n"
1563 << " };"
1564 << "\n"
1565 << "\n";
1566
1567
1568 outFile << "// Assignment to powers vector."
1569 << "\n"
1570 << "// The powers are stored row-wise, that is"
1571 << "\n"
1572 << "// p_ij = " << prefix << "gPower[i * NVariables + j];"
1573 << "\n";
1574 outFile << cv_qual << "int " << prefix << "gPower[] = {" << std::flush;
1575 for (i = 0; i < fNCoefficients; i++) {
1576 for (j = 0; j < fNVariables; j++) {
1577 if (j != 0)
1578 outFile << std::flush << " ";
1579 else
1580 outFile << "\n"
1581 << " ";
1582 outFile << fPowers[fPowerIndex[i] * fNVariables + j]
1583 << (i == fNCoefficients - 1 && j == fNVariables - 1 ? "" : ",") << std::flush;
1584 }
1585 }
1586 outFile << "\n"
1587 << "};"
1588 << "\n"
1589 << "\n";
1590
1591
1592
1593
1594 outFile << "// "
1595 << "\n"
1596 << "// The " << (isMethod ? "method " : "function ") << " double " << prefix << "MDF(double *x)"
1597 << "\n"
1598 << "// "
1599 << "\n";
1600 outFile << "double " << prefix << "MDF(double *x) {"
1601 << "\n"
1602 << " double returnValue = " << prefix << "gDMean;"
1603 << "\n"
1604 << " int i = 0, j = 0, k = 0;"
1605 << "\n"
1606 << " for (i = 0; i < " << prefix << "gNCoefficients ; i++) {"
1607 << "\n"
1608 << " // Evaluate the ith term in the expansion"
1609 << "\n"
1610 << " double term = " << prefix << "gCoefficient[i];"
1611 << "\n"
1612 << " for (j = 0; j < " << prefix << "gNVariables; j++) {"
1613 << "\n"
1614 << " // Evaluate the polynomial in the jth variable."
1615 << "\n"
1616 << " int power = " << prefix << "gPower[" << prefix << "gNVariables * i + j]; "
1617 << "\n"
1618 << " double p1 = 1, p2 = 0, p3 = 0, r = 0;"
1619 << "\n"
1620 << " double v = 1 + 2. / (" << prefix << "gXMax[j] - " << prefix << "gXMin[j]) * (x[j] - " << prefix
1621 << "gXMax[j]);"
1622 << "\n"
1623 << " // what is the power to use!"
1624 << "\n"
1625 << " switch(power) {"
1626 << "\n"
1627 << " case 1: r = 1; break; "
1628 << "\n"
1629 << " case 2: r = v; break; "
1630 << "\n"
1631 << " default: "
1632 << "\n"
1633 << " p2 = v; "
1634 << "\n"
1635 << " for (k = 3; k <= power; k++) { "
1636 << "\n"
1637 << " p3 = p2 * v;"
1638 << "\n";
1639 if (fPolyType == kLegendre)
1640 outFile << " p3 = ((2 * i - 3) * p2 * v - (i - 2) * p1)"
1641 << " / (i - 1);"
1642 << "\n";
1643 if (fPolyType == kChebyshev)
1644 outFile << " p3 = 2 * v * p2 - p1; "
1645 << "\n";
1646 outFile << " p1 = p2; p2 = p3; "
1647 << "\n"
1648 << " }"
1649 << "\n"
1650 << " r = p3;"
1651 << "\n"
1652 << " }"
1653 << "\n"
1654 << " // multiply this term by the poly in the jth var"
1655 << "\n"
1656 << " term *= r; "
1657 << "\n"
1658 << " }"
1659 << "\n"
1660 << " // Add this term to the final result"
1661 << "\n"
1662 << " returnValue += term;"
1663 << "\n"
1664 << " }"
1665 << "\n"
1666 << " return returnValue;"
1667 << "\n"
1668 << "}"
1669 << "\n"
1670 << "\n";
1671
1672
1673 outFile << "// EOF for " << filename << "\n";
1674
1675
1676 outFile.close();
1677
1678 if (fIsVerbose)
1679 edm::LogInfo("TMultiDimFet") << "done"
1680 << "\n";
1681 }
1682
1683
1684 void TMultiDimFet::Print(Option_t *option) const {
1685
1686
1687
1688
1689
1690
1691
1692
1693
1694
1695 Int_t i = 0;
1696 Int_t j = 0;
1697
1698 TString opt(option);
1699 opt.ToLower();
1700
1701 if (opt.Contains("p")) {
1702
1703 edm::LogInfo("TMultiDimFet") << "User parameters:"
1704 << "\n"
1705 << "----------------"
1706 << "\n"
1707 << " Variables: " << fNVariables << "\n"
1708 << " Data points: " << fSampleSize << "\n"
1709 << " Max Terms: " << fMaxTerms << "\n"
1710 << " Power Limit Parameter: " << fPowerLimit << "\n"
1711 << " Max functions: " << fMaxFunctions << "\n"
1712 << " Max functions to study: " << fMaxStudy << "\n"
1713 << " Max angle (optional): " << fMaxAngle << "\n"
1714 << " Min angle: " << fMinAngle << "\n"
1715 << " Relative Error accepted: " << fMinRelativeError << "\n"
1716 << " Maximum Powers: " << std::flush;
1717 for (i = 0; i < fNVariables; i++)
1718 edm::LogInfo("TMultiDimFet") << " " << fMaxPowers[i] - 1 << std::flush;
1719 edm::LogInfo("TMultiDimFet") << "\n"
1720 << "\n"
1721 << " Parameterisation will be done using " << std::flush;
1722 if (fPolyType == kChebyshev)
1723 edm::LogInfo("TMultiDimFet") << "Chebyshev polynomials"
1724 << "\n";
1725 else if (fPolyType == kLegendre)
1726 edm::LogInfo("TMultiDimFet") << "Legendre polynomials"
1727 << "\n";
1728 else
1729 edm::LogInfo("TMultiDimFet") << "Monomials"
1730 << "\n";
1731 edm::LogInfo("TMultiDimFet") << "\n";
1732 }
1733
1734 if (opt.Contains("s")) {
1735
1736 edm::LogInfo("TMultiDimFet") << "Sample statistics:"
1737 << "\n"
1738 << "------------------"
1739 << "\n"
1740 << " D" << std::flush;
1741 for (i = 0; i < fNVariables; i++)
1742 edm::LogInfo("TMultiDimFet") << " " << std::setw(10) << i + 1 << std::flush;
1743 edm::LogInfo("TMultiDimFet") << "\n"
1744 << " Max: " << std::setw(10) << std::setprecision(7) << fMaxQuantity << std::flush;
1745 for (i = 0; i < fNVariables; i++)
1746 edm::LogInfo("TMultiDimFet") << " " << std::setw(10) << std::setprecision(4) << fMaxVariables(i) << std::flush;
1747 edm::LogInfo("TMultiDimFet") << "\n"
1748 << " Min: " << std::setw(10) << std::setprecision(7) << fMinQuantity << std::flush;
1749 for (i = 0; i < fNVariables; i++)
1750 edm::LogInfo("TMultiDimFet") << " " << std::setw(10) << std::setprecision(4) << fMinVariables(i) << std::flush;
1751 edm::LogInfo("TMultiDimFet") << "\n"
1752 << " Mean: " << std::setw(10) << std::setprecision(7) << fMeanQuantity << std::flush;
1753 for (i = 0; i < fNVariables; i++)
1754 edm::LogInfo("TMultiDimFet") << " " << std::setw(10) << std::setprecision(4) << fMeanVariables(i) << std::flush;
1755 edm::LogInfo("TMultiDimFet") << "\n"
1756 << " Function Sum Squares: " << fSumSqQuantity << "\n"
1757 << "\n";
1758 }
1759
1760 if (opt.Contains("r")) {
1761 edm::LogInfo("TMultiDimFet") << "Results of Parameterisation:"
1762 << "\n"
1763 << "----------------------------"
1764 << "\n"
1765 << " Total reduction of square residuals " << fSumSqResidual << "\n"
1766 << " Relative precision obtained: " << fPrecision << "\n"
1767 << " Error obtained: " << fError << "\n"
1768 << " Multiple correlation coefficient: " << fCorrelationCoeff << "\n"
1769 << " Reduced Chi square over sample: " << fChi2 / (fSampleSize - fNCoefficients)
1770 << "\n"
1771 << " Maximum residual value: " << fMaxResidual << "\n"
1772 << " Minimum residual value: " << fMinResidual << "\n"
1773 << " Estimated root mean square: " << fRMS << "\n"
1774 << " Maximum powers used: " << std::flush;
1775 for (j = 0; j < fNVariables; j++)
1776 edm::LogInfo("TMultiDimFet") << fMaxPowersFinal[j] << " " << std::flush;
1777 edm::LogInfo("TMultiDimFet") << "\n"
1778 << " Function codes of candidate functions."
1779 << "\n"
1780 << " 1: considered,"
1781 << " 2: too little contribution,"
1782 << " 3: accepted." << std::flush;
1783 for (i = 0; i < fMaxFunctions; i++) {
1784 if (i % 60 == 0)
1785 edm::LogInfo("TMultiDimFet") << "\n"
1786 << " " << std::flush;
1787 else if (i % 10 == 0)
1788 edm::LogInfo("TMultiDimFet") << " " << std::flush;
1789 edm::LogInfo("TMultiDimFet") << fFunctionCodes[i];
1790 }
1791 edm::LogInfo("TMultiDimFet") << "\n"
1792 << " Loop over candidates stopped because " << std::flush;
1793 switch (fParameterisationCode) {
1794 case PARAM_MAXSTUDY:
1795 edm::LogInfo("TMultiDimFet") << "max allowed studies reached"
1796 << "\n";
1797 break;
1798 case PARAM_SEVERAL:
1799 edm::LogInfo("TMultiDimFet") << "all candidates considered several times"
1800 << "\n";
1801 break;
1802 case PARAM_RELERR:
1803 edm::LogInfo("TMultiDimFet") << "wanted relative error obtained"
1804 << "\n";
1805 break;
1806 case PARAM_MAXTERMS:
1807 edm::LogInfo("TMultiDimFet") << "max number of terms reached"
1808 << "\n";
1809 break;
1810 default:
1811 edm::LogInfo("TMultiDimFet") << "some unknown reason"
1812 << "\n";
1813 break;
1814 }
1815 edm::LogInfo("TMultiDimFet") << "\n";
1816 }
1817
1818 if (opt.Contains("f")) {
1819 edm::LogInfo("TMultiDimFet") << "Results of Fit:"
1820 << "\n"
1821 << "---------------"
1822 << "\n"
1823 << " Test sample size: " << fTestSampleSize << "\n"
1824 << " Multiple correlation coefficient: " << fTestCorrelationCoeff << "\n"
1825 << " Relative precision obtained: " << fTestPrecision << "\n"
1826 << " Error obtained: " << fTestError << "\n"
1827 << " Reduced Chi square over sample: " << fChi2 / (fSampleSize - fNCoefficients)
1828 << "\n"
1829 << "\n";
1830
1831
1832
1833
1834
1835
1836 }
1837
1838 if (opt.Contains("c")) {
1839 edm::LogInfo("TMultiDimFet") << "Coefficients:"
1840 << "\n"
1841 << "-------------"
1842 << "\n"
1843 << " # Value Error Powers"
1844 << "\n"
1845 << " ---------------------------------------"
1846 << "\n";
1847 for (i = 0; i < fNCoefficients; i++) {
1848 edm::LogInfo("TMultiDimFet") << " " << std::setw(3) << i << " " << std::setw(12) << fCoefficients(i) << " "
1849 << std::setw(12) << fCoefficientsRMS(i) << " " << std::flush;
1850 for (j = 0; j < fNVariables; j++)
1851 edm::LogInfo("TMultiDimFet") << " " << std::setw(3) << fPowers[fPowerIndex[i] * fNVariables + j] - 1
1852 << std::flush;
1853 edm::LogInfo("TMultiDimFet") << "\n";
1854 }
1855 edm::LogInfo("TMultiDimFet") << "\n";
1856 }
1857 if (opt.Contains("k") && fCorrelationMatrix.IsValid()) {
1858 edm::LogInfo("TMultiDimFet") << "Correlation Matrix:"
1859 << "\n"
1860 << "-------------------";
1861 fCorrelationMatrix.Print();
1862 }
1863
1864 if (opt.Contains("m")) {
1865 edm::LogInfo("TMultiDimFet") << std::setprecision(25);
1866 edm::LogInfo("TMultiDimFet") << "Parameterization:"
1867 << "\n"
1868 << "-----------------"
1869 << "\n"
1870 << " Normalised variables: "
1871 << "\n";
1872 for (i = 0; i < fNVariables; i++)
1873 edm::LogInfo("TMultiDimFet") << "\ty" << i << "\t:= 1 + 2 * (x" << i << " - " << fMaxVariables(i) << ") / ("
1874 << fMaxVariables(i) << " - " << fMinVariables(i) << ")"
1875 << "\n";
1876 edm::LogInfo("TMultiDimFet") << "\n"
1877 << " f[";
1878 for (i = 0; i < fNVariables; i++) {
1879 edm::LogInfo("TMultiDimFet") << "y" << i;
1880 if (i != fNVariables - 1)
1881 edm::LogInfo("TMultiDimFet") << ", ";
1882 }
1883 edm::LogInfo("TMultiDimFet") << "] := ";
1884 for (Int_t i = 0; i < fNCoefficients; i++) {
1885 if (i != 0)
1886 edm::LogInfo("TMultiDimFet") << " " << (fCoefficients(i) < 0 ? "- " : "+ ") << TMath::Abs(fCoefficients(i));
1887 else
1888 edm::LogInfo("TMultiDimFet") << fCoefficients(i);
1889 for (Int_t j = 0; j < fNVariables; j++) {
1890 Int_t p = fPowers[fPowerIndex[i] * fNVariables + j];
1891 switch (p) {
1892 case 1:
1893 break;
1894 case 2:
1895 edm::LogInfo("TMultiDimFet") << " * y" << j;
1896 break;
1897 default:
1898 switch (fPolyType) {
1899 case kLegendre:
1900 edm::LogInfo("TMultiDimFet") << " * L" << p - 1 << "(y" << j << ")";
1901 break;
1902 case kChebyshev:
1903 edm::LogInfo("TMultiDimFet") << " * C" << p - 1 << "(y" << j << ")";
1904 break;
1905 default:
1906 edm::LogInfo("TMultiDimFet") << " * y" << j << "^" << p - 1;
1907 break;
1908 }
1909 }
1910 }
1911 }
1912 edm::LogInfo("TMultiDimFet") << "\n";
1913 }
1914 }
1915
1916
1917 void TMultiDimFet::PrintPolynomialsSpecial(Option_t *option) const {
1918
1919
1920 Int_t i = 0;
1921
1922
1923 TString opt(option);
1924 opt.ToLower();
1925
1926 if (opt.Contains("m")) {
1927 edm::LogInfo("TMultiDimFet") << std::setprecision(25);
1928 edm::LogInfo("TMultiDimFet") << "Parameterization:"
1929 << "\n"
1930 << "-----------------"
1931 << "\n"
1932 << " Normalised variables: "
1933 << "\n";
1934 for (i = 0; i < fNVariables; i++)
1935 edm::LogInfo("TMultiDimFet") << "\tdouble y" << i << "\t=1+2*(x" << i << "-" << fMaxVariables(i) << ")/("
1936 << fMaxVariables(i) << "-" << fMinVariables(i) << ");"
1937 << "\n";
1938 edm::LogInfo("TMultiDimFet") << "\n"
1939 << " f[";
1940 for (i = 0; i < fNVariables; i++) {
1941 edm::LogInfo("TMultiDimFet") << "y" << i;
1942 if (i != fNVariables - 1)
1943 edm::LogInfo("TMultiDimFet") << ", ";
1944 }
1945 edm::LogInfo("TMultiDimFet") << "] := " << fMeanQuantity << " + ";
1946 for (Int_t i = 0; i < fNCoefficients; i++) {
1947 if (i != 0)
1948 edm::LogInfo("TMultiDimFet") << " " << (fCoefficients(i) < 0 ? "-" : "+") << TMath::Abs(fCoefficients(i));
1949 else
1950 edm::LogInfo("TMultiDimFet") << fCoefficients(i);
1951 for (Int_t j = 0; j < fNVariables; j++) {
1952 Int_t p = fPowers[fPowerIndex[i] * fNVariables + j];
1953 switch (p) {
1954 case 1:
1955 break;
1956 case 2:
1957 edm::LogInfo("TMultiDimFet") << "*y" << j;
1958 break;
1959 default:
1960 switch (fPolyType) {
1961 case kLegendre:
1962 edm::LogInfo("TMultiDimFet") << "*Leg(" << p - 1 << ",y" << j << ")";
1963 break;
1964 case kChebyshev:
1965 edm::LogInfo("TMultiDimFet") << "*C" << p - 1 << "(y" << j << ")";
1966 break;
1967 default:
1968 edm::LogInfo("TMultiDimFet") << "*y" << j << "**" << p - 1;
1969 break;
1970 }
1971 }
1972 }
1973 edm::LogInfo("TMultiDimFet") << "\n";
1974 }
1975 edm::LogInfo("TMultiDimFet") << "\n";
1976 }
1977 }
1978
1979
1980 Bool_t TMultiDimFet::Select(const Int_t *) {
1981
1982
1983
1984
1985
1986
1987
1988
1989
1990 return kTRUE;
1991 }
1992
1993
1994 void TMultiDimFet::SetMaxAngle(Double_t ang) {
1995
1996
1997
1998
1999
2000 if (ang >= 90 || ang < 0) {
2001 Warning("SetMaxAngle", "angle must be in [0,90)");
2002 return;
2003 }
2004
2005 fMaxAngle = ang;
2006 }
2007
2008
2009 void TMultiDimFet::SetMinAngle(Double_t ang) {
2010
2011
2012
2013
2014 if (ang > 90 || ang <= 0) {
2015 Warning("SetMinAngle", "angle must be in [0,90)");
2016 return;
2017 }
2018
2019 fMinAngle = ang;
2020 }
2021
2022
2023 void TMultiDimFet::SetPowers(const Int_t *powers, Int_t terms) {
2024
2025
2026
2027
2028
2029
2030 fIsUserFunction = kTRUE;
2031 fMaxFunctions = terms;
2032 fMaxTerms = terms;
2033 fMaxStudy = terms;
2034 fMaxFunctionsTimesNVariables = fMaxFunctions * fNVariables;
2035 fPowers.resize(fMaxFunctions * fNVariables);
2036 Int_t i, j;
2037 for (i = 0; i < fMaxFunctions; i++)
2038 for (j = 0; j < fNVariables; j++)
2039 fPowers[i * fNVariables + j] = powers[i * fNVariables + j] + 1;
2040 }
2041
2042
2043 void TMultiDimFet::SetPowerLimit(Double_t limit) {
2044
2045
2046
2047
2048 fPowerLimit = limit;
2049 }
2050
2051
2052 void TMultiDimFet::SetMaxPowers(const Int_t *powers) {
2053
2054
2055
2056 if (!powers)
2057 return;
2058
2059 for (Int_t i = 0; i < fNVariables; i++)
2060 fMaxPowers[i] = powers[i] + 1;
2061 }
2062
2063
2064 void TMultiDimFet::SetMinRelativeError(Double_t error) {
2065
2066
2067
2068
2069 fMinRelativeError = error;
2070 }
2071
2072
2073 Bool_t TMultiDimFet::TestFunction(Double_t squareResidual, Double_t dResidur) {
2074
2075
2076
2077
2078
2079 if (fNCoefficients != 0) {
2080
2081 if (fMaxAngle == 0) {
2082
2083 if (dResidur < squareResidual / (fMaxTerms - fNCoefficients + 1 + 1E-10)) {
2084 return kFALSE;
2085 }
2086 } else {
2087
2088
2089 if (TMath::Sqrt(dResidur / fSumSqAvgQuantity) < TMath::Cos(fMaxAngle * DEGRAD)) {
2090 return kFALSE;
2091 }
2092 }
2093 }
2094
2095 return kTRUE;
2096 }
2097
2098
2099
2100
2101
2102
2103
2104
2105
2106
2107
2108