English Français

Documentation / Developer's manual

Available modules

Calibration,  DataServer,  Launcher,  MetaModelOptim,  Modeler,  Optimizer,  ReLauncher,  Reliability,  ReOptimizer,  Sampler,  Sensitivity,  UncertModeler,  XmlProblem,   Uranie / Optimizer: TMultiGenCode.h Source File
Uranie / Optimizer  v4.10.0
/* @license-end */
TMultiGenCode.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 
47 #ifndef TMULTIGENCODE_H
48 #define TMULTIGENCODE_H
49 
50 // ROOT
51 #include "Math/IFunction.h"
52 #include "TH1.h"
53 #include "TF1.h"
54 #include "TFile.h"
55 #include "TNtupleD.h"
56 #include "TTreeFormula.h"
57 #include "TPad.h"
58 #include "TRandom3.h"
59 #include "TStyle.h"
60 
61 #include <vector>
62 #include <iostream>
63 
64 // URANIE
65 #include "Optimizer.h"
66 #include "TCode.h"
67 #include "TDataServer.h"
68 
69 namespace URANIE
70 {
71 namespace Optimizer
72 {
73 class TMultiGenCode: public ROOT::Math::IBaseFunctionMultiDim
74 {
75 
76 public:
77  //---------------------------------------------
81  TMultiGenCode(URANIE::DataServer::TDataServer *tds,
83  URANIE::Launcher::TCode *code, TString scost);
86 
87  //---------------------------------------------
91  void init();
94  void clean();
96 
97  //---------------------------------------------
101  TMultiGenCode * Clone() const
103  {
104  return new TMultiGenCode(_tds, _code, _sCost);
105  }
106 
108  unsigned int NDim() const
109  {
110  return _ninput;
111  }
113 
114  //---------------------------------------------
118  void setLog()
119  {
120  _blog = kTRUE;
121  }
122  void unsetLog()
123  {
124  _blog = kFALSE;
125  }
126 
127  void changeLog()
128  {
129  _blog = _blog ? kFALSE : kTRUE;
130  }
131  Bool_t getLog()
132  {
133  return _blog;
134  }
135 
137  {
138  _bintermedstep = true;
139  }
141 
142 
143 private:
144  //---------------------------------------------
148  double DoEval(const double * x) const;
151 
152  Bool_t _blog;
153  int _ninput;
154  TString _sinput;
155  double *_dinputValue;
156  double *_xin, *_codevalue;
157  int *_codeindex ;
158  int _noutput;
161  TString _soutput;
162  double * _doutputValue;
163  TString _sCost;
164  Int_t _nCost;
165  URANIE::DataServer::TDataServer *_tds;
166  URANIE::Launcher::TCode *_code;
167  TNtupleD *_informTree;
168  vector<TTreeFormula*> _vinForm;
170 
171  ClassDef(URANIE::Optimizer::TMultiGenCode, ID_OPTIMIZER)
172  //Classe de
173 
174 };
175 } // Fin du namespace Optimizer
176 } // Fin du namespace URANIE
177 
179 
180 #endif
Rosenbrock&#39;s function (n=2) with first and second order derivatives.
Definition: TBestEstimate.h:57
int _noutput
Definition: TMultiGenCode.h:158
unsigned int NDim() const
the NDIM method
Definition: TMultiGenCode.h:108
void unsetLog()
Definition: TMultiGenCode.h:122
void setLog()
Definition: TMultiGenCode.h:118
NEWMAT::ColumnVector & x
Definition: TOptimizerOpt.cxx:70
bool _bintermedstep
Definition: TMultiGenCode.h:169
void clean()
The clean method.
Definition: TMultiGenCode.h:73
TString _sinput
Name of the input attributes.
Definition: TMultiGenCode.h:154
TMultiGenCode(URANIE::DataServer::TDataServer *tds, URANIE::Launcher::TCode *code, TString scost)
Constructor with a dataserver.
TString _soutput
Name of the output attributes.
Definition: TMultiGenCode.h:161
URANIE::Launcher::TCode * _code
Pointeur vers un TDS.
Definition: TMultiGenCode.h:166
int _noutformula
Definition: TMultiGenCode.h:160
TString _sCost
The name of the selected cost (the number must be equal to 1)
Definition: TMultiGenCode.h:163
vector< TTreeFormula * > _vinForm
Definition: TMultiGenCode.h:168
int * _codeindex
Definition: TMultiGenCode.h:157
double * _codevalue
Definition: TMultiGenCode.h:156
double * _dinputValue
Definition: TMultiGenCode.h:155
int _ninput
Definition: TMultiGenCode.h:153
TNtupleD * _informTree
Pointeur vers un TCode.
Definition: TMultiGenCode.h:167
int _ninformula
Definition: TMultiGenCode.h:159
TMultiGenCode * Clone() const
The clone method.
Definition: TMultiGenCode.h:102
double * _xin
Definition: TMultiGenCode.h:156
R__EXTERN URANIE::Optimizer::TMultiGenCode * gTMGC
Definition: TMultiGenCode.h:178
void init()
The init method.
double DoEval(const double *x) const
the DoEval method
Bool_t _blog
Boolean for edit the log.
Definition: TMultiGenCode.h:152
Bool_t getLog()
Definition: TMultiGenCode.h:131
URANIE::DataServer::TDataServer * _tds
Definition: TMultiGenCode.h:165
void doNotSaveSteps()
Definition: TMultiGenCode.h:136
void changeLog()
Definition: TMultiGenCode.h:127
Int_t _nCost
the current position for the selected cost in the output list of attribute
Definition: TMultiGenCode.h:164
double * _doutputValue
The output attribute name.
Definition: TMultiGenCode.h:162