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.9.0
/* @license-end */
THSIC.h
Go to the documentation of this file.
1
2// 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
65#ifndef THSIC_H
66#define THSIC_H
67
69#include "TVector.h"
70#include "TMatrixDSymfwd.h"
71
73#include "TSensitivityTest.h"
74
76#define NDefaultPermutation 500
77#define NDefaultThresholdPermutation 7
78#define NDefaultThresholdGamma 7
79#define NMinimumPatterns 6
80#define dEPSILONSigma 1.e-12
81
82namespace URANIE {
83namespace Sensitivity {
84class THSIC: public TSensitivityTest {
85public:
92 //---------------------------------------------
100
109 THSIC(URANIE::DataServer::TDataServer *tds, const char *sInputs, const char *sOutputs, Option_t * option = "");
110
114 virtual ~THSIC();
116
117 //---------------------------------------------
129 virtual void generateSample(Option_t * option = "");
130
135 void computeIndexes(Option_t * option = "");
140 void evaluateIndexes(Option_t * option = "");
141
149 void drawIndexes(TString sTitle, const char *select = "", Option_t * option = "");
151
152
154
167
175
181 return _nPermutationSample;
182 }
183
188 void setNPermutationSample(Int_t ns);
189
201
207 void setThresholdGamma(Int_t n);
208
215 return _nThresholdGamma;
216 }
218
219 //---------------------------------------------
231 void generateHMatrix(Option_t *option = "");
232
240 Double_t computeTrace(TMatrixD mat, TMatrixD crtMatHLH);
241
248 TMatrixDSym doubleCentering(const TMatrixDSym& mat);
249
256 TMatrixDSym getMatrixKGaussian(const char *sAtt, Double_t dvariance);
257
264 Double_t getSumOfMatrixBk(TMatrixD matHLkH, TMatrixD matHLH);
265
272 Double_t estimateSigma(const char *sAtt, Option_t *option = "");
273
288 Double_t unBiasedHSIC( TMatrixD matKx, TMatrixD matKy);
289
297 void permutedPermMatrix(TMatrixDSym& theMatK, TList *tlstIndex, TList *tlstMatrix);
299
300 //---------------------------------------------
314 virtual void parseOption(Option_t *option = "");
316
317 //---------------------------------------------
326 virtual void printLog(Option_t *option = "");
327
335 void printSummary(Option_t * option = "");
337
338 ClassDef(URANIE::Sensitivity::THSIC, ID_SENSITIVITY)
339
340 //---------------------------------------------
345public:
346 TMatrixDSym _matH;
347private:
353
354};
355
356} // End Of namespace Sensitivity
357} // End Of namespace URANIE
358#endif
Interface of class URANIE::Sensitivity::TSensitivityTest.
Description of the class THSIC.
Definition THSIC.h:84
void setThresholdPermutation(Int_t n)
Int_t _nThresholdPermutation
The number of permutation sample ( permutation ) (default NDefaultPermutation)
Definition THSIC.h:350
void generateHMatrix(Option_t *option="")
generate the matrix
void printSummary(Option_t *option="")
Print the summary of the HSIC method.
THSIC()
Default constructor.
Int_t getNPermutationSample()
getNPermutationSample
Definition THSIC.h:180
void setEstimatedVariance(EEstimatedVariance nMethod)
Set the method used to estimate the variance of the gaussian kernel ( kUnknown | kStdEmpirical | kMed...
virtual void generateSample(Option_t *option="")
The generate Sample method.
Int_t getThresholdGamma()
Get the Threshold (nS/nX) to use the Gamma approximation to cumpute the p-value (default NDefaultThre...
Definition THSIC.h:214
TMatrixDSym doubleCentering(const TMatrixDSym &mat)
compute the matrix matrix
virtual void parseOption(Option_t *option="")
Parse the option given by the user.
void evaluateIndexes(Option_t *option="")
evaluateIndexes
void drawIndexes(TString sTitle, const char *select="", Option_t *option="")
Draw the indexes computed by the method.
Double_t unBiasedHSIC(TMatrixD matKx, TMatrixD matKy)
Use the unbiased formula to compute .
void setNPermutationSample(Int_t ns)
Set the Permutation sample permutation Size (default NDefaultPermutation)
Int_t getThresholdPermutation()
Get the Threshold (nS/nX) to use the permutation permutation test to cumpute the p-value (default NDe...
Definition THSIC.h:198
TMatrixDSym getMatrixKGaussian(const char *sAtt, Double_t dvariance)
getMatrixKGaussian
Int_t _nPermutationSample
the method to compute the variance for the gaussian kernel
Definition THSIC.h:349
virtual void printLog(Option_t *option="")
Print the Log.
void permutedPermMatrix(TMatrixDSym &theMatK, TList *tlstIndex, TList *tlstMatrix)
Create the list of permuted Permutation sample.
virtual ~THSIC()
Default destructor.
Double_t computeTrace(TMatrixD mat, TMatrixD crtMatHLH)
compute the trace of the product matrix
Int_t _nThresholdGamma
The Threshold (nS/nX) to use the permutation permutation test to cumpute the p-value (default 5)
Definition THSIC.h:351
Double_t getSumOfMatrixBk(TMatrixD matHLkH, TMatrixD matHLH)
getSumOfMatrixBk
void setThresholdGamma(Int_t n)
Set the Threshold (nS/nX) to use the Gamma approximation to cumpute the p-value (default NDefaultThre...
EEstimatedVariance getEstimatedVariance()
Which method used to estimate the variance of the gaussian kernel ( kUnknown | kStdEmpirical | kMedia...
Definition THSIC.h:164
THSIC(URANIE::DataServer::TDataServer *tds, const char *sInputs, const char *sOutputs, Option_t *option="")
THSIC constructor with data contined in the TDS and the inputs and outputs attributes.
Double_t estimateSigma(const char *sAtt, Option_t *option="")
estimateSigma
EEstimatedVariance
Definition THSIC.h:87
@ kMedianeDelta
Definition THSIC.h:90
@ kStdEmpirical
Definition THSIC.h:89
@ kUnknown
Definition THSIC.h:88
void computeIndexes(Option_t *option="")
computeIndexes
EEstimatedVariance _nEstimatedVariance
< The H matrix defined by H := ( Kronecker(i,j) - 1/nS)
Definition THSIC.h:348
Description of the class TSensitivityTest.
Definition TSensitivityTest.h:53
ROOT.
Definition TCMN.h:45