Documentation / Manuel développeur
Modules disponibles
Calibration,  DataServer,  Launcher,  MetaModelOptim,  Modeler,  Optimizer,  ReLauncher,  Reliability,  ReOptimizer,  Sampler,  Sensitivity,  UncertModeler,  XmlProblem,  
Uranie / Sensitivity v4.9.0
|
Description of the class TMorris. More...
#include <TMorris.h>
Public Member Functions | |
Constructor and Destructor | |
TMorris (URANIE::DataServer::TDataServer *tds, const char *fcn, Int_t nReplique, Int_t level, Double_t delta=0, TString sout="") | |
Default constructor with the name of a function. | |
TMorris (URANIE::DataServer::TDataServer *tds, void(*fcn)(Double_t *, Double_t *), TString sin, TString sout, Int_t nReplique, Int_t level, Double_t delta=0) | |
Default constructor with a function. | |
TMorris (URANIE::DataServer::TDataServer *tds, const char *varexpinput, const char *varexpoutput, Option_t *option="") | |
Default constructor with a TDataServer filling. | |
TMorris (URANIE::DataServer::TDataServer *tds, URANIE::Launcher::TCode *fcode, Int_t nReplique, Int_t level, Double_t delta=0) | |
Default constructor with a tcode. | |
TMorris (URANIE::DataServer::TDataServer *tds, URANIE::Relauncher::TRun *run, Int_t nReplique, Int_t level, Double_t delta=0) | |
Default constructor with a trun. | |
~TMorris () | |
Default destructor. | |
Setting and Getting attributs | |
URANIE::DataServer::TDSNtupleD * | getMorrisResults () |
Get Morris results. | |
TTree * | getResultTuple (bool commonresulttuple=false) |
Get the result ntuple (default parameter unused but for Morris method) | |
void | setLevel (Int_t val) |
Sets the level. | |
Int_t | getLevel () |
Gets the level. | |
void | setDelta (Double_t val) |
Sets the delta. | |
Double_t | getDelta () |
Gets the delta. | |
Int_t | getNReplica () |
Gets the number of replica. | |
Generation of the sample | |
virtual void | generateSample (Option_t *option="") |
void | computeOrientationMatrixOld () |
Computes the orientation matrix . | |
void | computeOrientationMatrix () |
Evaluate the indexes | |
void | evaluateIndexes (Option_t *option="") |
Evaluates the index from a Specific TDataServer. | |
void | preTreatment () |
PreTreatment for every output. | |
visualization and Scan the data | |
virtual void | createTuple () |
virtual void | drawIndexes (TString sTitre, const char *select="", Option_t *option="", bool internalCall=kFALSE) |
Draws the indexes. | |
void | drawSample (TString svar="", Int_t nreplique=0, Option_t *option="") |
void | drawIndexes (Option_t *option="", bool internalCall=kFALSE) |
void | drawIndexesAll (Option_t *option="") |
void | pyDrawIndexes (Option_t *option="") |
Printing Log | |
void | printLog (Option_t *option="") |
Public Member Functions inherited from URANIE::Sensitivity::TScreeming | |
TScreeming (URANIE::DataServer::TDataServer *tds, const char *fcn, Int_t ns=-1, const char *varexpinput="", const char *varexpoutput="") | |
Default constructor with the name of a function. | |
TScreeming (URANIE::DataServer::TDataServer *tds, void(*fcn)(Double_t *, Double_t *), const char *varexpinput, const char *varexpoutput, Int_t ns=-1) | |
TScreeming (URANIE::DataServer::TDataServer *tds, const char *varexpinput, const char *varexpoutput, Option_t *option="") | |
Default constructor with a TDataServer filling. | |
TScreeming (URANIE::DataServer::TDataServer *tds, URANIE::Launcher::TCode *code, Int_t ns=100) | |
Default constructor with TCode arg. | |
TScreeming (URANIE::DataServer::TDataServer *tds, URANIE::Relauncher::TRun *code, Int_t ns=100) | |
Default constructor with TRun arg. | |
virtual | ~TScreeming () |
Default destructor. | |
Public Member Functions inherited from URANIE::Sensitivity::TSensitivity | |
TSensitivity () | |
Default constructor. | |
TSensitivity (URANIE::DataServer::TDataServer *tds, const char *fcn, Int_t ns, const char *varexpinput="", const char *varexpoutput="", Option_t *option="") | |
Default constructor with the name of a function. | |
TSensitivity (URANIE::DataServer::TDataServer *tds, void(*fcn)(Double_t *, Double_t *), const char *varexpinput, const char *varexpoutput, Int_t ns, Option_t *option="") | |
TSensitivity (URANIE::DataServer::TDataServer *tds, const char *varexpinput, const char *varexpoutput, Option_t *option="") | |
Default constructor. | |
TSensitivity (URANIE::DataServer::TDataServer *tds, URANIE::Launcher::TCode *code, Int_t ns, Option_t *option="") | |
Default constructor with TCode arg. | |
TSensitivity (URANIE::DataServer::TDataServer *tds, URANIE::Relauncher::TRun *run, Int_t ns, Option_t *option="") | |
Default constructor with TRun arg. | |
virtual | ~TSensitivity () |
Default destructor. | |
Int_t | getID () |
double | getValue (const char *sorder="", const char *sinputname="", const char *sselect="") |
vector< int > * | getAttributeElements (string str) |
void | setFunction (void(*fct)(Double_t *, Double_t *), Int_t nx=-1, Int_t ny=1) |
TString | getFunctionName () |
void | setSeed (UInt_t nval) |
UInt_t | getSeed () |
virtual void | setMethodName (const char *str) |
Set the Method name. | |
const char * | getMethodName () |
Get the method name. | |
Bool_t | getNoIntermediateSaved () |
Get the noIntermediateSaved flag. | |
const char * | getIteratorName () |
Get the name of the iterator attribut of the method. | |
void | setSensitivityIteratorName (const char *str) |
Set the iterator name devoted to compute sensitivity indexes. | |
void | setTimeName (TString sname) |
Set the name of the time attribute (only one) | |
TString | getTimeName () |
Get the name of the time attribute. | |
virtual void | setDrawProgressBar (Bool_t bbool=kTRUE) |
Set the "draw progress bar" flag. | |
void | setUsingErrors (bool thebool) |
Set the "using error results anyway" option. | |
Bool_t | getDrawProgressBar () |
Get the clean flag. | |
Bool_t | isInputCorrelated () |
TMatrixD | getMatrixInputCorrelation () |
Int_t | getNInput () |
Get the number of input attributes. | |
Int_t | getNOutput () |
Get the number of output attributes. | |
void | setInputCorrelationMatrix (TMatrixD Corr) |
virtual void | parseOption (Option_t *option="") |
Read the possible options. | |
void | checkOutputRequested (string attlist, bool fromoption=false) |
Check the output list requested by the user. | |
void | computeIndexes (Option_t *option="") |
Compute the Sensitivity Indexes. | |
void | fillIndex (const char *sinputname, const char *sorder, Double_t dval, const char *algo="", Double_t dvalCILower=-1.0, Double_t dvalCIUpper=-1.0) |
Method to fill in the tree the value of Sensitivity indexes for an input attribute. | |
virtual void | createTuple (Option_t *option="") |
virtual void | drawIndexes (TString sTitre, const char *select="", Option_t *option="") |
Draws the indexes. | |
virtual void | setLog () |
void | unsetLog () |
void | changeLog () |
Bool_t | getLog () |
void | setNLauncher (ELauncher codeLauncher) |
Private Member Functions | |
void | stdInit () |
constructor helpers | |
Private Attributes | |
Double_t | _delta |
The Delta parameter. | |
Int_t | _level |
The level parameter (known as p) | |
TMatrixD | _matB |
TVectorD | _vecDirection |
Matrix of . | |
TNtupleD * | _ntdatamorris |
Vector of direction. | |
URANIE::DataServer::TDSNtupleD * | _ntmorrisresult |
The tntuple of result . | |
TRandom2 * | _rand |
Additional Inherited Members | |
Public Types inherited from URANIE::Sensitivity::TSensitivity | |
enum | ELauncher { kCode , kCodeRemote , kFunction , kRun , kUnknown } |
Public Attributes inherited from URANIE::Sensitivity::TSensitivity | |
URANIE::DataServer::TDataServer * | _tds |
Pointeur vers un TDS. | |
TList * | _listOfInputAttributes |
List of the input branches. | |
TList * | _listOfOutputAttributes |
List of the input branches. | |
TString | _sTimeAttribute |
The name of the Time attribute. | |
Int_t | _nS |
The number of simulation or other information depend on the used method. | |
Int_t | _nX |
Dimension of the input. | |
UInt_t | _nY |
Dimension of the target. | |
UInt_t | _nElY |
Number of element for one selected output. | |
Int_t | _nbOut |
Total number of Output to be considered. | |
Int_t | _iOut |
counter for output | |
unsigned int | _iy |
iterator over number of element | |
unsigned int | _iely |
iterator over number of element | |
ELauncher | _nLauncher |
The type of launcher. | |
TString | _sFunctionName |
The Name of the evaluatuor. | |
URANIE::Launcher::TCode * | _code |
The tcode. | |
URANIE::Relauncher::TRun * | _run |
TObjArray * | _drawingGarbageCollector |
Garbage collector for prints. | |
Int_t | _nSeed |
The seed of the random generator. | |
Bool_t | _bChosenOutputs |
Fact that the input list is provided or not. | |
Bool_t | _blog |
Boolean for edit the log. | |
Bool_t | _bdrawProgressBar |
Boolean to know if the progress bar has to be drawn. | |
Bool_t | _bnoIntermediateSaved |
Boolean to know if the progress bar has to be drawn. | |
TString | _sIteratorName |
The specific iterator attribute for the method. | |
TString | _sMethodName |
The method name. | |
TString | _sSelectedOutput |
The output. | |
TString | _sSelectedInput |
The input. | |
map< string, unsigned int > | _mAttributeSize |
Map of size of element for attribute;. | |
map< string, vector< int > > | _mAttributeElements |
Map of Elements number to run (if vector subselection is requested) | |
vector< string > | _vOutputNames |
Name of the output. | |
TCanvas * | _canvas |
Canvas object to deal with. | |
void(* | _pFunction )(double *, double *) |
TTree * | _ntresult |
The TTree of results. | |
Protected Member Functions inherited from URANIE::Sensitivity::TSensitivity | |
void | checkCanvasCreation (bool newcan) |
Create a canvas if needed. | |
void | drawIndexesHistogram (TString sTitre, const char *select="", Option_t *option="") |
Draws indexes with an histogram. | |
void | drawIndexesPie (TString sTitre, const char *select="", Option_t *option="") |
Draws indexes with an pie chart. | |
virtual void | postTreatment () |
PostTreatment for every output. | |
void | setNoIntermediateSaved (Bool_t bbool=kTRUE) |
Set the "only final file" flag. | |
Protected Attributes inherited from URANIE::Sensitivity::TSensitivity | |
Char_t | _sOutputAttribute [MAXLENGTHSTRING] |
The name of the output attribute. | |
Char_t | _sInputAttribute [MAXLENGTHSTRING] |
The name of the input attribute. | |
Char_t | _sOrder [MAXLENGTHSTRING] |
The order of sensitivity indexes. | |
Char_t | _sMethod [MAXLENGTHSTRING] |
The name of the method. | |
Char_t | _sAlgorithm [MAXLENGTHSTRING] |
The name of the algorithm to compute the index. | |
Double_t | _valSobolCrt |
The value of sensitivity indexes. | |
Double_t | _valSobolCILower |
The value of lower Condidence Interval (95) | |
Double_t | _valSobolCIUpper |
The value of upper Condidence Interval (95) | |
TMatrixD | _minputCorr |
Input correlation matrix if sample needs to be correlated. | |
Bool_t | _bisInputCorrelated |
State whether the input correlation matrix is set. | |
Bool_t | _bgoingThroughError |
State whether the error must not block the computation. | |
Detailed Description
Description of the class TMorris.
The method recommended by Morris (1991) is valid for p even, and the value for is given by the relation
where is the level.
The number for the code evaluation is given by the relation
Using TMorris class inside a script
This is the use of TMorris class. For a generic documentation, please refer to the section Using the TSensitivity class.
TMorris * scmo = new TMorris(tds, function, replique, level); TMorris * scmo = new TMorris(tds, TCode, replique, level); scmo->generateSample(); scmo->computeIndexes(); scmo->drawIndexes("ddpie", "", "nonewcanv,pie,first"); ntsampl->Print();
Constructor & Destructor Documentation
◆ TMorris() [1/5]
URANIE::Sensitivity::TMorris::TMorris | ( | URANIE::DataServer::TDataServer * | tds, |
const char * | fcn, | ||
Int_t | nReplique, | ||
Int_t | level, | ||
Double_t | delta = 0 , |
||
TString | sout = "" |
||
) |
Default constructor with the name of a function.
Referenced by ClassImp().
◆ TMorris() [2/5]
URANIE::Sensitivity::TMorris::TMorris | ( | URANIE::DataServer::TDataServer * | tds, |
void(*)(Double_t *, Double_t *) | fcn, | ||
TString | sin, | ||
TString | sout, | ||
Int_t | nReplique, | ||
Int_t | level, | ||
Double_t | delta = 0 |
||
) |
Default constructor with a function.
◆ TMorris() [3/5]
URANIE::Sensitivity::TMorris::TMorris | ( | URANIE::DataServer::TDataServer * | tds, |
const char * | varexpinput, | ||
const char * | varexpoutput, | ||
Option_t * | option = "" |
||
) |
Default constructor with a TDataServer filling.
◆ TMorris() [4/5]
URANIE::Sensitivity::TMorris::TMorris | ( | URANIE::DataServer::TDataServer * | tds, |
URANIE::Launcher::TCode * | fcode, | ||
Int_t | nReplique, | ||
Int_t | level, | ||
Double_t | delta = 0 |
||
) |
Default constructor with a tcode.
◆ TMorris() [5/5]
URANIE::Sensitivity::TMorris::TMorris | ( | URANIE::DataServer::TDataServer * | tds, |
URANIE::Relauncher::TRun * | run, | ||
Int_t | nReplique, | ||
Int_t | level, | ||
Double_t | delta = 0 |
||
) |
Default constructor with a trun.
◆ ~TMorris()
URANIE::Sensitivity::TMorris::~TMorris | ( | ) |
Default destructor.
Referenced by ClassImp().
Member Function Documentation
◆ computeOrientationMatrix()
void URANIE::Sensitivity::TMorris::computeOrientationMatrix | ( | ) |
Referenced by ClassImp().
◆ computeOrientationMatrixOld()
void URANIE::Sensitivity::TMorris::computeOrientationMatrixOld | ( | ) |
Computes the orientation matrix .
algorithme
The matrix d�fined in the paper of Morris is computed
- Warning
- the author underlines the fact that the column of the matrix colonne is either the matrix of B, or the column of B by replacing the 0 by 1 or reciprocally the 1 by 0. Thus, we can even compute .
Referenced by ClassImp().
◆ createTuple()
|
virtual |
Referenced by ClassImp().
◆ drawIndexes() [1/2]
void URANIE::Sensitivity::TMorris::drawIndexes | ( | Option_t * | option = "" , |
bool | internalCall = kFALSE |
||
) |
◆ drawIndexes() [2/2]
|
virtual |
Draws the indexes.
Referenced by ClassImp().
◆ drawIndexesAll()
void URANIE::Sensitivity::TMorris::drawIndexesAll | ( | Option_t * | option = "" | ) |
Referenced by ClassImp().
◆ drawSample()
void URANIE::Sensitivity::TMorris::drawSample | ( | TString | svar = "" , |
Int_t | nreplique = 0 , |
||
Option_t * | option = "" |
||
) |
Referenced by ClassImp().
◆ evaluateIndexes()
|
virtual |
Evaluates the index from a Specific TDataServer.
The TDS must contain an attribute Name + "n__iter__sobol" with values
- -1 : for the M matrix
- 0 : for the N matrix
- ix+1 : for the matrix Mi
- 100*_nX+ ix+1 : for the matrix Ni
- Parameters
-
option (Option_t *) option to pass [""]
Implements URANIE::Sensitivity::TSensitivity.
Referenced by ClassImp().
◆ generateSample()
|
virtual |
Implements URANIE::Sensitivity::TSensitivity.
Referenced by ClassImp().
◆ getDelta()
|
inline |
Gets the delta.
References _delta.
◆ getLevel()
|
inline |
Gets the level.
References _level.
◆ getMorrisResults()
|
inline |
Get Morris results.
References _ntmorrisresult.
◆ getNReplica()
|
inline |
Gets the number of replica.
References URANIE::Sensitivity::TSensitivity::_nS.
◆ getResultTuple()
|
virtual |
Get the result ntuple (default parameter unused but for Morris method)
Reimplemented from URANIE::Sensitivity::TSensitivity.
Referenced by ClassImp().
◆ preTreatment()
|
virtual |
PreTreatment for every output.
Doing a pre-treatment when running over outputs.
Reimplemented from URANIE::Sensitivity::TSensitivity.
Referenced by ClassImp().
◆ printLog()
|
virtual |
Reimplemented from URANIE::Sensitivity::TScreeming.
Referenced by ClassImp().
◆ pyDrawIndexes()
void URANIE::Sensitivity::TMorris::pyDrawIndexes | ( | Option_t * | option = "" | ) |
Referenced by ClassImp().
◆ setDelta()
|
inline |
Sets the delta.
References _delta.
◆ setLevel()
|
inline |
Sets the level.
References _level.
◆ stdInit()
|
private |
constructor helpers
Referenced by ClassImp().
Member Data Documentation
◆ _delta
|
private |
The Delta parameter.
Referenced by ClassImp(), getDelta(), and setDelta().
◆ _level
|
private |
The level parameter (known as p)
Referenced by ClassImp(), getLevel(), and setLevel().
◆ _matB
|
private |
Referenced by ClassImp().
◆ _ntdatamorris
|
private |
◆ _ntmorrisresult
|
private |
The tntuple of result .
Referenced by ClassImp(), and getMorrisResults().
◆ _rand
|
private |
Referenced by ClassImp().
◆ _vecDirection
|
private |
Matrix of .
Referenced by ClassImp().