English Français

Documentation / Developer's manual

Available modules

Calibration,  DataServer,  Launcher,  MetaModelOptim,  Modeler,  Optimizer,  ReLauncher,  Reliability,  ReOptimizer,  Sampler,  Sensitivity,  UncertModeler,  XmlProblem,   Uranie / Sensitivity: THSIC.h Source File
Uranie / Sensitivity  v4.10.0
/* @license-end */
THSIC.h
Go to the documentation of this file.
1 // Copyright (C) 2013-2024 CEA/DES
3 //
4 // This program is free software: you can redistribute it and/or modify
5 // it under the terms of the GNU Lesser General Public License as published
6 // by the Free Software Foundation, either version 3 of the License, or any
7 // later version.
8 //
9 // This program is distributed in the hope that it will be useful,
10 // but WITHOUT ANY WARRANTY; without even the implied warranty of
11 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
12 // GNU Lesser General Public License for more details.
13 //
14 // You should have received a copy of the GNU Lesser General Public License
15 // along with this program. If not, see <http://www.gnu.org/licenses/>.
18 // $Id$
19 // $Author$
20 // $Date$
21 // $Revision$
22 // $State$
24 
66 #ifndef THSIC_H
67 #define THSIC_H
68 
70 #include "TVector.h"
71 #include "TMatrixDSymfwd.h"
72 
74 #include "TSensitivityTest.h"
75 
77 #define NDefaultPermutation 500
78 #define NDefaultThresholdPermutation 7
79 #define NDefaultThresholdGamma 7
80 #define NMinimumPatterns 6
81 #define dEPSILONSigma 1.e-12
82 
83 namespace URANIE {
84 namespace Sensitivity {
85 class THSIC: public TSensitivityTest {
86 public:
87 
89  {
93  };
94 
95  //---------------------------------------------
99 
105  static TString convertRule(EEstimatedVariance evRule);
106 
108 
109 
110 
111  //---------------------------------------------
115 
119  THSIC();
120 
129  THSIC(URANIE::DataServer::TDataServer* tds, const char* sInputs, const char* sOutputs, Option_t* option="");
130 
134  virtual ~THSIC();
135 
137 
138 
139 
140  //---------------------------------------------
144 
153  return _nEstimatedVariance;
154  }
155 
163 
169  return _nPermutationSample;
170  }
171 
176  void setNPermutationSample(Int_t ns);
177 
180  void setThresholdPermutation(Int_t n);
187  return _nThresholdPermutation;
188  }
189 
195  void setThresholdGamma(Int_t n);
196 
203  return _nThresholdGamma;
204  }
205 
207 
208 
209 
210  //---------------------------------------------
214 
228  virtual void parseOption(Option_t *option = "");
229 
231 
232 
233 
234  //---------------------------------------------
238 
247  void printSummary();
248 
254  virtual void printLog(Option_t* option="");
255 
257 
258 
259 
260  //---------------------------------------------
265 
273  void generateHMatrix(Option_t *option = "");
274 
286  TMatrixD doubleCenteringRect(const TMatrixD& matK);
287 
299  TMatrixDSym doubleCentering(const TMatrixDSym& matK);
300 
314  Double_t estimateSigma(const char* sAtt, Option_t* option="");
315 
327  TMatrixDSym getMatrixKGaussian(const char* sAtt, Double_t dTheta);
328 
348  Double_t biasedHSIC(const TMatrixDSym& matKx, const TMatrixDSym& matKy);
349 
373  Double_t unbiasedHSIC(const TMatrixDSym& matKx, const TMatrixDSym& matKy);
374 
392  Double_t getProductOfNorms(const TMatrixDSym& matKx, const TMatrixDSym& matKy);
393 
411  void permuteKMatrix(const TMatrixDSym& theMatK, TList* tlstPerm, TList* tlstPermMat);
412 
414 
415 
416 
417  //---------------------------------------------
421 
429  virtual void generateSample(Option_t* option="");
430 
436  void computeIndexes(Option_t* option="");
437 
463  void evaluateIndexes(Option_t* option);
464 
472  void drawIndexes(TString sTitle, const char* select="", Option_t* option="");
473 
475 
476 
477 
478  ClassDef(URANIE::Sensitivity::THSIC, ID_SENSITIVITY)
479 
480  //---------------------------------------------
485 
486  public:
487  TMatrixDSym _matH;
488  private:
493 
495 
496 
497 
498 };
499 
500 } // End Of namespace Sensitivity
501 } // End Of namespace URANIE
502 
503 #endif
ROOT.
Definition: TCMN.h:45
virtual ~THSIC()
Default destructor.
TMatrixD doubleCenteringRect(const TMatrixD &matK)
Centers the rows and columns of the provided matrix.
void setThresholdGamma(Int_t n)
Set the Threshold (nS/nX) to use the Gamma approximation to cumpute the p-value (default NDefaultThre...
Double_t estimateSigma(const char *sAtt, Option_t *option="")
Computes the bandwidth parameter of the specified attribute.
Int_t getThresholdPermutation()
Get the Threshold (nS/nX) to use the permutation permutation test to compute the p-value (default NDe...
Definition: THSIC.h:186
virtual void generateSample(Option_t *option="")
No action is performed, but a warning is returned.
Int_t _nPermutationSample
the method to compute the variance for the gaussian kernel
Definition: THSIC.h:490
TMatrixDSym doubleCentering(const TMatrixDSym &matK)
Centers the rows and columns of the provided matrix.
void setEstimatedVariance(EEstimatedVariance nMethod)
Set the method used to estimate the variance of the gaussian kernel ( kUnknown | kStdEmpirical | kMed...
TMatrixDSym getMatrixKGaussian(const char *sAtt, Double_t dTheta)
Constructs the Gram matrix associated with the specified attribute.
static TString convertRule(EEstimatedVariance evRule)
Outputs the string representation of any enumerator of EEstimatedVariance.
THSIC()
Default constructor.
EEstimatedVariance
Definition: THSIC.h:88
void permuteKMatrix(const TMatrixDSym &theMatK, TList *tlstPerm, TList *tlstPermMat)
Given a Gram matrix and a list of permutations, constructs the permuted matrices and store them in a ...
Int_t getNPermutationSample()
getNPermutationSample
Definition: THSIC.h:168
Int_t _nThresholdPermutation
The number of permutation sample ( permutation ) (default NDefaultPermutation)
Definition: THSIC.h:491
void evaluateIndexes(Option_t *option)
Applies the HSIC methodology.
virtual void parseOption(Option_t *option="")
Parse the option given by the user.
void printSummary()
Prints a concise summary of results to the terminal.
void computeIndexes(Option_t *option="")
Calls the method evaluateIndexes.
EEstimatedVariance getEstimatedVariance()
Which method used to estimate the variance of the gaussian kernel ( kUnknown | kStdEmpirical | kMedia...
Definition: THSIC.h:152
Int_t _nThresholdGamma
The Threshold (nS/nX) to use the permutation permutation test to cumpute the p-value (default 5) ...
Definition: THSIC.h:492
Description of the class TSensitivityTest.
Definition: TSensitivityTest.h:53
EEstimatedVariance _nEstimatedVariance
< The H matrix defined by H := ( Kronecker(i,j) - 1/nS)
Definition: THSIC.h:489
Interface of class URANIE::Sensitivity::TSensitivityTest.
void setThresholdPermutation(Int_t n)
Int_t getThresholdGamma()
Get the Threshold (nS/nX) to use the Gamma approximation to cumpute the p-value (default NDefaultThre...
Definition: THSIC.h:202
void setNPermutationSample(Int_t ns)
Set the Permutation sample permutation Size (default NDefaultPermutation)
void generateHMatrix(Option_t *option="")
generate the matrix
void drawIndexes(TString sTitle, const char *select="", Option_t *option="")
Draw the indexes computed by the method.
virtual void printLog(Option_t *option="")
Prints all attributes of the class for the current instance.
Double_t unbiasedHSIC(const TMatrixDSym &matKx, const TMatrixDSym &matKy)
Computes the U-statistic estimate of . The computation is based on the Gram matrices and which are ...
Double_t biasedHSIC(const TMatrixDSym &matKx, const TMatrixDSym &matKy)
Computes the V-statistic estimate of where . The computation is based on the Gram matrices and whi...
Double_t getProductOfNorms(const TMatrixDSym &matKx, const TMatrixDSym &matKy)
Estimates the product of the squared Hilbert-Schmidt (HS) norms of the cross-covariance operators an...
Description of the class THSIC.
Definition: THSIC.h:85