Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-10-04 22:55:07

0001 /****************************************************************************
0002 *
0003 * This is a part of TOTEM ofFileline software.
0004 *
0005 * Based on TMultiDimFet.cc (from ROOT) by Christian Holm Christense.
0006 *
0007 * Authors:
0008 *   Hubert Niewiadomski
0009 *   Jan Kašpar (jan.kaspar@gmail.com)
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 //static void mdfHelper(int&, double*, double&, double*, int);
0042 
0043 //____________________________________________________________________
0044 ClassImp(TMultiDimFet);
0045 
0046 //____________________________________________________________________
0047 // Static instance. Used with mdfHelper and TMinuit
0048 //TMultiDimFet* TMultiDimFet::fgInstance = nullptr;
0049 
0050 //____________________________________________________________________
0051 TMultiDimFet::TMultiDimFet() {
0052   // Empty CTOR. Do not use
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   //fFitter                  = nullptr;
0080   //fgInstance               = nullptr;
0081 }
0082 
0083 const TMultiDimFet &TMultiDimFet::operator=(const TMultiDimFet &in) {
0084   if (this == &in) {
0085     return *this;
0086   }
0087 
0088   fMeanQuantity = in.fMeanQuantity;  // Mean of dependent quantity
0089 
0090   fMaxQuantity = 0.0;       //! Max value of dependent quantity
0091   fMinQuantity = 0.0;       //! Min value of dependent quantity
0092   fSumSqQuantity = 0.0;     //! SumSquare of dependent quantity
0093   fSumSqAvgQuantity = 0.0;  //! Sum of squares away from mean
0094 
0095   fNVariables = in.fNVariables;  // Number of independent variables
0096 
0097   fMaxVariables.ResizeTo(in.fMaxVariables.GetLwb(), in.fMaxVariables.GetUpb());
0098   fMaxVariables = in.fMaxVariables;  // max value of independent variables
0099 
0100   fMinVariables.ResizeTo(in.fMinVariables.GetLwb(), in.fMinVariables.GetUpb());
0101   fMinVariables = in.fMinVariables;  // min value of independent variables
0102 
0103   fSampleSize = 0;      //! Size of training sample
0104   fTestSampleSize = 0;  //! Size of test sample
0105   fMinAngle = 1;        //! Min angle for acepting new function
0106   fMaxAngle = 0.0;      //! Max angle for acepting new function
0107 
0108   fMaxTerms = in.fMaxTerms;  // Max terms expected in final expr.
0109 
0110   fMinRelativeError = 0.0;  //! Min relative error accepted
0111 
0112   fMaxPowers.clear();  //! [fNVariables] maximum powers
0113   fPowerLimit = 1;     //! Control parameter
0114 
0115   fMaxFunctions = in.fMaxFunctions;  // max number of functions
0116 
0117   fFunctionCodes.clear();  //! [fMaxFunctions] acceptance code
0118   fMaxStudy = 0;           //! max functions to study
0119 
0120   fMaxPowersFinal.clear();  //! [fNVariables] maximum powers from fit;
0121 
0122   fMaxFunctionsTimesNVariables = in.fMaxFunctionsTimesNVariables;  // fMaxFunctionsTimesNVariables
0123   fPowers = in.fPowers;
0124 
0125   fPowerIndex = in.fPowerIndex;  // [fMaxTerms] Index of accepted powers
0126 
0127   fMaxResidual = 0.0;    //! Max redsidual value
0128   fMinResidual = 0.0;    //! Min redsidual value
0129   fMaxResidualRow = 0;   //! Row giving max residual
0130   fMinResidualRow = 0;   //! Row giving min residual
0131   fSumSqResidual = 0.0;  //! Sum of Square residuals
0132 
0133   fNCoefficients = in.fNCoefficients;  // Dimension of model coefficients
0134 
0135   fCoefficients.ResizeTo(in.fCoefficients.GetLwb(), in.fCoefficients.GetUpb());
0136   fCoefficients = in.fCoefficients;  // Vector of the final coefficients
0137 
0138   fRMS = 0.0;                   //! Root mean square of fit
0139   fChi2 = 0.0;                  //! Chi square of fit
0140   fParameterisationCode = 0;    //! Exit code of parameterisation
0141   fError = 0.0;                 //! Error from parameterization
0142   fTestError = 0.0;             //! Error from test
0143   fPrecision = 0.0;             //! Relative precision of param
0144   fTestPrecision = 0.0;         //! Relative precision of test
0145   fCorrelationCoeff = 0.0;      //! Multi Correlation coefficient
0146   fTestCorrelationCoeff = 0.0;  //! Multi Correlation coefficient
0147   fHistograms = nullptr;        //! List of histograms
0148   fHistogramMask = 0;           //! Bit pattern of hisograms used
0149   //fFitter = nullptr;            //! Fit object (MINUIT)
0150 
0151   fPolyType = in.fPolyType;                // Type of polynomials to use
0152   fShowCorrelation = in.fShowCorrelation;  // print correlation matrix
0153   fIsUserFunction = in.fIsUserFunction;    // Flag for user defined function
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   // Constructor
0168   // Second argument is the type of polynomials to use in
0169   // parameterisation, one of:
0170   //      TMultiDimFet::kMonomials
0171   //      TMultiDimFet::kChebyshev
0172   //      TMultiDimFet::kLegendre
0173   //
0174   // Options:
0175   //   K      Compute (k)correlation matrix
0176   //   V      Be verbose
0177   //
0178   // Default is no options.
0179   //
0180 
0181   //fgInstance = this;
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   //fFitter                 = nullptr;
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   // Add a row consisting of fNVariables independent variables, the
0238   // known, dependent quantity, and optionally, the square error in
0239   // the dependent quantity, to the training sample to be used for the
0240   // parameterization.
0241   // The mean of the variables and quantity is calculated on the fly,
0242   // as outlined in TPrincipal::AddRow.
0243   // This sample should be representive of the problem at hand.
0244   // Please note, that if no error is given Poisson statistics is
0245   // assumed and the square error is set to the value of dependent
0246   // quantity.  See also the
0247   // Begin_Html<a href="#TMultiDimFet:description">class description</a>End_Html
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   // If the vector isn't big enough to hold the new data, then
0267   // expand the vector by half it's size.
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   // Store the value
0275   fQuantity(fSampleSize - 1) = D;
0276   fSqError(fSampleSize - 1) = (E == 0 ? D : E);
0277 
0278   // Store data point in internal vector
0279   // If the vector isn't big enough to hold the new data, then
0280   // expand the vector by half it's size
0281   size = fVariables.GetNrows();
0282   if (fSampleSize * fNVariables > size)
0283     fVariables.ResizeTo(size + size / 2);
0284 
0285   // Increment the data point counter
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       // Update the maximum value for this component
0297       if (x[i] >= fMaxVariables(i))
0298         fMaxVariables(i) = x[i];
0299 
0300       // Update the minimum value for this component
0301       if (x[i] <= fMinVariables(i))
0302         fMinVariables(i) = x[i];
0303     }
0304 
0305     // Store the data.
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   // Add a row consisting of fNVariables independent variables, the
0314   // known, dependent quantity, and optionally, the square error in
0315   // the dependent quantity, to the test sample to be used for the
0316   // test of the parameterization.
0317   // This sample needn't be representive of the problem at hand.
0318   // Please note, that if no error is given Poisson statistics is
0319   // assumed and the square error is set to the value of dependent
0320   // quantity.  See also the
0321   // Begin_Html<a href="#TMultiDimFet:description">class description</a>End_Html
0322   if (fTestSampleSize++ == 0) {
0323     fTestQuantity.ResizeTo(fNVariables);
0324     fTestSqError.ResizeTo(fNVariables);
0325     fTestVariables.ResizeTo(fNVariables * 100);
0326   }
0327 
0328   // If the vector isn't big enough to hold the new data, then
0329   // expand the vector by half it's size.
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   // Store the value
0337   fTestQuantity(fTestSampleSize - 1) = D;
0338   fTestSqError(fTestSampleSize - 1) = (E == 0 ? D : E);
0339 
0340   // Store data point in internal vector
0341   // If the vector isn't big enough to hold the new data, then
0342   // expand the vector by half it's size
0343   size = fTestVariables.GetNrows();
0344   if (fTestSampleSize * fNVariables > size)
0345     fTestVariables.ResizeTo(size + size / 2);
0346 
0347   // Increment the data point counter
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   // Clear internal structures and variables
0363   Int_t i, j, n = fNVariables, m = fMaxFunctions;
0364 
0365   // Training sample, dependent quantity
0366   fQuantity.Zero();
0367   fSqError.Zero();
0368   fMeanQuantity = 0;
0369   fMaxQuantity = 0;
0370   fMinQuantity = 0;
0371   fSumSqQuantity = 0;
0372   fSumSqAvgQuantity = 0;
0373 
0374   // Training sample, independent variables
0375   fVariables.Zero();
0376   fNVariables = 0;
0377   fSampleSize = 0;
0378   fMeanVariables.Zero();
0379   fMaxVariables.Zero();
0380   fMinVariables.Zero();
0381 
0382   // Test sample
0383   fTestQuantity.Zero();
0384   fTestSqError.Zero();
0385   fTestVariables.Zero();
0386   fTestSampleSize = 0;
0387 
0388   // Functions
0389   fFunctions.Zero();
0390   fMaxFunctions = 0;
0391   fMaxStudy = 0;
0392   fMaxFunctionsTimesNVariables = 0;
0393   fOrthFunctions.Zero();
0394   fOrthFunctionNorms.Zero();
0395 
0396   // Control parameters
0397   fMinRelativeError = 0;
0398   fMinAngle = 0;
0399   fMaxAngle = 0;
0400   fMaxTerms = 0;
0401 
0402   // Powers
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   // Residuals
0412   fMaxResidual = 0;
0413   fMinResidual = 0;
0414   fMaxResidualRow = 0;
0415   fMinResidualRow = 0;
0416   fSumSqResidual = 0;
0417 
0418   // Fit
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   // Coefficients
0430   fCoefficients.Zero();
0431   fCoefficientsRMS.Zero();
0432   fResiduals.Zero();
0433   fHistograms->Clear(option);
0434 
0435   // Options
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   // Evaluate parameterization at point x. Optional argument coeff is
0444   // a vector of coefficients for the parameterisation, fNCoefficients
0445   // elements long.
0446   //   int fMaxFunctionsTimesNVariables = fMaxFunctions * fNVariables;
0447   Double_t returnValue = fMeanQuantity;
0448   Double_t term = 0;
0449   Int_t i, j;
0450 
0451   for (i = 0; i < fNCoefficients; i++) {
0452     // Evaluate the ith term in the expansion
0453     term = (coeff ? coeff[i] : fCoefficients(i));
0454     for (j = 0; j < fNVariables; j++) {
0455       // Evaluate the factor (polynomial) in the j-th variable.
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     // Add this term to the final result
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   // PRIVATE METHOD:
0517   // Calculate the control parameter from the passed powers
0518   Double_t s = 0;
0519   Double_t epsilon = 1e-6;  // a small number
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   // PRIVATE METHOD:
0530   // Evaluate function with power p at variable value x
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   // Find the parameterization
0564   //
0565   // Options:
0566   //     None so far
0567   //
0568   // For detailed description of what this entails, please refer to the
0569   // Begin_Html<a href="#TMultiDimFet:description">class description</a>End_Html
0570   MakeNormalized();
0571   MakeCandidates();
0572   MakeParameterization();
0573   MakeCoefficients();
0574   ReducePolynomial(precision);
0575 }
0576 
0577 //____________________________________________________________________
0578 /*
0579 void TMultiDimFet::Fit(Option_t *option)
0580 {
0581    // Try to fit the found parameterisation to the test sample.
0582    //
0583    // Options
0584    //     M     use Minuit to improve coefficients
0585    //
0586    // Also, refer to
0587    // Begin_Html<a href="#TMultiDimFet:description">class description</a>End_Html
0588 
0589    //fgInstance = this;
0590 
0591    Int_t i, j;
0592    Double_t*      x    = new Double_t[fNVariables];
0593    Double_t  sumSqD    = 0;
0594    Double_t    sumD    = 0;
0595    Double_t  sumSqR    = 0;
0596    Double_t    sumR    = 0;
0597 
0598    // Calculate the residuals over the test sample
0599    for (i = 0; i < fTestSampleSize; i++) {
0600       for (j = 0; j < fNVariables; j++)
0601          x[j] = fTestVariables(i * fNVariables + j);
0602       Double_t res =  fTestQuantity(i) - Eval(x);
0603       sumD         += fTestQuantity(i);
0604       sumSqD       += fTestQuantity(i) * fTestQuantity(i);
0605       sumR         += res;
0606       sumSqR       += res * res;
0607       if (TESTBIT(fHistogramMask,HIST_RTEST))
0608          ((TH1D*)fHistograms->FindObject("res_test"))->Fill(res);
0609    }
0610    Double_t dAvg         = sumSqD - (sumD * sumD) / fTestSampleSize;
0611    Double_t rAvg         = sumSqR - (sumR * sumR) / fTestSampleSize;
0612    fTestCorrelationCoeff = (dAvg - rAvg) / dAvg;
0613    fTestError            = sumSqR;
0614    fTestPrecision        = sumSqR / sumSqD;
0615 
0616    TString opt(option);
0617    opt.ToLower();
0618 
0619    if (!opt.Contains("m"))
0620       MakeChi2();
0621 
0622    if (fNCoefficients * 50 > fTestSampleSize)
0623       Warning("Fit", "test sample is very small");
0624 
0625    if (!opt.Contains("m"))
0626       return;
0627 
0628    fFitter = TVirtualFitter::Fitter(nullptr,fNCoefficients);
0629    fFitter->SetFCN(mdfHelper);
0630 
0631    const Int_t  maxArgs = 16;
0632    Int_t           args = 1;
0633    Double_t*   arglist  = new Double_t[maxArgs];
0634    arglist[0]           = -1;
0635    fFitter->ExecuteCommand("SET PRINT",arglist,args);
0636 
0637    for (i = 0; i < fNCoefficients; i++) {
0638       Double_t startVal = fCoefficients(i);
0639       Double_t startErr = fCoefficientsRMS(i);
0640       fFitter->SetParameter(i, Form("coeff%02d",i),
0641          startVal, startErr, 0, 0);
0642    }
0643 
0644    args                 = 1;
0645    fFitter->ExecuteCommand("MIGRAD",arglist,args);
0646 
0647    for (i = 0; i < fNCoefficients; i++) {
0648       Double_t val = 0, err = 0, low = 0, high = 0;
0649       fFitter->GetParameter(i, Form("coeff%02d",i),
0650          val, err, low, high);
0651       fCoefficients(i)    = val;
0652       fCoefficientsRMS(i) = err;
0653    }
0654 }
0655 */
0656 
0657 //____________________________________________________________________
0658 void TMultiDimFet::MakeCandidates() {
0659   // PRIVATE METHOD:
0660   // Create list of candidate functions for the parameterisation. See
0661   // also
0662   // Begin_Html<a href="#TMultiDimFet:description">class description</a>End_Html
0663   Int_t i = 0;
0664   Int_t j = 0;
0665   Int_t k = 0;
0666 
0667   // The temporary array to store the powers in. We don't need to
0668   // initialize this array however.
0669   assert(fNVariables > 0);
0670   Int_t *powers = new Int_t[fNVariables * fMaxFunctions];
0671 
0672   // store of `control variables'
0673   Double_t *control = new Double_t[fMaxFunctions];
0674 
0675   // We've better initialize the variables
0676   Int_t *iv = new Int_t[fNVariables];
0677   for (i = 0; i < fNVariables; i++)
0678     iv[i] = 1;
0679 
0680   if (!fIsUserFunction) {
0681     // Number of funcs selected
0682     Int_t numberFunctions = 0;
0683 
0684     while (kTRUE) {
0685       // Get the control value for this function
0686       Double_t s = EvalControl(iv);
0687 
0688       if (s <= fPowerLimit) {
0689         // Call over-loadable method Select, as to allow the user to
0690         // interfere with the selection of functions.
0691         if (Select(iv)) {
0692           numberFunctions++;
0693 
0694           // If we've reached the user defined limit of how many
0695           // functions we can consider, break out of the loop
0696           if (numberFunctions > fMaxFunctions)
0697             break;
0698 
0699           // Store the control value, so we can sort array of powers
0700           // later on
0701           control[numberFunctions - 1] = Int_t(1.0e+6 * s);
0702 
0703           // Store the powers in powers array.
0704           for (i = 0; i < fNVariables; i++) {
0705             j = (numberFunctions - 1) * fNVariables + i;
0706             powers[j] = iv[i];
0707           }
0708         }  // if (Select())
0709       }  // if (s <= fPowerLimit)
0710 
0711       for (i = 0; i < fNVariables; i++)
0712         if (iv[i] < fMaxPowers[i])
0713           break;
0714 
0715       // If all variables have reached their maximum power, then we
0716       // break out of the loop
0717       if (i == fNVariables) {
0718         fMaxFunctions = numberFunctions;
0719         fMaxFunctionsTimesNVariables = fMaxFunctions * fNVariables;
0720         break;
0721       }
0722 
0723       // Next power in variable i
0724       iv[i]++;
0725 
0726       for (j = 0; j < i; j++)
0727         iv[j] = 1;
0728     }  // while (kTRUE)
0729   } else {
0730     // In case the user gave an explicit function
0731     for (i = 0; i < fMaxFunctions; i++) {
0732       // Copy the powers to working arrays
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   // Now we need to sort the powers according to least `control
0743   // variable'
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   // Calculate Chi square over either the test sample. The optional
0782   // argument coeff is a vector of coefficients to use in the
0783   // evaluation of the parameterisation. If coeff == 0, then the found
0784   // coefficients is used.
0785   // Used my MINUIT for fit (see TMultDimFit::Fit)
0786   fChi2 = 0;
0787   Int_t i, j;
0788   Double_t *x = new Double_t[fNVariables];
0789   for (i = 0; i < fTestSampleSize; i++) {
0790     // Get the stored point
0791     for (j = 0; j < fNVariables; j++)
0792       x[j] = fTestVariables(i * fNVariables + j);
0793 
0794     // Evaluate function. Scale to shifted values
0795     Double_t f = Eval(x, coeff);
0796 
0797     // Calculate contribution to Chic square
0798     fChi2 += 1. / TMath::Max(fTestSqError(i), 1e-20) * (fTestQuantity(i) - f) * (fTestQuantity(i) - f);
0799   }
0800 
0801   // Clean up
0802   delete[] x;
0803 
0804   return fChi2;
0805 }
0806 
0807 //____________________________________________________________________
0808 void TMultiDimFet::MakeCode(const char *filename, Option_t *option) {
0809   // Generate the file <filename> with .C appended if argument doesn't
0810   // end in .cxx or .C. The contains the implementation of the
0811   // function:
0812   //
0813   //   Double_t <funcname>(Double_t *x)
0814   //
0815   // which does the same as TMultiDimFet::Eval. Please refer to this
0816   // method.
0817   //
0818   // Further, the static variables:
0819   //
0820   //     Int_t    gNVariables
0821   //     Int_t    gNCoefficients
0822   //     Double_t gDMean
0823   //     Double_t gXMean[]
0824   //     Double_t gXMin[]
0825   //     Double_t gXMax[]
0826   //     Double_t gCoefficient[]
0827   //     Int_t    gPower[]
0828   //
0829   // are initialized. The only ROOT header file needed is Rtypes.h
0830   //
0831   // See TMultiDimFet::MakeRealCode for a list of options
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   // PRIVATE METHOD:
0842   // Compute the errors on the coefficients. For this to be done, the
0843   // curvature matrix of the non-orthogonal functions, is computed.
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   // Build the curvature matrix
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   // Calculate Chi Square
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   // Invert the curvature matrix
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   // PRIVATE METHOD:
0891   // Invert the model matrix B, and compute final coefficients. For a
0892   // more thorough discussion of what this means, please refer to the
0893   // Begin_Html<a href="#TMultiDimFet:description">class description</a>End_Html
0894   //
0895   // First we invert the lower triangle matrix fOrthCurvatureMatrix
0896   // and store the inverted matrix in the upper triangle.
0897 
0898   Int_t i = 0, j = 0;
0899   Int_t col = 0, row = 0;
0900 
0901   // Invert the B matrix
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   // Compute the final coefficients
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   // Compute the final residuals
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   // Compute the max and minimum, and squared sum of the evaluated
0930   // residuals
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   // If we use histograms, fill some more
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   }  // If histograms
0962 }
0963 
0964 //____________________________________________________________________
0965 void TMultiDimFet::MakeCorrelation() {
0966   // PRIVATE METHOD:
0967   // Compute the correlation matrix
0968   if (!fShowCorrelation)
0969     return;
0970 
0971   fCorrelationMatrix.ResizeTo(fNVariables, fNVariables + 1);
0972 
0973   Double_t d2 = 0;
0974   Double_t ddotXi = 0;   // G.Q. needs to be reinitialized in the loop over i fNVariables
0975   Double_t xiNorm = 0;   // G.Q. needs to be reinitialized in the loop over i fNVariables
0976   Double_t xidotXj = 0;  // G.Q. needs to be reinitialized in the loop over j fNVariables
0977   Double_t xjNorm = 0;   // G.Q. needs to be reinitialized in the loop over j fNVariables
0978 
0979   Int_t i, j, k, l, m;  // G.Q. added m variable
0980   for (i = 0; i < fSampleSize; i++)
0981     d2 += fQuantity(i) * fQuantity(i);
0982 
0983   for (i = 0; i < fNVariables; i++) {
0984     ddotXi = 0.;  // G.Q. reinitialisation
0985     xiNorm = 0.;  // G.Q. reinitialisation
0986     for (j = 0; j < fSampleSize; j++) {
0987       // Index of sample j of variable i
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.;  // G.Q. reinitialisation
0996       xjNorm = 0.;   // G.Q. reinitialisation
0997       for (k = 0; k < fSampleSize; k++) {
0998         // Index of sample j of variable i
0999         // l =  j * fNVariables + k;  // G.Q.
1000         l = k * fNVariables + j;  // G.Q.
1001         m = k * fNVariables + i;  // G.Q.
1002         // G.Q.        xidotXj += (fVariables(i) - fMeanVariables(i))
1003         // G.Q.          * (fVariables(l) - fMeanVariables(j));
1004         xidotXj +=
1005             (fVariables(m) - fMeanVariables(i)) * (fVariables(l) - fMeanVariables(j));  // G.Q. modified index for Xi
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   // PRIVATE METHOD:
1016   // Make Gram-Schmidt orthogonalisation. The class description gives
1017   // a thorough account of this algorithm, as well as
1018   // references. Please refer to the
1019   // Begin_Html<a href="#TMultiDimFet:description">class description</a>End_Html
1020 
1021   // calculate w_i, that is, evaluate the current function at data
1022   // point i
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     // First, however, we need to calculate f_fNCoefficients
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     // Calculate f dot f in f2
1040     f2 += fFunctions(fNCoefficients, j) * fFunctions(fNCoefficients, j);
1041     // Assign to w_fNCoefficients f_fNCoefficients
1042     fOrthFunctions(fNCoefficients, j) = fFunctions(fNCoefficients, j);
1043   }
1044 
1045   // the first column of w is equal to f
1046   for (j = 0; j < fNCoefficients; j++) {
1047     Double_t fdw = 0;
1048     // Calculate (f_fNCoefficients dot w_j) / w_j^2
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     // and subtract it from the current value of w_ij
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     // calculate squared length of w_fNCoefficients
1061     fOrthFunctionNorms(fNCoefficients) += fOrthFunctions(fNCoefficients, j) * fOrthFunctions(fNCoefficients, j);
1062 
1063     // calculate D dot w_fNCoefficients in A
1064     fOrthCoefficients(fNCoefficients) += fQuantity(j) * fOrthFunctions(fNCoefficients, j);
1065   }
1066 
1067   // First test, but only if didn't user specify
1068   if (!fIsUserFunction)
1069     if (TMath::Sqrt(fOrthFunctionNorms(fNCoefficients) / (f2 + 1e-10)) < TMath::Sin(fMinAngle * DEGRAD))
1070       return 0;
1071 
1072   // The result found by this code for the first residual is always
1073   // much less then the one found be MUDIFI. That's because it's
1074   // supposed to be. The cause is the improved precision of Double_t
1075   // over DOUBLE PRECISION!
1076   fOrthCurvatureMatrix(fNCoefficients, fNCoefficients) = 1;
1077   Double_t b = fOrthCoefficients(fNCoefficients);
1078   fOrthCoefficients(fNCoefficients) /= fOrthFunctionNorms(fNCoefficients);
1079 
1080   // Calculate the residual from including this fNCoefficients.
1081   Double_t dResidur = fOrthCoefficients(fNCoefficients) * b;
1082 
1083   return dResidur;
1084 }
1085 
1086 //____________________________________________________________________
1087 void TMultiDimFet::MakeHistograms(Option_t *option) {
1088   // Make histograms of the result of the analysis. This message
1089   // should be sent after having read all data points, but before
1090   // finding the parameterization
1091   //
1092   // Options:
1093   //     A         All the below
1094   //     X         Original independent variables
1095   //     D         Original dependent variables
1096   //     N         Normalised independent variables
1097   //     S         Shifted dependent variables
1098   //     R1        Residuals versus normalised independent variables
1099   //     R2        Residuals versus dependent variable
1100   //     R3        Residuals computed on training sample
1101   //     R4        Residuals computed on test sample
1102   //
1103   // For a description of these quantities, refer to
1104   // Begin_Html<a href="#TMultiDimFet:description">class description</a>End_Html
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   // Counter variable
1115   Int_t i = 0;
1116 
1117   // Histogram of original variables
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   // Histogram of original dependent variable
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   // Histograms of normalized variables
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   // Histogram of shifted dependent variable
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   // Residual from training sample versus independent variables
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   // Residual from training sample versus. dependent variable
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   // Residual from training sample
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   // Generate the file <classname>MDF.cxx which contains the
1202   // implementation of the method:
1203   //
1204   //   Double_t <classname>::MDF(Double_t *x)
1205   //
1206   // which does the same as  TMultiDimFet::Eval. Please refer to this
1207   // method.
1208   //
1209   // Further, the public static members:
1210   //
1211   //   Int_t    <classname>::fgNVariables
1212   //   Int_t    <classname>::fgNCoefficients
1213   //   Double_t <classname>::fgDMean
1214   //   Double_t <classname>::fgXMean[]       //[fgNVariables]
1215   //   Double_t <classname>::fgXMin[]        //[fgNVariables]
1216   //   Double_t <classname>::fgXMax[]        //[fgNVariables]
1217   //   Double_t <classname>::fgCoefficient[] //[fgNCoeffficents]
1218   //   Int_t    <classname>::fgPower[]       //[fgNCoeffficents*fgNVariables]
1219   //
1220   // are initialized, and assumed to exist. The class declaration is
1221   // assumed to be in <classname>.h and assumed to be provided by the
1222   // user.
1223   //
1224   // See TMultiDimFet::MakeRealCode for a list of options
1225   //
1226   // The minimal class definition is:
1227   //
1228   //   class <classname> {
1229   //   public:
1230   //     Int_t    <classname>::fgNVariables;     // Number of variables
1231   //     Int_t    <classname>::fgNCoefficients;  // Number of terms
1232   //     Double_t <classname>::fgDMean;          // Mean from training sample
1233   //     Double_t <classname>::fgXMean[];        // Mean from training sample
1234   //     Double_t <classname>::fgXMin[];         // Min from training sample
1235   //     Double_t <classname>::fgXMax[];         // Max from training sample
1236   //     Double_t <classname>::fgCoefficient[];  // Coefficients
1237   //     Int_t    <classname>::fgPower[];        // Function powers
1238   //
1239   //     Double_t Eval(Double_t *x);
1240   //   };
1241   //
1242   // Whether the method <classname>::Eval should be static or not, is
1243   // up to the user.
1244 
1245   MakeRealCode(Form("%sMDF.cxx", classname), classname, option);
1246 }
1247 
1248 //____________________________________________________________________
1249 void TMultiDimFet::MakeNormalized() {
1250   // PRIVATE METHOD:
1251   // Normalize data to the interval [-1;1]. This is needed for the
1252   // classes method to work.
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       // Fill histograms of original independent variables
1273       if (TESTBIT(fHistogramMask, HIST_XORIG))
1274         ((TH1D *)fHistograms->FindObject(Form("x_%d_orig", j)))->Fill(fVariables(k));
1275 
1276       // Normalise independent variables
1277       fVariables(k) = 1 + 2 * range * (fVariables(k) - fMaxVariables(j));
1278 
1279       // Fill histograms of normalised independent variables
1280       if (TESTBIT(fHistogramMask, HIST_XNORM))
1281         ((TH1D *)fHistograms->FindObject(Form("x_%d_norm", j)))->Fill(fVariables(k));
1282     }
1283   }
1284   // Shift min and max of dependent variable
1285   fMaxQuantity -= fMeanQuantity;
1286   fMinQuantity -= fMeanQuantity;
1287 
1288   // Shift mean of independent variables
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   // PRIVATE METHOD:
1298   // Find the parameterization over the training sample. A full account
1299   // of the algorithm is given in the
1300   // Begin_Html<a href="#TMultiDimFet:description">class description</a>End_Html
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   // Loop over the number of functions we want to study.
1331   // increment inspection counter
1332   while (kTRUE) {
1333     // Reach user defined limit of studies
1334     if (studied++ >= fMaxStudy) {
1335       fParameterisationCode = PARAM_MAXSTUDY;
1336       break;
1337     }
1338 
1339     // Considered all functions several times
1340     if (k >= maxPass) {
1341       fParameterisationCode = PARAM_SEVERAL;
1342       break;
1343     }
1344 
1345     // increment function counter
1346     i++;
1347 
1348     // If we've reached the end of the functions, restart pass
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     // Print a happy message
1363     if (fIsVerbose && studied == 1)
1364       edm::LogInfo("TMultiDimFet") << "Coeff   SumSqRes    Contrib   Angle      QM   Func"
1365                                    << "     Value        W^2  Powers"
1366                                    << "\n";
1367 
1368     // Make the Gram-Schmidt
1369     Double_t dResidur = MakeGramSchmidt(i);
1370 
1371     if (dResidur == 0) {
1372       // This function is no good!
1373       // First test is in MakeGramSchmidt
1374       fFunctionCodes[i] = 1;
1375       continue;
1376     }
1377 
1378     // If user specified function, assume she/he knows what he's doing
1379     if (!fIsUserFunction) {
1380       // Flag this function as considered
1381       fFunctionCodes[i] = 2;
1382 
1383       // Test if this function contributes to the fit
1384       if (!TestFunction(squareResidual, dResidur)) {
1385         fFunctionCodes[i] = 1;
1386         continue;
1387       }
1388     }
1389 
1390     // If we get to here, the function currently considered is
1391     // fNCoefficients, so we increment the counter
1392     // Flag this function as OK, and store and the number in the
1393     // index.
1394     fFunctionCodes[i] = 3;
1395     fPowerIndex[fNCoefficients] = i;
1396     fNCoefficients++;
1397 
1398     // We add the current contribution to the sum of square of
1399     // residuals;
1400     squareResidual -= dResidur;
1401 
1402     // Calculate control parameter from this function
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     // Print the statistics about this function
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 /* && fIsVerbose */) {
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   // PRIVATE METHOD:
1443   // This is the method that actually generates the code for the
1444   // evaluation the parameterization on some point.
1445   // It's called by TMultiDimFet::MakeCode and TMultiDimFet::MakeMethod.
1446   //
1447   // The options are: NONE so far
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   // Write header of file
1464   //
1465   // Emacs mode line ;-)
1466   outFile << "// -*- mode: c++ -*-"
1467           << "\n";
1468   // Info about creator
1469   outFile << "// "
1470           << "\n"
1471           << "// File " << filename << " generated by TMultiDimFet::MakeRealCode"
1472           << "\n";
1473   // Time stamp
1474   TDatime date;
1475   outFile << "// on " << date.AsString() << "\n";
1476   // ROOT version info
1477   outFile << "// ROOT version " << gROOT->GetVersion() << "\n"
1478           << "//"
1479           << "\n";
1480   // General information on the code
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   // Header files
1501   if (isMethod)
1502     // If these are methods, we need the class header
1503     outFile << "#include \"" << classname << ".h\""
1504             << "\n";
1505 
1506   //
1507   // Now for the data
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   // Assignment to mean vector.
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   // Assignment to minimum vector.
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   // Assignment to maximum vector.
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   // Assignment to coefficients vector.
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   // Assignment to powers vector.
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   // Finally we reach the function itself
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   // EOF
1673   outFile << "// EOF for " << filename << "\n";
1674 
1675   // Close the file
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   // Print statistics etc.
1686   // Options are
1687   //   P        Parameters
1688   //   S        Statistics
1689   //   C        Coefficients
1690   //   R        Result of parameterisation
1691   //   F        Result of fit
1692   //   K        Correlation Matrix
1693   //   M        Pretty print formula
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     // Print basic parameters for this object
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     // Print statistics for read data
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       if (fFitter) {
1832          fFitter->PrintResults(1,1);
1833          edm::LogInfo("TMultiDimFet") << "\n";
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   //   M        Pretty print formula
1919   //
1920   Int_t i = 0;
1921   //   Int_t j = 0;
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   // Selection method. User can override this method for specialized
1982   // selection of acceptable functions in fit. Default is to select
1983   // all. This message is sent during the build-up of the function
1984   // candidates table once for each set of powers in
1985   // variables. Notice, that the argument array contains the powers
1986   // PLUS ONE. For example, to De select the function
1987   //     f = x1^2 * x2^4 * x3^5,
1988   // this method should return kFALSE if given the argument
1989   //     { 3, 4, 6 }
1990   return kTRUE;
1991 }
1992 
1993 //____________________________________________________________________
1994 void TMultiDimFet::SetMaxAngle(Double_t ang) {
1995   // Set the max angle (in degrees) between the initial data vector to
1996   // be fitted, and the new candidate function to be included in the
1997   // fit.  By default it is 0, which automatically chooses another
1998   // selection criteria. See also
1999   // Begin_Html<a href="#TMultiDimFet:description">class description</a>End_Html
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   // Set the min angle (in degrees) between a new candidate function
2011   // and the subspace spanned by the previously accepted
2012   // functions. See also
2013   // Begin_Html<a href="#TMultiDimFet:description">class description</a>End_Html
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   // Define a user function. The input array must be of the form
2025   // (p11, ..., p1N, ... ,pL1, ..., pLN)
2026   // Where N is the dimension of the data sample, L is the number of
2027   // terms (given in terms) and the first number, labels the term, the
2028   // second the variable.  More information is given in the
2029   // Begin_Html<a href="#TMultiDimFet:description">class description</a>End_Html
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   // Set the user parameter for the function selection. The bigger the
2045   // limit, the more functions are used. The meaning of this variable
2046   // is defined in the
2047   // Begin_Html<a href="#TMultiDimFet:description">class description</a>End_Html
2048   fPowerLimit = limit;
2049 }
2050 
2051 //____________________________________________________________________
2052 void TMultiDimFet::SetMaxPowers(const Int_t *powers) {
2053   // Set the maximum power to be considered in the fit for each
2054   // variable. See also
2055   // Begin_Html<a href="#TMultiDimFet:description">class description</a>End_Html
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   // Set the acceptable relative error for when sum of square
2066   // residuals is considered minimized. For a full account, refer to
2067   // the
2068   // Begin_Html<a href="#TMultiDimFet:description">class description</a>End_Html
2069   fMinRelativeError = error;
2070 }
2071 
2072 //____________________________________________________________________
2073 Bool_t TMultiDimFet::TestFunction(Double_t squareResidual, Double_t dResidur) {
2074   // PRIVATE METHOD:
2075   // Test whether the currently considered function contributes to the
2076   // fit. See also
2077   // Begin_Html<a href="#TMultiDimFet:description">class description</a>End_Html
2078 
2079   if (fNCoefficients != 0) {
2080     // Now for the second test:
2081     if (fMaxAngle == 0) {
2082       // If the user hasn't supplied a max angle do the test as,
2083       if (dResidur < squareResidual / (fMaxTerms - fNCoefficients + 1 + 1E-10)) {
2084         return kFALSE;
2085       }
2086     } else {
2087       // If the user has provided a max angle, test if the calculated
2088       // angle is less then the max angle.
2089       if (TMath::Sqrt(dResidur / fSumSqAvgQuantity) < TMath::Cos(fMaxAngle * DEGRAD)) {
2090         return kFALSE;
2091       }
2092     }
2093   }
2094   // If we get here, the function is OK
2095   return kTRUE;
2096 }
2097 
2098 //____________________________________________________________________
2099 //void mdfHelper(int& /*npar*/, double* /*divs*/, double& chi2,
2100 //               double* coeffs, int /*flag*/)
2101 /*{
2102    // Helper function for doing the minimisation of Chi2 using Minuit
2103 
2104    // Get pointer  to current TMultiDimFet object.
2105   // TMultiDimFet* mdf = TMultiDimFet::Instance();
2106    //chi2     = mdf->MakeChi2(coeffs);
2107 }
2108 */