English Français

Documentation / Developer's manual

Available modules

Calibration,  DataServer,  Launcher,  MetaModelOptim,  Modeler,  Optimizer,  ReLauncher,  Reliability,  ReOptimizer,  Sampler,  Sensitivity,  UncertModeler,  XmlProblem,   Uranie / Sampler: TSampling.h Source File
Uranie / Sampler v4.9.0
/* @license-end */
TSampling.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
32#ifndef TSAMPLING_H
33#define TSAMPLING_H
34
35#ifndef ROOT_TMatrixD
36#include "TMatrixD.h"
37#endif
38#ifndef ROOT_TVectorD
39#include "TVectorD.h"
40#endif
41#ifndef ROOT_TNtuple
42#include "TNtuple.h"
43#endif
44#ifndef ROOT_TString
45#include "TString.h"
46#endif
47#ifndef ROOT_TNtupleD
48#include "TNtupleD.h"
49#endif
50#ifndef ROOT_TNamed
51#include "TNamed.h"
52#endif
53
54// Uranie
55#include "TSamplerStochastic.h"
56#include "TAttribute.h"
57#include "TStochAttribut.h"
58#include "TDataServer.h"
59#include "LHSInput.H"
60
61namespace URANIE
62{
63namespace Sampler
64{
66{
67 // Associations
68
69 // Attributes
70public:
71 enum EType
72 {
74 };
75 LHSInput* _fLhsInput;
76 TMatrixD _corrmatrix;
79 bool _buseSVD;
81
82 // Operations
83private:
84 //---------------------------------------------
88
94 void computeCorrelationMatrix(TMatrixD matOrg);
96
97public:
98 //---------------------------------------------
102
114 TSampling(URANIE::DataServer::TDataServer *tds, Option_t *option = "lhs",
115 Int_t nCalcul = 1000);
117 virtual ~TSampling();
119
120 //---------------------------------------------
129
138 void generateSample(Option_t *option = "");
140
144
145 //---------------------------------------------
150
158 void setUserCorrelation(Int_t indx, Int_t indy, double value);
159
161
167 void setUserCorrelation(TString xname, TString yname, double value);
168
170
177 void setUserCorrelation(URANIE::DataServer::TAttribute *x,
178 URANIE::DataServer::TAttribute *y, Double_t value);
179
181
186 void setCorrelationMatrix(TMatrixD corrMat);
187
189
191
195 void setSaveRang(Int_t rang = 0);
196 void setTypeSampling(Option_t *option);
197
198 //---------------------------------------------
202
203 virtual void printLog(Option_t *option = "");
205
206 ClassDef(URANIE::Sampler::TSampling, ID_SAMPLER)
207 //Definition d'un plan d'exp�rience
208};
209
210} // Fin du namespace Sampler
211} // Fin du namespace URANIE
212#endif
213// fin du fichier $RCSfile$.
Creation of the abstract class TSamplerStochastic.
Definition TSamplerStochastic.h:44
Definition TSampling.h:66
void CorrelationFactorization()
Decomposes the correlation matrix.
void generateSample(Option_t *option="")
Generates the sample.
TMatrixD _corrmatrix
The correlation matrix.
Definition TSampling.h:76
EType _ntype
The type of sampling (SRS, LHS)
Definition TSampling.h:78
void generateSampleByDakota()
Generates the sample by the Dakota module.
TMatrixD _corrIndexPosition
The permutation for giving the rank correlation matrix given by the user.
Definition TSampling.h:77
bool _buseSVD
useSVD instead of Cholesky to factorize the input correlation matrix
Definition TSampling.h:79
void setTypeSampling(Option_t *option)
EType
Definition TSampling.h:72
@ kLHS
Definition TSampling.h:73
@ kUnknown
Definition TSampling.h:73
@ kSRS
Definition TSampling.h:73
TSampling(URANIE::DataServer::TDataServer *tds, Option_t *option="lhs", Int_t nCalcul=1000)
Constructor from a TDataServer, the name of the method and the size of the sample.
void setCorrelationMatrix(TMatrixD corrMat)
Set the correlation matrix.
void setUserCorrelation(TString xname, TString yname, double value)
Defines a correlation between two attributes given by their names.
void computeCorrelationMatrix(TMatrixD matOrg)
Computes the correlation matrix for a TMatrixD.
void setSaveRang(Int_t rang=0)
EMPTY:: TSampling::setSaveRang.
LHSInput * _fLhsInput
Definition TSampling.h:75
void setUserCorrelation(Int_t indx, Int_t indy, double value)
Defines a correlation between two attributes given by their indexes.
virtual ~TSampling()
Default destructor.
virtual void printLog(Option_t *option="")
Prints the log.
bool _buseSVDquiet
useSVD instead of Cholesky to factorize the input correlation matrix
Definition TSampling.h:80
void setUserCorrelation(URANIE::DataServer::TAttribute *x, URANIE::DataServer::TAttribute *y, Double_t value)
Defines a correlation between two attributes.
Definition TAMHCopula.h:60